aboutsummaryrefslogtreecommitdiff
path: root/src/SetTmunu.F
blob: e75b14a6f85e7ff8d1e1a312158b74769e6b6376 (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
#include "cctk.h"
#include "cctk_Arguments.h"
#include "cctk_Functions.h"
#include "cctk_Parameters.h"


      
c     Calculate the contribution to the stress energy tensor T_munu
c     which are calcualated via the CalcTmunu.inc mechanism.
c     Then make a copy of that into the T2munu variables.

      subroutine TmunuBase_SetTmunu (CCTK_ARGUMENTS)
      implicit none
      DECLARE_CCTK_ARGUMENTS
      DECLARE_CCTK_FUNCTIONS
      DECLARE_CCTK_PARAMETERS
      
#define TMUNUBASE_SETTMUNU

#include "EinsteinBase/ADMMacros/src/macro/STRESSENERGY_declare.h"
      
      integer i, j, k
      
      do k = 1, cctk_lsh(3)
         do j = 1, cctk_lsh(2)
            do i = 1, cctk_lsh(1)
               
#include "EinsteinBase/ADMMacros/src/macro/STRESSENERGY_guts.h"

               eTtt(i,j,k) = Ttt
               
               eTtx(i,j,k) = Ttx
               eTty(i,j,k) = Tty
               eTtz(i,j,k) = Ttz
               
               eTxx(i,j,k) = Txx
               eTxy(i,j,k) = Txy
               eTxz(i,j,k) = Txz
               eTyy(i,j,k) = Tyy
               eTyz(i,j,k) = Tyz
               eTzz(i,j,k) = Tzz
               
#include "EinsteinBase/ADMMacros/src/macro/STRESSENERGY_undefine.h"
               
            end do
         end do
      end do
      
#undef TMUNUBASE_SETTMUNU
      
      eT2tt = eTtt
      
      eT2tx = eTtx
      eT2ty = eTty
      eT2tz = eTtz
      
      eT2xx = eTxx
      eT2xy = eTxy
      eT2xz = eTxz
      eT2yy = eTyy
      eT2yz = eTyz
      eT2zz = eTzz
      
      end