aboutsummaryrefslogtreecommitdiff
path: root/src/ADMConstraints.F
blob: c97828d320e6fc340a4b0e15d55fa8d818b3e839 (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
cc/*@@
c  @file      ADMConstraints.F
c  @date      Aug 98
c  @desc 
c     Calculate the ADM Constraints for output:
c    
c     Hamiltonian Constraint is:
c
c       R - K^i_j K^j_i + trK^2 - 16 Pi rho
c
c     Momentum Constraints are:
c
c       Del_j K_i^j - Del_i trK - 8 Pi j_i
c
c  @enddesc 
c@@*/

#include "cctk.h"
#include "cctk_Parameters.h"
#include "cctk_Arguments.h"
 
#include "CactusEinstein/Einstein/src/Einstein.h"

c/*@@
c  @routine    ADMConstraints
c  @date       Aug 98
c  @desc 
c    Calculate the ADM Constraints
c  @enddesc 
c  @calls     
c  @calledby   

c  @history 
c
c  @endhistory 
c@@*/

      subroutine ADMConstraints(CCTK_ARGUMENTS)
      
      implicit none

      DECLARE_CCTK_ARGUMENTS
      DECLARE_CCTK_PARAMETERS
      DECLARE_CCTK_FUNCTIONS

      integer :: i,j,k
      integer :: nx,ny,nz

c     Stencil width used for calculating constraints
c     (for outer boundary condition)

      integer, dimension(3),parameter :: sw = 1

c     Return code from Cactus boundary conditions.

      integer ierr

c     Various real variables.

      CCTK_REAL :: dx,dy,dz
      CCTK_REAL :: det,uxx,uxy,uxz,uyy,uyz,uzz
      CCTK_REAL :: m_rho,m_sx,m_sy,m_sz
      CCTK_REAL :: pi,ialp,ialp2

c     Temporaries for the Stress-Energy tensor.

      CCTK_REAL :: Ttt,Ttx,Tty,Ttz,Txx,Txy,Txz,Tyy,Tyz,Tzz

c     Matter declarations.

#include "CalcTmunu_temps.inc"

c     Macros from Standard Einstein.

#include "CactusEinstein/Einstein/src/macro/HAMADM_declare.h"
#include "CactusEinstein/Einstein/src/macro/MOMXADM_declare.h"
#include "CactusEinstein/Einstein/src/macro/MOMYADM_declare.h"
#include "CactusEinstein/Einstein/src/macro/MOMZADM_declare.h"
#include "CactusEinstein/Einstein/src/macro/DETG_declare.h"
#include "CactusEinstein/Einstein/src/macro/UPPERMET_declare.h"

c --------------------------------------------------------------

c     Grid parameters.

      dx = cctk_delta_space(1)
      dy = cctk_delta_space(2)
      dz = cctk_delta_space(3)

      nx = cctk_lsh(1)
      ny = cctk_lsh(2)
      nz = cctk_lsh(3)

c     Fill with zeros.

      ham = 0.0D0

      momx = 0.0D0
      momy = 0.0D0
      momz = 0.0D0

c     Calculate constraints.

      pi = acos(-1.0D0)

      do k=2,nz-1
         do j=2,ny-1
            do i=2,nx-1

               ialp  = 1.0D0/alp(i,j,k)
               ialp2 = ialp**2

c              Calculate the stress energy tensor at this point
c              ------------------------------------------------

c              Fill with zeroes.

               Ttt = 0.0D0; Ttx = 0.0D0; Tty = 0.0D0; Ttz = 0.0D0

               Txx = 0.0D0; Txy = 0.0D0; Txz = 0.0D0
               Tyy = 0.0D0; Tyz = 0.0D0; Tzz = 0.0D0

c              Inverse metric may be needed for CalcTmunu.

#include "CactusEinstein/Einstein/src/macro/DETG_guts.h"
#include "CactusEinstein/Einstein/src/macro/UPPERMET_guts.h"

               det = DETG_DETCG

               uxx = UPPERMET_UXX; uxy = UPPERMET_UXY; uxz = UPPERMET_UXZ
               uyy = UPPERMET_UYY; uyz = UPPERMET_UYZ; uzz = UPPERMET_UZZ

c              Call macro for stress energy tensor.

#include "CalcTmunu.inc"

c              Calculate the hamiltonian constraint
c              ------------------------------------

c              Geometric piece.

#include "CactusEinstein/Einstein/src/macro/HAMADM_guts.h"

c              Add matter terms: - 16*pi*rho
c
c              with rho defined as:
c
c              rho = n_a n_b T^{ab} = n^a n^b T_{ab}
c                  = (T_00 - 2 beta^i T_{i0} + beta^i beta^j T_{ij})/alpha^2

               m_rho = ialp2*Ttt


               if (shift_state == SHIFT_ACTIVE) then

                  m_rho = m_rho + ialp2
     &               *(betax(i,j,k)**2*Txx
     &               + betay(i,j,k)**2*Tyy
     &               + betaz(i,j,k)**2*Tzz
     &               +(betax(i,j,k)*betay(i,j,k)*Txy
     &               + betax(i,j,k)*betaz(i,j,k)*Txz
     &               + betay(i,j,k)*betaz(i,j,k)*Tyz)*2.0D0
     &               -(betax(i,j,k)**2*Ttx
     &               + betay(i,j,k)**2*Tty
     &               + betaz(i,j,k)**2*Ttz)*2.0D0)


               end if

               ham(i,j,k) = HAMADM_HAMADM - 16.0D0*pi*m_rho

c              From thorn MAHC (this overwrites all we just did).

#ifdef THORN_MAHC
               if ((some_hydro.eq.MAHCHYDRO).or.
     .             (some_hydro.eq.MAHCHYDRODUST)) then
                  ham(i,j,k) = ADM_ham(i,j,k) - 16.0d0*pi*
     .                 ((rho(i,j,k) + rho(i,j,k)*eps(i,j,k) + press(i,j,k))*
     .                 w_lorentz(i,j,k)**2 - press(i,j,k))
               end if
#endif

c              Calculate the Momentum constraints
c              ----------------------------------

c              Geometric piece.

#include "CactusEinstein/Einstein/src/macro/MOMXADM_guts.h"
#include "CactusEinstein/Einstein/src/macro/MOMYADM_guts.h"
#include "CactusEinstein/Einstein/src/macro/MOMZADM_guts.h"

c              Add matter terms: - 8*pi*S_i
c
c              with S_i defined as:
c
c              S_i = - g_{ia} n_b T^{ab} = - g_i^a n^b T_{ab} 
c                  = - (T_{i0} - beta^j T_{ij})/alpha

               m_sx = - ialp*Ttx

               m_sy = - ialp*Tty
               m_sz = - ialp*Ttz

               if (shift_state == SHIFT_ACTIVE) then

                  m_sx = m_sx + ialp
     &                 *(betax(i,j,k)*Txx
     &                 + betay(i,j,k)*Txy
     &                 + betaz(i,j,k)*Txz)

                  m_sy = m_sy + ialp
     &                 *(betax(i,j,k)*Txy
     &                 + betay(i,j,k)*Tyy
     &                 + betaz(i,j,k)*Tyz)

                  m_sz = m_sz + ialp
     &                 *(betax(i,j,k)*Txz
     &                 + betay(i,j,k)*Tyz
     &                 + betaz(i,j,k)*Tzz)

               end if

               momx(i,j,k) = MOMXADM_MOMXADM - 8.0D0*pi*m_sx
               momy(i,j,k) = MOMYADM_MOMYADM - 8.0D0*pi*m_sy
               momz(i,j,k) = MOMZADM_MOMZADM - 8.0D0*pi*m_sz

c              From thorn MAHC (this overwrites all we just did).

#ifdef THORN_MAHC
               if ((some_hydro.eq.MAHCHYDRO).or.
     .             (some_hydro.eq.MAHCHYDRODUST)) then
                  momx(i,j,k) = ADM_momx(i,j,k)
     .                 - 8.0d0*pi*sx(i,j,k)/sqrt(det)
                  momy(i,j,k) = ADM_momy(i,j,k)
     .                 - 8.0d0*pi*sy(i,j,k)/sqrt(det)
                  momz(i,j,k) = ADM_momz(i,j,k)
     .                 - 8.0d0*pi*sz(i,j,k)/sqrt(det)
               end if
#endif

            end do
         end do 
      end do


#include "CactusEinstein/Einstein/src/macro/DETG_undefine.h"
#include "CactusEinstein/Einstein/src/macro/UPPERMET_undefine.h"
#include "CactusEinstein/Einstein/src/macro/HAMADM_undefine.h"
#include "CactusEinstein/Einstein/src/macro/MOMXADM_undefine.h"
#include "CactusEinstein/Einstein/src/macro/MOMYADM_undefine.h"
#include "CactusEinstein/Einstein/src/macro/MOMZADM_undefine.h"

c     Apply symmetry boundary conditions.

      call CartSymGN(ierr,cctkGH,"admconstraints::admconstraints")

c     Apply flat boundary conditions at outer boundaries.

      if (CCTK_Equals(bound,"flat") == 1) then
         call BndFlatGN(ierr,cctkGH,sw,"admconstraints::admconstraints")
      end if

c     Synchronize.

      if (constraint_communication.eq.1) then
         call CCTK_SyncGroup(cctkGH,"admconstraints::admconstraints")
      end if

c     End

      end subroutine ADMConstraints