aboutsummaryrefslogtreecommitdiff
path: root/src/ADMConstraints.F
blob: c999936996a0cf8e652ede34703bfe8d5e5d5d32 (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
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@@*/

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

c CCTK Header files
#include "cctk.h"
#include "cctk_Parameters.h"
#include "cctk_Arguments.h"
 
#include "CactusEinstein/Einstein/src/Einstein.h"

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

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,i1,i2,j1,j2,k1,k2

c     Stencil width used for calculating constraints (for outer boundary condition)
      integer, dimension(3),parameter :: sw = 1

c     Return code from Cactus boundary conditions
      integer ierr

      CCTK_REAL :: dx,dy,dz

      CCTK_REAL :: pi,det,uxx,uxy,uxz,uyy,uyz,uzz

c     Temporaries for the Stress-Energy tensor
      CCTK_REAL :: Ttt,Ttx,Tty,Ttz,Txx,Txy,Txz,Tyy,Tyz,Tzz


#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"

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

      i1 = 2
      i2 = cctk_lsh(1)-1
      j1 = 2
      j2 = cctk_lsh(2)-1
      k1 = 2
      k2 = cctk_lsh(3)-1
  
      do i=i1,i2
        do j=j1,j2
          do k=k1,k2
            ham(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

      Pi =  ACOS(-1D0)

      do k = k1, k2
         do j = j1, j2
            do i = i1, i2 
               
c              Calculate the stress energy tensor at this point
c              ------------------------------------------------

               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              This may be needed for CalcTmunu etc
#include "CactusEinstein/Einstein/src/macro/DETG_guts.h"
               det = DETG_DETCG
#include "CactusEinstein/Einstein/src/macro/UPPERMET_guts.h"
               uxx = UPPERMET_UXX; uxy = UPPERMET_UXY; uxz = UPPERMET_UXZ
               uyy = UPPERMET_UYY; uyz = UPPERMET_UYZ; uzz = UPPERMET_UZZ
#include "CalcTmunu.inc"
              

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

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

               ham(i,j,k) = HAMADM_HAMADM - 16D0*Pi*Ttt/alp(i,j,k)**2

#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))
               endif
#endif


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

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

               momx(i,j,k) = MOMXADM_MOMXADM - 8D0*Pi*Ttx/alp(i,j,k)**2
               momy(i,j,k) = MOMYADM_MOMYADM - 8D0*Pi*Tty/alp(i,j,k)**2
               momz(i,j,k) = MOMZADM_MOMZADM - 8D0*Pi*Ttz/alp(i,j,k)**2


#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)
               endif
#endif
               
            enddo
         enddo
      enddo

#include "CactusEinstein/Einstein/src/macro/DETG_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"
#include "CactusEinstein/Einstein/src/macro/UPPERMET_undefine.h"


c     Synchronize and apply flat boundary conditions

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

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

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

      end subroutine ADMConstraints