aboutsummaryrefslogtreecommitdiff
path: root/src/include/GRHydro_EMTensor.inc
blob: 4318c16ac6d86852dd6b77e9246506e273801cdf (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
 /*@@
   @file      GRHydro_EMTensor.inc
   @date      Sat Jan 26 02:03:56 2002
   @author    Ian Hawke
   @desc 
     The calculation of the stress energy tensor.
     The version used here was worked out by Miguel Alcubierre. I
     think it was an extension of the routine from GR3D, written
     by Mark Miller.
     C version added by Ian Hawke.

     Lower components of the stress-energy tensor obtained from
     the hydro variables.  The components are given by:

     T      =  rho h  u   u    +  P  g
      mu nu            mu  nu         mu nu

     where rho is the energy density of the fluid, h the enthalpy
     and P the pressure.  The enthalpy is given in terms of the
     basic variables as:

     h  =  1  +  e  +  P/rho

     with e the internal energy (eps here).

     In the expresion for T_{mu,nu} we also have the four-velocity
     of the fluid given by (v_i is the 3-velocity field):

                                 i
     u  =  W ( - alpha +  v  beta  )
      0                    i

     u  =  W v
      i       i
                                                i  -1/2
     with W the Lorentz factor:   W = ( 1 -  v v  )
                                              i

     and where alpha and beta are the lapse and shift vector.

     Finally, the 4 metric is given by

                   2             i
     g   =  - alpha  + beta  beta
      00                   i

     g   =  beta
      0i        i


     g   =  gamma      (the spatial metric)
      ij        ij


   @enddesc 
 @@*/

#ifdef FCODE

/*       if (SpaceMask_CheckStateBitsF90(space_mask, i,j,k, \
                                       atmosphere_field_descriptor, \
                                       atmosphere_normal_descriptor)) then*/
!!$      call SpaceMask_GetTypeBits(type_bits, "Hydro_Atmosphere")
!!$      call SpaceMask_GetStateBits(atmosphere, "Hydro_Atmosphere", "in_atmosphere")
       if (conformal_state .eq. 0) then
         velxlow = gxx(i,j,k)*velx(i,j,k) + gxy(i,j,k)*vely(i,j,k) + gxz(i,j,k)*velz(i,j,k)
         velylow = gxy(i,j,k)*velx(i,j,k) + gyy(i,j,k)*vely(i,j,k) + gyz(i,j,k)*velz(i,j,k)
         velzlow = gxz(i,j,k)*velx(i,j,k) + gyz(i,j,k)*vely(i,j,k) + gzz(i,j,k)*velz(i,j,k)
       else
         velxlow = psi(i,j,k)**4 * (gxx(i,j,k)*velx(i,j,k) + gxy(i,j,k)*vely(i,j,k) + gxz(i,j,k)*velz(i,j,k) )
         velylow = psi(i,j,k)**4 * (gxy(i,j,k)*velx(i,j,k) + gyy(i,j,k)*vely(i,j,k) + gyz(i,j,k)*velz(i,j,k) )
         velzlow = psi(i,j,k)**4 * (gxz(i,j,k)*velx(i,j,k) + gyz(i,j,k)*vely(i,j,k) + gzz(i,j,k)*velz(i,j,k) )
       end if

!!$       Calculate lower components and square of shift vector.

       if (shift_state .ne. 0) then

         if (conformal_state .eq. 0) then
           betaxlow = gxx(i,j,k)*betax(i,j,k) + gxy(i,j,k)*betay(i,j,k) + gxz(i,j,k)*betaz(i,j,k)
           betaylow = gxy(i,j,k)*betax(i,j,k) + gyy(i,j,k)*betay(i,j,k) + gyz(i,j,k)*betaz(i,j,k)
           betazlow = gxz(i,j,k)*betax(i,j,k) + gyz(i,j,k)*betay(i,j,k) + gzz(i,j,k)*betaz(i,j,k)
         else
           betaxlow = psi(i,j,k)**4 * (gxx(i,j,k)*betax(i,j,k) + gxy(i,j,k)*betay(i,j,k) + gxz(i,j,k)*betaz(i,j,k) )
           betaylow = psi(i,j,k)**4 * (gxy(i,j,k)*betax(i,j,k) + gyy(i,j,k)*betay(i,j,k) + gyz(i,j,k)*betaz(i,j,k) )
           betazlow = psi(i,j,k)**4 * (gxz(i,j,k)*betax(i,j,k) + gyz(i,j,k)*betay(i,j,k) + gzz(i,j,k)*betaz(i,j,k) )
         end if

          beta2 = betax(i,j,k)*betaxlow + betay(i,j,k)*betaylow + betaz(i,j,k)*betazlow 

       else

          betaxlow = 0.0D0
          betaylow = 0.0D0
          betazlow = 0.0D0

          beta2 = 0.0D0

       end if

!!$       Calculate the specific relativistic enthalpy times rho times the
!!$       square of the lorentz factor.

       rhoenthalpy = w_lorentz(i,j,k)**2*(rho(i,j,k)*(1.0D0 + eps(i,j,k)) + press(i,j,k))

!!$       Calculate lower components of 4-velocity (without the Lorent factor).

       utlow = (-alp(i,j,k) + velx(i,j,k)*betaxlow + vely(i,j,k)*betaylow + velz(i,j,k)*betazlow)

       uxlow = velxlow
       uylow = velylow
       uzlow = velzlow

!!$       Calculate Tmunu (the lower components!).

!!$      if ( .not.(SpaceMask_CheckStateBitsF90(space_mask, i, j, k, type_bits, atmosphere) )) then

         Ttt = Ttt + rhoenthalpy*utlow**2 + press(i,j,k)*(beta2 - alp(i,j,k)**2)

         Ttx = Ttx + rhoenthalpy*utlow*uxlow + press(i,j,k)*betaxlow
         Tty = Tty + rhoenthalpy*utlow*uylow + press(i,j,k)*betaylow
         Ttz = Ttz + rhoenthalpy*utlow*uzlow + press(i,j,k)*betazlow

         if (conformal_state .eq. 0) then

           Txx = Txx + rhoenthalpy*uxlow**2 + press(i,j,k)*gxx(i,j,k)
           Tyy = Tyy + rhoenthalpy*uylow**2 + press(i,j,k)*gyy(i,j,k)
           Tzz = Tzz + rhoenthalpy*uzlow**2 + press(i,j,k)*gzz(i,j,k)

           Txy = Txy + rhoenthalpy*uxlow*uylow + press(i,j,k)*gxy(i,j,k)
           Txz = Txz + rhoenthalpy*uxlow*uzlow + press(i,j,k)*gxz(i,j,k)
           Tyz = Tyz + rhoenthalpy*uylow*uzlow + press(i,j,k)*gyz(i,j,k)

         else

           Txx = Txx + rhoenthalpy*uxlow**2 + psi(i,j,k)**4 * press(i,j,k)*gxx(i,j,k)
           Tyy = Tyy + rhoenthalpy*uylow**2 + psi(i,j,k)**4 * press(i,j,k)*gyy(i,j,k)
           Tzz = Tzz + rhoenthalpy*uzlow**2 + psi(i,j,k)**4 * press(i,j,k)*gzz(i,j,k)

           Txy = Txy + rhoenthalpy*uxlow*uylow + psi(i,j,k)**4 * press(i,j,k)*gxy(i,j,k)
           Txz = Txz + rhoenthalpy*uxlow*uzlow + psi(i,j,k)**4 * press(i,j,k)*gxz(i,j,k)
           Tyz = Tyz + rhoenthalpy*uylow*uzlow + psi(i,j,k)**4 * press(i,j,k)*gyz(i,j,k)

         end if
!!$      end if
/*       end if*/

#endif

#ifdef CCODE

       if (!*conformal_state)
       {
         velxlow = gxx[ijk]*velx[ijk] + gxy[ijk]*vely[ijk] + 
                   gxz[ijk]*velz[ijk];
         velylow = gxy[ijk]*velx[ijk] + gyy[ijk]*vely[ijk] + 
                   gyz[ijk]*velz[ijk];
         velzlow = gxz[ijk]*velx[ijk] + gyz[ijk]*vely[ijk] + 
                   gzz[ijk]*velz[ijk];
       }
       else
       {
         velxlow = pow(psi[ijk],4) * (gxx[ijk]*velx[ijk] + 
                                      gxy[ijk]*vely[ijk] + 
                                      gxz[ijk]*velz[ijk]);
         velylow = pow(psi[ijk],4) * (gxy[ijk]*velx[ijk] + 
                                      gyy[ijk]*vely[ijk] + 
                                      gyz[ijk]*velz[ijk]);
         velzlow = pow(psi[ijk],4) * (gxz[ijk]*velx[ijk] + 
                                      gyz[ijk]*vely[ijk] + 
                                      gzz[ijk]*velz[ijk]);
       }

/*
 *        Calculate lower components and square of shift vector.
 */

       if (*shift_state)
       {
         if (!*conformal_state)
         {
           betaxlow = gxx[ijk]*betax[ijk] + gxy[ijk]*betay[ijk] 
                    + gxz[ijk]*betaz[ijk];
           betaylow = gxy[ijk]*betax[ijk] + gyy[ijk]*betay[ijk]
                    + gyz[ijk]*betaz[ijk];
           betazlow = gxz[ijk]*betax[ijk] + gyz[ijk]*betay[ijk]
                    + gzz[ijk]*betaz[ijk];
         }
         else
         {
           betaxlow = pow(psi[ijk],4) * (gxx[ijk]*betax[ijk] + 
                                         gxy[ijk]*betay[ijk] +
                                         gxz[ijk]*betaz[ijk]);
           betaylow = pow(psi[ijk],4) * (gxy[ijk]*betax[ijk] + 
                                         gyy[ijk]*betay[ijk] +
                                         gyz[ijk]*betaz[ijk]);
           betazlow = pow(psi[ijk],4) * (gxz[ijk]*betax[ijk] + 
                                         gyz[ijk]*betay[ijk] +
                                         gzz[ijk]*betaz[ijk]);
         }
         beta2 = betax[ijk]*betaxlow + betay[ijk]*betaylow
                + betaz[ijk]*betazlow ;
       }
       else
       {
          betaxlow = 0.0;
          betaylow = 0.0;
          betazlow = 0.0;

          beta2 = 0.0;
       }

/*
 *        Calculate the specific relativistic enthalpy times rho times the
 *        square of the lorentz factor.
 */

       rhoenthalpy = w_lorentz[ijk]*w_lorentz[ijk]*
                      (rho[ijk]*(1.0 + eps[ijk]) + press[ijk]);

/*
 *        Calculate lower components of 4-velocity (without the Lorent factor).
 */
       utlow = (-alp[ijk] + velx[ijk]*betaxlow + vely[ijk]*betaylow
                          + velz[ijk]*betazlow);

       uxlow = velxlow;
       uylow = velylow;
       uzlow = velzlow;

/*
 *        Calculate Tmunu (the lower components!).
 */

       Ttt = Ttt + rhoenthalpy*utlow*utlow +
             press[ijk]*(beta2 - alp[ijk]*alp[ijk]);

       Ttx = Ttx + rhoenthalpy*utlow*uxlow + press[ijk]*betaxlow;
       Tty = Tty + rhoenthalpy*utlow*uylow + press[ijk]*betaylow;
       Ttz = Ttz + rhoenthalpy*utlow*uzlow + press[ijk]*betazlow;

       if (!*conformal_state)
       {
         Txx = Txx + rhoenthalpy*uxlow*uxlow + press[ijk]*gxx[ijk];
         Tyy = Tyy + rhoenthalpy*uylow*uylow + press[ijk]*gyy[ijk];
         Tzz = Tzz + rhoenthalpy*uzlow*uzlow + press[ijk]*gzz[ijk];

         Txy = Txy + rhoenthalpy*uxlow*uylow + press[ijk]*gxy[ijk];
         Txz = Txz + rhoenthalpy*uxlow*uzlow + press[ijk]*gxz[ijk];
         Tyz = Tyz + rhoenthalpy*uylow*uzlow + press[ijk]*gyz[ijk];
       }
       else
       {
         Txx = Txx + rhoenthalpy*uxlow*uxlow + 
                     pow(psi[ijk],4) * press[ijk]*gxx[ijk];
         Tyy = Tyy + rhoenthalpy*uylow*uylow + 
                     pow(psi[ijk],4) * press[ijk]*gyy[ijk];
         Tzz = Tzz + rhoenthalpy*uzlow*uzlow + 
                     pow(psi[ijk],4) * press[ijk]*gzz[ijk];

         Txy = Txy + rhoenthalpy*uxlow*uylow + 
                     pow(psi[ijk],4) * press[ijk]*gxy[ijk];
         Txz = Txz + rhoenthalpy*uxlow*uzlow + 
                     pow(psi[ijk],4) * press[ijk]*gxz[ijk];
         Tyz = Tyz + rhoenthalpy*uylow*uzlow + 
                     pow(psi[ijk],4) * press[ijk]*gyz[ijk];
       }

#endif