aboutsummaryrefslogtreecommitdiff
path: root/src/ADMConstraints.F
blob: 090da3cfa36e1f29e355b895add945d1c9fdc05c (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
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
/*@@
  @file      ADMConstraints.F77
  @date      August 98
  @desc 
             Calculate the ADM Constraints for output:
    
             Hamiltonian Constraint is:

             H = R - K^i_j K^j_i + trK^2 - 16 Pi rho

             Momentum Constraints are:

             M_i = Del_j K_i^j - Del_i trK - 8 Pi S_i

  @enddesc
  @version $Header$
@@*/

#include "cctk.h"

#include "cctk_Arguments.h"
#include "cctk_Functions.h"
#include "cctk_Parameters.h"

#include "SpaceMask.h"

      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 sw(3)
      logical docalc

c     flags for excision

      integer ex_field, ex_type_excised

c     Various real variables.

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

#include "EinsteinBase/ADMMacros/src/macro/ADM_Spacing_declare.h"

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 "EinsteinBase/ADMMacros/src/macro/HAMADM_declare.h"
#include "EinsteinBase/ADMMacros/src/macro/MOMXADM_declare.h"
#include "EinsteinBase/ADMMacros/src/macro/MOMYADM_declare.h"
#include "EinsteinBase/ADMMacros/src/macro/MOMZADM_declare.h"
#include "EinsteinBase/ADMMacros/src/macro/DETG_declare.h"
#include "EinsteinBase/ADMMacros/src/macro/UPPERMET_declare.h"

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

      if (check_excision_bitmask .ne. 0) then
          call SpaceMask_GetTypeBits(ex_field, excision_mask_name)
          call SpaceMask_GetStateBits(ex_type_excised, excision_mask_name, \
                                      excision_type_excised)
      end if

      sw(1) = 1
      sw(2) = 1
      sw(3) = 1

c     Grid parameters.

#include "EinsteinBase/ADMMacros/src/macro/ADM_Spacing.h"

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

c     Fill with zeros.

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

              ham(i,j,k) = 0.0D0
              hamnormalized(i,j,k) = 0.0D0
        
              momx(i,j,k) = 0.0D0
              momy(i,j,k) = 0.0D0
              momz(i,j,k) = 0.0D0       

            end do
         end do
      end do

c     Calculate constraints.

      pi = acos(-1.0D0)

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

              if ( (i<3).or.(i>cctk_lsh(1)-2).or.
     .             (j<3).or.(j>cctk_lsh(2)-2).or.
     .             (k<3).or.(k>cctk_lsh(3)-2) ) then
                local_spatial_order = 2
              else
                local_spatial_order = spatial_order
              end if

               docalc = .TRUE.

               if (use_mask .eq. 1) then
                  if (abs(emask(i,j,k)-1) > 0.001) then
                     docalc = .FALSE.
                  elseif (check_excision_bitmask .ne. 0) then
                     if (SpaceMask_CheckStateBitsF90(space_mask, i, j, k,\
                          ex_field, ex_type_excised)) then
                        docalc = .FALSE.
                     end if
                  end if
               end if
                  
               if (docalc) then

                  ialp  = 1.0D0/alp(i,j,k)
                  ialp2 = ialp**2
                  
c                 Calculate the stress energy tensor at this point
c                 ------------------------------------------------

c                 This may be needed for CalcTmunu

#include "EinsteinBase/ADMMacros/src/macro/DETG_guts.h"

                  det = DETG_DETCG

#include "EinsteinBase/ADMMacros/src/macro/UPPERMET_guts.h"

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

c                 Initialize stress-energy tensor components.

                  if (stress_energy_state .ne. 0) then

                     Ttt = eTtt(i,j,k)

                     Ttx = eTtx(i,j,k); Tty = eTty(i,j,k); Ttz = eTtz(i,j,k)

                     Txx = eTxx(i,j,k); Tyy = eTyy(i,j,k); Tzz = eTzz(i,j,k)

                     Txy = eTxy(i,j,k); Txz = eTxz(i,j,k); Tyz = eTyz(i,j,k)

                  else

                     Ttt = 0.0D0

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

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

c                    Include macro for stress energy tensor.

#include "CalcTmunu.inc"
                  
                  end if

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

c                 Geometric piece.

#include "EinsteinBase/ADMMacros/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 .eq. 1) 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)*Ttx
     &                  + betay(i,j,k)*Tty
     &                  + betaz(i,j,k)*Ttz)*2.0D0)

                  end if

                  ham(i,j,k) = HAMADM_HAMADM - 16.0D0*pi*m_rho
                  if ((HAMADM_HAMADMABS + abs(16.0D0*pi*m_rho))==0) then
                     hamnormalized(i,j,k) = abs(ham(i,j,k))
                  else
                     hamnormalized(i,j,k) = abs(ham(i,j,k))/
     &              (HAMADM_HAMADMABS + abs(16.0D0*pi*m_rho))
                  end if


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

c                 Geometric piece.

#include "EinsteinBase/ADMMacros/src/macro/MOMXADM_guts.h"
#include "EinsteinBase/ADMMacros/src/macro/MOMYADM_guts.h"
#include "EinsteinBase/ADMMacros/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 .eq. 1) 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

               else

                  ham(i,j,k) = excised_value
                  momx(i,j,k) = excised_value
                  momy(i,j,k) = excised_value
                  momz(i,j,k) = excised_value

               end if

            end do
         end do 
      end do

#include "EinsteinBase/ADMMacros/src/macro/DETG_undefine.h"
#include "EinsteinBase/ADMMacros/src/macro/UPPERMET_undefine.h"
#include "EinsteinBase/ADMMacros/src/macro/HAMADM_undefine.h"
#include "EinsteinBase/ADMMacros/src/macro/MOMXADM_undefine.h"
#include "EinsteinBase/ADMMacros/src/macro/MOMYADM_undefine.h"
#include "EinsteinBase/ADMMacros/src/macro/MOMZADM_undefine.h"

      return
      end

      

      subroutine ADMConstraints_Boundaries(CCTK_ARGUMENTS)
      
      implicit none

      DECLARE_CCTK_ARGUMENTS
      DECLARE_CCTK_PARAMETERS
      DECLARE_CCTK_FUNCTIONS

      CCTK_INT, parameter :: izero = 0
      integer,  parameter :: ik = kind (izero)

c     Return code from Cactus sync routine and boundary conditions.

      integer ierr

c     Apply flat boundary conditions at outer boundaries.

      if (CCTK_EQUALS(bound,"flat")) then
         ierr = Boundary_SelectGroupForBC(cctkGH, int(CCTK_ALL_FACES,ik), 1_ik, -1_ik,
     $        "admconstraints::hamiltonian", "Flat")
         ierr = Boundary_SelectGroupForBC(cctkGH, int(CCTK_ALL_FACES,ik), 1_ik, -1_ik,
     $        "admconstraints::normalized_hamiltonian", "Flat")
         ierr = Boundary_SelectGroupForBC(cctkGH, int(CCTK_ALL_FACES,ik), 1_ik, -1_ik,
     $        "admconstraints::momentum", "Flat")
      else
         ierr = Boundary_SelectGroupForBC(cctkGH, int(CCTK_ALL_FACES,ik), 1_ik, -1_ik,
     $        "admconstraints::hamiltonian", "None")
         ierr = Boundary_SelectGroupForBC(cctkGH, int(CCTK_ALL_FACES,ik), 1_ik, -1_ik,
     $        "admconstraints::normalized_hamiltonian", "None")
         ierr = Boundary_SelectGroupForBC(cctkGH, int(CCTK_ALL_FACES,ik), 1_ik, -1_ik,
     $        "admconstraints::momentum", "None")
      end if

      return
      end