aboutsummaryrefslogtreecommitdiff
path: root/src/Utils.F90
blob: 9f670e25b9bb05726c591b7f81571d91049e3c3e (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
 /*@@
   @file      Utils.F
   @date      Sat Jan 26 02:28:46 2002
   @author    
   @desc 
   Utility functions for other thorns. Calculation of the determinant
   of the spatial metric and the upper metric.
   @enddesc 
 @@*/

#include "cctk.h"
#include "cctk_Parameters.h"
#include "cctk_Arguments.h"
#include "GRHydro_Macros.h"

 /*@@
   @routine    GRHydro_RefinementLevel
   @date       July 2005
   @author     
   @desc 
   Calculates the current refinement level from flesh data
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 
@@*/

subroutine GRHydro_Debug(CCTK_ARGUMENTS)

  implicit none
  DECLARE_CCTK_ARGUMENTS
  integer i,j,k
  integer nx, ny, nz

  GRHydro_reflevel = int(log10(dble(cctk_levfac(1)))/log10(2.0d0))
  
  nx = cctk_lsh(1)
  ny = cctk_lsh(2)
  nz = cctk_lsh(3)

  if (GRHydro_reflevel .lt. 3) return
 
  do i=4,nx
     do j=4,ny
        do k=4,4
           if( r(i,j,k)-4.0d0 .lt. 96.0d0 .and. &
                r(i,j,k)+4.0d0 .gt. 96.0d0) then
              write(6,"(4i4,1P10E15.6)") GRHydro_reflevel, i,j,k,&
                   r(i,j,k),rho(i,j,k),eps(i,j,k),y_e(i,j,k),temperature(i,j,k)
           endif
        enddo
     enddo
  enddo

!  call flush(6)
!  if(GRHydro_reflevel.eq.4.and.temperature(1,1,1).lt.0.1d0) then
!     call CCTK_WARN(0,"stop")
!  endif

end subroutine GRHydro_Debug


subroutine GRHydro_RefinementLevel(CCTK_ARGUMENTS)

  implicit none
  DECLARE_CCTK_ARGUMENTS

  GRHydro_reflevel = int(log10(dble(cctk_levfac(1)))/log10(2.0d0))

end subroutine GRHydro_RefinementLevel

subroutine GRHydro_SqrtSpatialDeterminant(CCTK_ARGUMENTS)

  implicit none
  DECLARE_CCTK_ARGUMENTS
  integer i,j,k
  integer nx, ny, nz

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

  !$OMP PARALLEL DO
  do k=1,nz
     do j=1,ny
        do i=1,nx
           sdetg(i,j,k) = -(gxz(i,j,k)**2)*gyy(i,j,k) + & 
                2.0d0*gxy(i,j,k)*gxz(i,j,k)*gyz(i,j,k) - &
                gxx(i,j,k)*(gyz(i,j,k)**2) - &
                (gxy(i,j,k)**2)*gzz(i,j,k)  + &
                (gxx(i,j,k)*gyy(i,j,k))*gzz(i,j,k)
           sdetg(i,j,k) = sqrt(sdetg(i,j,k))
        enddo
     enddo
  enddo
  !$OMP END PARALLEL DO

end subroutine GRHydro_SqrtSpatialDeterminant


subroutine SpatialDeterminant(gxx,gxy,gxz,gyy,gyz,gzz,det)

  implicit none
  
  CCTK_REAL :: det
  CCTK_REAL :: gxx,gxy,gxz,gyy,gyz,gzz

  det = -gxz**2*gyy + 2*gxy*gxz*gyz - gxx*gyz**2 - gxy**2*gzz + gxx*gyy*gzz


!!$  Why is this weird order necessary? Search me. It just seemed 
!!$  to make a really odd NaN go away.

!!$  det =  -gxz**2*gyy + 2*gxy*gxz*gyz 
!!$  det = det - gxx*gyz**2 - gxy**2*gzz
!!$  det = det + gxx*gyy*gzz
  
  return
  
end subroutine SpatialDeterminant



/*@@
@routine    UpperMetric
@date       Sat Jan 26 02:32:26 2002
@author     
@desc 
Calculates the upper metric. The determinant is given, not
calculated.
@enddesc 
@calls     
@calledby   
@history 

@endhistory 
@@*/

subroutine UpperMetric(uxx, uxy, uxz, uyy, uyz, uzz, &
     det, gxx, gxy, gxz, gyy, gyz, gzz)
  
  implicit none
  
  CCTK_REAL :: uxx, uxy, uxz, uyy, uyz, uzz, det, &
       gxx, gxy, gxz, gyy, gyz, gzz, invdet
  
  invdet = 1.d0 / det
  uxx=(-gyz**2 + gyy*gzz)*invdet
  uxy=(gxz*gyz - gxy*gzz)*invdet
  uyy=(-gxz**2 + gxx*gzz)*invdet
  uxz=(-gxz*gyy + gxy*gyz)*invdet
  uyz=(gxy*gxz - gxx*gyz)*invdet
  uzz=(-gxy**2 + gxx*gyy)*invdet
  
  return
  
end subroutine UpperMetric