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
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
|
/*@@
@file GRHydro_UpdateMask.F90
@date Wed Mar 13 14:18:38 2002
@author
@desc
Alter the update terms if inside the atmosphere or excision region
@enddesc
@@*/
#include "cctk.h"
#include "cctk_Parameters.h"
#include "cctk_Arguments.h"
#include "GRHydro_Macros.h"
#include "SpaceMask.h"
#define velx(i,j,k) vel(i,j,k,1)
#define vely(i,j,k) vel(i,j,k,2)
#define velz(i,j,k) vel(i,j,k,3)
#define velx_p(i,j,k) vel_p(i,j,k,1)
#define vely_p(i,j,k) vel_p(i,j,k,2)
#define velz_p(i,j,k) vel_p(i,j,k,3)
#define velx_p_p(i,j,k) vel_p_p(i,j,k,1)
#define vely_p_p(i,j,k) vel_p_p(i,j,k,2)
#define velz_p_p(i,j,k) vel_p_p(i,j,k,3)
subroutine GRHydroUpdateAtmosphereMask(CCTK_ARGUMENTS)
implicit none
DECLARE_CCTK_ARGUMENTS
DECLARE_CCTK_PARAMETERS
CCTK_INT :: i,j,k
CCTK_REAL :: frac
CCTK_INT :: type_bits, atmosphere
call SpaceMask_GetTypeBits(type_bits, "Hydro_Atmosphere")
call SpaceMask_GetStateBits(atmosphere, "Hydro_Atmosphere", "in_atmosphere")
frac = CCTK_DELTA_TIME
!$OMP PARALLEL DO PRIVATE(j,i)
do k = 1, cctk_lsh(3)
do j = 1, cctk_lsh(2)
do i = 1, cctk_lsh(1)
if ( GRHydro_enable_internal_excision /= 0 .and. (hydro_excision_mask(i,j,k) .ne. 0) .or. &
(SpaceMask_CheckStateBitsF90(space_mask, i, j, k, \
type_bits, atmosphere)) .or. &
(tau(i,j,k) + frac * taurhs(i,j,k) .le. 0.d0) .or. &
(dens(i,j,k) + frac * densrhs(i,j,k) .le. 0.d0) ) then
densrhs(i,j,k) = 0.0d0
srhs(i,j,k,:) = 0.0d0
taurhs(i,j,k) = 0.0d0
atmosphere_mask(i,j,k) = 1
SpaceMask_SetStateBitsF90(space_mask, i, j, k, type_bits, atmosphere)
end if
end do
end do
end do
!$OMP END PARALLEL DO
end subroutine GRHydroUpdateAtmosphereMask
/*@@
@routine GRHydro_SetupMask
@date Thu Jun 20 13:27:28 2002
@author Ian Hawke
@desc
Initialize the mask to be zero.
@enddesc
@calls
@calledby
@history
@endhistory
@@*/
subroutine GRHydro_SetupMask(CCTK_ARGUMENTS)
implicit none
DECLARE_CCTK_ARGUMENTS
DECLARE_CCTK_PARAMETERS
CCTK_INT :: type_bits, not_atmosphere, i, j, k
! Initialize all rhs variables and the mask.
! The former vars need to be initialized since there is
! no rhs computation in CCTK_INITIAL or POSTINITIAL.
densrhs = 0.0d0
taurhs = 0.0d0
srhs = 0.0d0
if (evolve_tracer .ne. 0) then
cons_tracerrhs = 0.0d0
endif
atmosphere_mask = 0
call SpaceMask_GetTypeBits(type_bits, "Hydro_Atmosphere")
call SpaceMask_GetStateBits(not_atmosphere, &
&"Hydro_Atmosphere", "not_in_atmosphere")
do k = 1, cctk_lsh(3)
do j = 1, cctk_lsh(2)
do i = 1, cctk_lsh(1)
SpaceMask_SetStateBitsF90(space_mask, i, j, k, type_bits, not_atmosphere)
end do
end do
end do
call CCTK_INFO("Setting up the atmosphere mask: all points are not_atmosphere")
end subroutine GRHydro_SetupMask
/*@@
@routine GRHydro_InitAtmosMask
@date Thu Jun 20 13:27:28 2002
@author Ian Hawke
@desc
Initialize the mask based on rho_min. This is used only if wk_atmosphere=yes.
@enddesc
@calls
@calledby
@history
@endhistory
@@*/
subroutine GRHydro_InitAtmosMask(CCTK_ARGUMENTS)
implicit none
DECLARE_CCTK_ARGUMENTS
CCTK_INT :: type_bits, atmosphere, i, j, k
call SpaceMask_GetTypeBits(type_bits, "Hydro_Atmosphere")
call SpaceMask_GetStateBits(atmosphere, &
&"Hydro_Atmosphere", "in_atmosphere")
do k = 1, cctk_lsh(3)
do j = 1, cctk_lsh(2)
do i = 1, cctk_lsh(1)
if (rho(i,j,k) .le. GRHydro_rho_min) then
SpaceMask_SetStateBitsF90(space_mask, i, j, k, type_bits, atmosphere)
atmosphere_mask(i,j,k) = 1
end if
end do
end do
end do
call CCTK_INFO("Setting up the atmosphere mask: points with rho<rho_min are set to be atmosphere")
end subroutine GRHydro_InitAtmosMask
/*@@
@routine GRHydro_AtmosphereReset
@date Thu Jun 20 13:30:51 2002
@author Ian Hawke
@desc
After MoL has evolved, if a point is supposed to be reset then do so.
@enddesc
@calls
@calledby
@history
@endhistory
@@*/
subroutine GRHydro_AtmosphereReset(CCTK_ARGUMENTS)
implicit none
DECLARE_CCTK_ARGUMENTS
DECLARE_CCTK_PARAMETERS
CCTK_INT :: i, j, k
CCTK_REAL :: det, psi4pt
CCTK_INT :: type_bits, atmosphere, not_atmosphere
call SpaceMask_GetTypeBits(type_bits, "Hydro_Atmosphere")
call SpaceMask_GetStateBits(atmosphere, "Hydro_Atmosphere",&
"in_atmosphere")
call SpaceMask_GetStateBits(not_atmosphere, "Hydro_Atmosphere",&
"not_in_atmosphere")
do k = 1, cctk_lsh(3)
do j = 1, cctk_lsh(2)
do i = 1, cctk_lsh(1)
if ( (atmosphere_mask(i, j, k) .eq. 1) &
&.or. (SpaceMask_CheckStateBitsF90(space_mask,i, j, k, type_bits,\
atmosphere)) &
&) then
!!$ write(*,*) 'Resetting at ',i,j,k, atmosphere_mask(i, j, k), &
!!$ & (SpaceMask_CheckStateBitsF90(space_mask,i, j, k, type_bits,\
!!$ atmosphere))
rho(i,j,k) = GRHydro_rho_min
velx(i,j,k) = 0.0d0
vely(i,j,k) = 0.0d0
velz(i,j,k) = 0.0d0
det = SPATIAL_DETERMINANT(gxx(i,j,k), gxy(i,j,k), gxz(i,j,k), \
gyy(i,j,k), gyz(i,j,k), gzz(i,j,k))
call prim2conpolytype(GRHydro_polytrope_handle, &
gxx(i,j,k), gxy(i,j,k), gxz(i,j,k), &
gyy(i,j,k), gyz(i,j,k), gzz(i,j,k), det, &
dens(i,j,k), scon(i,j,k,1), scon(i,j,k,2), scon(i,j,k,3), &
tau(i,j,k), rho(i,j,k), velx(i,j,k), vely(i,j,k), &
velz(i,j,k), eps(i,j,k), press(i,j,k), w_lorentz(i,j,k))
if (wk_atmosphere .eq. 0) then
atmosphere_mask(i, j, k) = 0
SpaceMask_SetStateBitsF90(space_mask, i, j, k, type_bits,\
not_atmosphere)
end if
end if
end do
end do
end do
!!$ call GRHydro_Boundaries(CCTK_PASS_FTOF)
end subroutine GRHydro_AtmosphereReset
subroutine GRHydro_InitialAtmosphereReset(CCTK_ARGUMENTS)
implicit none
DECLARE_CCTK_ARGUMENTS
DECLARE_CCTK_PARAMETERS
CCTK_INT :: i, j, k
CCTK_REAL :: det, psi4pt
CCTK_INT :: type_bits, atmosphere, not_atmosphere
CCTK_INT :: eos_handle
#if !USE_EOS_OMNI
#ifdef _EOS_BASE_INC_
#undef _EOS_BASE_INC_
#endif
#include "EOS_Base.inc"
#endif
#if USE_EOS_OMNI
! begin EOS Omni vars
integer :: n = 1
integer :: keytemp = 0
integer :: anyerr = 0
integer :: keyerr(1) = 0
real*8 :: xpress = 0.0d0
real*8 :: xeps = 0.0d0
real*8 :: xtemp = 0.0d0
real*8 :: xye = 0.0d0
! end EOS Omni vars
#endif
eos_handle = GRHydro_polytrope_handle
call SpaceMask_GetTypeBits(type_bits, "Hydro_Atmosphere")
call SpaceMask_GetStateBits(atmosphere, "Hydro_Atmosphere",&
"in_atmosphere")
call SpaceMask_GetStateBits(not_atmosphere, "Hydro_Atmosphere",&
"not_in_atmosphere")
do k = 1, cctk_lsh(3)
do j = 1, cctk_lsh(2)
do i = 1, cctk_lsh(1)
if (rho(i,j,k) .le. GRHydro_rho_min) then
rho(i,j,k) = GRHydro_rho_min
velx(i,j,k) = 0.0d0
vely(i,j,k) = 0.0d0
velz(i,j,k) = 0.0d0
#if USE_EOS_OMNI
call EOS_Omni_press(eos_handle,keytemp,GRHydro_eos_rf_prec,n,&
rho(i,j,k),eps(i,j,k),xtemp,xye,press(i,j,k),keyerr,anyerr)
call EOS_Omni_EpsFromPress(eos_handle,keytemp,GRHydro_eos_rf_prec,n,&
rho(i,j,k),xeps,xtemp,xye,press(i,j,k),eps(i,j,k),keyerr,anyerr)
#else
press(i,j,k) = EOS_Pressure(eos_handle, GRHydro_rho_min, eps(i,j,k))
eps(i,j,k) = EOS_SpecificIntEnergy(eos_handle, GRHydro_rho_min, press(i,j,k))
#endif
det = SPATIAL_DETERMINANT(gxx(i,j,k), gxy(i,j,k), gxz(i,j,k), \
gyy(i,j,k), gyz(i,j,k), gzz(i,j,k))
call prim2conpolytype(eos_handle, &
gxx(i,j,k), gxy(i,j,k), gxz(i,j,k), &
gyy(i,j,k), gyz(i,j,k), gzz(i,j,k), det, &
dens(i,j,k), scon(i,j,k,1), scon(i,j,k,2), scon(i,j,k,3), &
tau(i,j,k), rho(i,j,k), velx(i,j,k), vely(i,j,k), &
velz(i,j,k), eps(i,j,k), press(i,j,k), w_lorentz(i,j,k))
end if
if (timelevels .gt. 1) then
if (rho_p(i,j,k) .le. GRHydro_rho_min) then
rho_p(i,j,k) = GRHydro_rho_min
velx_p(i,j,k) = 0.0d0
vely_p(i,j,k) = 0.0d0
velz_p(i,j,k) = 0.0d0
#if USE_EOS_OMNI
call EOS_Omni_press(eos_handle,keytemp,GRHydro_eos_rf_prec,n,&
rho_p(i,j,k),eps_p(i,j,k),xtemp,xye,press_p(i,j,k),keyerr,anyerr)
call EOS_Omni_EpsFromPress(eos_handle,keytemp,GRHydro_eos_rf_prec,n,&
rho_p(i,j,k),xeps,xtemp,xye,press_p(i,j,k),eps_p(i,j,k),keyerr,anyerr)
#else
press_p(i,j,k) = EOS_Pressure(eos_handle, GRHydro_rho_min, eps_p(i,j,k))
eps_p(i,j,k) = EOS_SpecificIntEnergy(eos_handle, GRHydro_rho_min, press_p(i,j,k))
#endif
det = SPATIAL_DETERMINANT(gxx_p(i,j,k), gxy_p(i,j,k), gxz_p(i,j,k), \
gyy_p(i,j,k), gyz_p(i,j,k), gzz_p(i,j,k))
call prim2conpolytype(eos_handle, &
gxx_p(i,j,k), gxy_p(i,j,k), gxz_p(i,j,k), &
gyy_p(i,j,k), gyz_p(i,j,k), gzz_p(i,j,k), det, &
dens_p(i,j,k), scon_p(i,j,k,1), scon_p(i,j,k,2), scon_p(i,j,k,3), &
tau_p(i,j,k), rho_p(i,j,k), velx_p(i,j,k), vely_p(i,j,k), &
velz_p(i,j,k), eps_p(i,j,k), press_p(i,j,k), w_lorentz_p(i,j,k))
endif
end if
if (timelevels .gt. 2) then
if (rho_p_p(i,j,k) .le. GRHydro_rho_min) then
rho_p_p(i,j,k) = GRHydro_rho_min
velx_p_p(i,j,k) = 0.0d0
vely_p_p(i,j,k) = 0.0d0
velz_p_p(i,j,k) = 0.0d0
#if USE_EOS_OMNI
call EOS_Omni_press(eos_handle,keytemp,GRHydro_eos_rf_prec,n,&
rho_p_p(i,j,k),eps_p_p(i,j,k),xtemp,xye,press_p_p(i,j,k),keyerr,anyerr)
call EOS_Omni_EpsFromPress(eos_handle,keytemp,GRHydro_eos_rf_prec,n,&
rho_p_p(i,j,k),xeps,xtemp,xye,press_p_p(i,j,k),eps_p_p(i,j,k),keyerr,anyerr)
#else
press_p_p(i,j,k) = EOS_Pressure(eos_handle, GRHydro_rho_min, eps_p_p(i,j,k))
eps_p_p(i,j,k) = EOS_SpecificIntEnergy(eos_handle, GRHydro_rho_min, press_p_p(i,j,k))
#endif
det = SPATIAL_DETERMINANT(gxx_p_p(i,j,k), gxy_p_p(i,j,k), gxz_p_p(i,j,k), \
gyy_p_p(i,j,k), gyz_p_p(i,j,k), gzz_p_p(i,j,k))
call prim2conpolytype(eos_handle, &
gxx_p_p(i,j,k), gxy_p_p(i,j,k), gxz_p_p(i,j,k), &
gyy_p_p(i,j,k), gyz_p_p(i,j,k), gzz_p_p(i,j,k), det, &
dens_p_p(i,j,k), scon_p_p(i,j,k,1), scon_p_p(i,j,k,2), scon_p_p(i,j,k,3), &
tau_p_p(i,j,k), rho_p_p(i,j,k), velx_p_p(i,j,k), vely_p_p(i,j,k), &
velz_p_p(i,j,k), eps_p_p(i,j,k), press_p_p(i,j,k), w_lorentz_p_p(i,j,k))
endif
endif
end do
end do
end do
write(*,*) " GRHydro_InitialAtmosphereReset"
!!$ call GRHydro_Boundaries(CCTK_PASS_FTOF)
end subroutine GRHydro_InitialAtmosphereReset
|