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
|
! Calculation of the sources for the level set function.
! $Header$
#include "cctk.h"
#include "cctk_Parameters.h"
#include "cctk_Arguments.h"
subroutine EHFinder_Sources(CCTK_ARGUMENTS)
use EHFinder_mod
implicit none
DECLARE_CCTK_PARAMETERS
DECLARE_CCTK_ARGUMENTS
DECLARE_CCTK_FUNCTIONS
CCTK_INT :: i, j, k
CCTK_REAL :: idx, idy, idz, mdelta
CCTK_REAL :: a, b, c
CCTK_REAL :: g4tt, g4tx, g4ty, g4tz, g4xx, g4xy, g4xz, g4yy, g4yz, g4zz
CCTK_REAL :: g3xx, g3xy, g3xz, g3yy, g3yz, g3zz
CCTK_REAL :: gxxc, gxyc, gxzc, gyyc, gyzc, gzzc, psito4
CCTK_REAL :: idetg, alp2, tmp1, tmp2, tmp3
CCTK_REAL :: ratio, cfactor
CCTK_REAL, dimension(3) :: maxpos, cdx, dfup
CCTK_REAL :: al, ar, bl, br, cl, cr
CCTK_REAL :: alminus, alplus, blminus, blplus, clminus, clplus
CCTK_REAL :: arminus, arplus, brminus, brplus, crminus, crplus
! calculate 1/(2*delta) in each direction
#include "include/physical_part.h"
idx = half / cctk_delta_space(1)
idy = half / cctk_delta_space(2)
idz = half / cctk_delta_space(3)
mdelta = maxval ( cctk_delta_space )
if ( CCTK_EQUALS ( upwind_type, 'intrinsic' ) ) then
do k = kzl, kzr
do j = jyl, jyr
do i = ixl, ixr
if ( f(i,j,k) .gt. fmin_bound - three*mdelta ) then
#include "include/centered_second2.h"
else
#include "include/upwind_second2.h"
end if
end do
end do
end do
end if
if ( CCTK_EQUALS ( upwind_type, 'shift' ) ) then
do k = kzl, kzr
do j = jyl, jyr
do i = ixl, ixr
#include "include/upwind_shift_second2.h"
end do
end do
end do
end if
if ( CCTK_EQUALS ( upwind_type, 'characteristic' ) ) then
do k = kzl, kzr
do j = jyl, jyr
do i = ixl, ixr
#include "include/metric.h"
#include "include/upwind_characteristic_second2.h"
end do
end do
end do
end if
if ( CCTK_EQUALS ( metric_type, 'physical' ) ) then
do k = kzl, kzr
do j = jyl, jyr
do i = ixl, ixr
if ( eh_mask(i,j,k) .ge. 0 ) then
alp2 = alp(i,j,k)**2
g4tt = -one
g4tx = betax(i,j,k)
g4ty = betay(i,j,k)
g4tz = betaz(i,j,k)
tmp1 = gyy(i,j,k)*gzz(i,j,k) - gyz(i,j,k)**2
tmp2 = gxz(i,j,k)*gyz(i,j,k) - gxy(i,j,k)*gzz(i,j,k)
tmp3 = gxy(i,j,k)*gyz(i,j,k) - gxz(i,j,k)*gyy(i,j,k)
idetg = one / ( gxx(i,j,k)*tmp1 + gxy(i,j,k)*tmp2 + gxz(i,j,k)*tmp3 )
g4xx = tmp1 * idetg
g4xy = tmp2 * idetg
g4xz = tmp3 * idetg
g4yy = ( gxx(i,j,k)*gzz(i,j,k) - gxz(i,j,k)**2 ) * idetg
g4yz = ( gxy(i,j,k)*gxz(i,j,k) - gxx(i,j,k)*gyz(i,j,k) ) * idetg
g4zz = ( gxx(i,j,k)*gyy(i,j,k) - gxy(i,j,k)**2 ) * idetg
g4xx = alp2*g4xx - betax(i,j,k)**2
g4xy = alp2*g4xy - betax(i,j,k)*betay(i,j,k)
g4xz = alp2*g4xz - betax(i,j,k)*betaz(i,j,k)
g4yy = alp2*g4yy - betay(i,j,k)**2
g4yz = alp2*g4yz - betay(i,j,k)*betaz(i,j,k)
g4zz = alp2*g4zz - betaz(i,j,k)**2
a = g4tt
b = two * ( g4tx * dfx(i,j,k) + &
g4ty * dfy(i,j,k) + &
g4tz * dfz(i,j,k) )
c = g4xx * dfx(i,j,k)**2 + &
g4yy * dfy(i,j,k)**2 + &
g4zz * dfz(i,j,k)**2 + &
two * ( g4xy * dfx(i,j,k) * dfy(i,j,k) + &
g4xz * dfx(i,j,k) * dfz(i,j,k) + &
g4yz * dfy(i,j,k) * dfz(i,j,k) )
if ( b .lt. zero ) then
sf(i,j,k) = ( -b + sqrt(b**2 - four*a*c) ) / ( two*a )
else
sf(i,j,k) = two*c / ( -b - sqrt(b**2 - 4*a*c) )
end if
else
sf(i,j,k) = zero
end if
end do
end do
end do
end if
if ( CCTK_EQUALS ( metric_type, 'static conformal' ) ) then
do k = kzl, kzr
do j = jyl, jyr
do i = ixl, ixr
if ( eh_mask(i,j,k) .ge. 0 ) then
alp2 = alp(i,j,k)**2
g4tt = -one
g4tx = betax(i,j,k)
g4ty = betay(i,j,k)
g4tz = betaz(i,j,k)
psito4 = psi(i,j,k)**4
gxxc = gxx(i,j,k) * psito4
gxyc = gxy(i,j,k) * psito4
gxzc = gxz(i,j,k) * psito4
gyyc = gyy(i,j,k) * psito4
gyzc = gyz(i,j,k) * psito4
gzzc = gzz(i,j,k) * psito4
tmp1 = gyyc*gzzc - gyzc**2
tmp2 = gxzc*gyzc - gxyc*gzzc
tmp3 = gxyc*gyzc - gxzc*gyyc
idetg = one / ( gxxc*tmp1 + gxyc*tmp2 + gxzc*tmp3 )
g4xx = tmp1 * idetg
g4xy = tmp2 * idetg
g4xz = tmp3 * idetg
g4yy = ( gxxc*gzzc - gxzc**2 ) * idetg
g4yz = ( gxyc*gxzc - gxxc*gyzc ) * idetg
g4zz = ( gxxc*gyyc - gxyc**2 ) * idetg
g4xx = alp2*g4xx - betax(i,j,k)**2
g4xy = alp2*g4xy - betax(i,j,k)*betay(i,j,k)
g4xz = alp2*g4xz - betax(i,j,k)*betaz(i,j,k)
g4yy = alp2*g4yy - betay(i,j,k)**2
g4yz = alp2*g4yz - betay(i,j,k)*betaz(i,j,k)
g4zz = alp2*g4zz - betaz(i,j,k)**2
a = g4tt
b = two * ( g4tx * dfx(i,j,k) + &
g4ty * dfy(i,j,k) + &
g4tz * dfz(i,j,k) )
c = g4xx * dfx(i,j,k)**2 + &
g4yy * dfy(i,j,k)**2 + &
g4zz * dfz(i,j,k)**2 + &
two * ( g4xy * dfx(i,j,k) * dfy(i,j,k) + &
g4xz * dfx(i,j,k) * dfz(i,j,k) + &
g4yz * dfy(i,j,k) * dfz(i,j,k) )
if ( b .le. zero ) then
sf(i,j,k) = ( -b + sqrt(b**2 - four*a*c) ) / ( two*a )
else
sf(i,j,k) = two*c / ( -b - sqrt(b**2 - 4*a*c) )
end if
sf(i,j,k) = sf(i,j,k)*sign(one,alp(i,j,k))
! if ( i .eq. 41 .and. j .eq. 20 .and. k .eq. 25 ) then
! print*,i,j,k
! print*,f(i,j,k),eh_mask(i,j,k)
! print*,f(i-1,j,k),f(i+1,j,k),eh_mask(i-1,j,k),eh_mask(i+1,j,k)
! print*,f(i,j-1,k),f(i,j+1,k),eh_mask(i,j-1,k),eh_mask(i,j+1,k)
! print*,f(i,j,k-1),f(i,j,k+1),eh_mask(i,j,k-1),eh_mask(i,j,k+1)
! print*,dfx(i,j,k),dfy(i,j,k),dfz(i,j,k)
! print*
! print*,x(i,j,k),y(i,j,k),z(i,j,k)
! print*,alp(i,j,k)
! print*,betax(i,j,k),betay(i,j,k),betaz(i,j,k)
! print*,gxx(i,j,k),gxy(i,j,k),gxz(i,j,k)
! print*,gyy(i,j,k),gyz(i,j,k),gzz(i,j,k)
! print*,psi(i,j,k)
! print*
! print*,tmp1,tmp2,tmp3
! print*,gxxc*tmp1 + gxyc*tmp2 + gxzc*tmp3
! print*,idetg
! print*
! print*,g4xx,g4xy,g4xz
! print*,g4yy,g4yz,g4zz
! print*
! print*,a,b,c
! print*,sf(i,j,k)
! print*
! print*
! end if
else
sf(i,j,k) = zero
end if
end do
end do
end do
end if
return
end subroutine EHFinder_Sources
|