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
|
/*@@
@file GRHydro_UpdateMaskM.F90
@date Sep 2, 2010
@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"
!!$ We don't need to adapt GRHydroUpdateAtmosphereMask, GRHydro_SetupMask, or
!!$ since we need to evolve Bvec in the atmosphere
!!$ In GRHydro_AtmosphereResetM, we just need to switch the P2C calls to MHD
/*@@
@routine GRHydro_AtmosphereResetM
@date Sep 2, 2010
@author Joshua Faber, Scott Noble, Bruno Mundim, 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_AtmosphereResetM(CCTK_ARGUMENTS)
implicit none
DECLARE_CCTK_ARGUMENTS
DECLARE_CCTK_PARAMETERS
CCTK_INT :: i, j, k
CCTK_REAL :: det, psi4pt
! begin EOS Omni vars
integer :: n,keytemp,anyerr,keyerr(1)
real*8 :: xpress,xeps,xtemp,xye
n=1;keytemp=0;anyerr=0;keyerr(1)=0
xpress=0.0d0;xye=0.0d0;xeps=0.0d0;xtemp=0.0d0
! end EOS Omni vars
do k = 1, cctk_lsh(3)
do j = 1, cctk_lsh(2)
do i = 1, cctk_lsh(1)
if (atmosphere_mask(i, j, k) .ne. 0) then
rho(i,j,k) = GRHydro_rho_min
vel(i,j,k,1) = 0.0d0
vel(i,j,k,2) = 0.0d0
vel(i,j,k,3) = 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))
if(evolve_temper.ne.0) then
! ! set the temperature to be relatively low
temperature(i,j,k) = grhydro_hot_atmo_temp
y_e(i,j,k) = grhydro_hot_atmo_Y_e
keytemp = 1
call EOS_Omni_press(GRHydro_eos_handle,keytemp,GRHydro_eos_rf_prec,n,&
rho(i,j,k),eps(i,j,k),temperature(i,j,k),y_e(i,j,k),&
press(i,j,k),keyerr,anyerr)
call prim2conM_hot(GRHydro_eos_handle, GRHydro_reflevel,&
i,j,k,x(i,j,k),y(i,j,k),z(i,j,k),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),Bcons(i,j,k,1),Bcons(i,j,k,2),Bcons(i,j,k,3),&
rho(i,j,k),vel(i,j,k,1),vel(i,j,k,2),vel(i,j,k,3),&
eps(i,j,k),press(i,j,k),Bvec(i,j,k,1), &
Bvec(i,j,k,2), Bvec(i,j,k,3), w_lorentz(i,j,k),&
temperature(i,j,k),y_e(i,j,k))
y_e_con(i,j,k) = dens(i,j,k) * y_e(i,j,k)
else
call prim2conpolytypeM(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), Bcons(i,j,k,1),Bcons(i,j,k,2),Bcons(i,j,k,3),&
rho(i,j,k), vel(i,j,k,1), vel(i,j,k,2), &
vel(i,j,k,3), eps(i,j,k), press(i,j,k), &
Bvec(i,j,k,1),Bvec(i,j,k,2),Bvec(i,j,k,3),w_lorentz(i,j,k))
if (wk_atmosphere .eq. 0) then
atmosphere_mask(i, j, k) = 0
end if
end if
endif
end do
end do
end do
!!$ call GRHydro_BoundariesM(CCTK_PASS_FTOF)
end subroutine GRHydro_AtmosphereResetM
subroutine GRHydro_InitialAtmosphereResetM(CCTK_ARGUMENTS)
implicit none
DECLARE_CCTK_ARGUMENTS
DECLARE_CCTK_PARAMETERS
CCTK_INT :: i, j, k
CCTK_REAL :: det, psi4pt
CCTK_REAL :: rho_min
CCTK_INT :: eos_handle
! begin EOS Omni vars
integer :: n,keytemp,anyerr,keyerr(1)
real*8 :: xpress,xeps,xtemp,xye
n=1;keytemp=0;anyerr=0;keyerr(1)=0
xpress=0.0d0;xye=0.0d0;xeps=0.0d0;xtemp=0.0d0
! end EOS Omni vars
eos_handle = GRHydro_polytrope_handle
rho_min = GRHydro_rho_min
if (initial_atmosphere_factor .gt. 0) then
rho_min = rho_min * initial_atmosphere_factor
endif
do k = 1, cctk_lsh(3)
do j = 1, cctk_lsh(2)
do i = 1, cctk_lsh(1)
if (rho(i,j,k) .le. rho_min) then
rho(i,j,k) = rho_min
vel(i,j,k,1) = 0.0d0
vel(i,j,k,2) = 0.0d0
vel(i,j,k,3) = 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))
if(evolve_temper.ne.0) then
! ! set the temperature to be relatively low
temperature(i,j,k) = grhydro_hot_atmo_temp
y_e(i,j,k) = grhydro_hot_atmo_Y_e
keytemp = 1
call EOS_Omni_press(GRHydro_eos_handle,keytemp,GRHydro_eos_rf_prec,n,&
rho(i,j,k),eps(i,j,k),temperature(i,j,k),y_e(i,j,k),&
press(i,j,k),keyerr,anyerr)
call prim2conM_hot(GRHydro_eos_handle, GRHydro_reflevel,&
i,j,k,x(i,j,k),y(i,j,k),z(i,j,k),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),Bcons(i,j,k,1),Bcons(i,j,k,2),Bcons(i,j,k,3),&
rho(i,j,k),vel(i,j,k,1),vel(i,j,k,2),vel(i,j,k,3),&
eps(i,j,k),press(i,j,k),Bvec(i,j,k,1), &
Bvec(i,j,k,2), Bvec(i,j,k,3), w_lorentz(i,j,k),&
temperature(i,j,k),y_e(i,j,k))
y_e_con(i,j,k) = dens(i,j,k) * y_e(i,j,k)
else
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)
call prim2conpolytypeM(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), Bcons(i,j,k,1),Bcons(i,j,k,2),Bcons(i,j,k,3),&
rho(i,j,k), vel(i,j,k,1), vel(i,j,k,2), &
vel(i,j,k,3), eps(i,j,k), press(i,j,k), &
Bvec(i,j,k,1),Bvec(i,j,k,2),Bvec(i,j,k,3),w_lorentz(i,j,k))
end if
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
vel_p(i,j,k,1) = 0.0d0
vel_p(i,j,k,2) = 0.0d0
vel_p(i,j,k,3) = 0.0d0
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))
if(evolve_temper.ne.0) then
! ! set the temperature to be relatively low
temperature_p(i,j,k) = grhydro_hot_atmo_temp
y_e_p(i,j,k) = grhydro_hot_atmo_Y_e
keytemp = 1
call EOS_Omni_press(GRHydro_eos_handle,keytemp,GRHydro_eos_rf_prec,n,&
rho_p(i,j,k),eps_p(i,j,k),temperature_p(i,j,k),y_e_p(i,j,k),&
press_p(i,j,k),keyerr,anyerr)
call prim2conM_hot(GRHydro_eos_handle, GRHydro_reflevel,&
i,j,k,x(i,j,k),y(i,j,k),z(i,j,k),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),Bcons_p(i,j,k,1),Bcons_p(i,j,k,2),Bcons_p(i,j,k,3),&
rho_p(i,j,k),vel_p(i,j,k,1),vel_p(i,j,k,2),vel_p(i,j,k,3),&
eps_p(i,j,k),press_p(i,j,k),Bvec_p(i,j,k,1), &
Bvec_p(i,j,k,2), Bvec_p(i,j,k,3), w_lorentz_p(i,j,k),&
temperature_p(i,j,k),y_e_p(i,j,k))
y_e_con_p(i,j,k) = dens_p(i,j,k) * y_e_p(i,j,k)
else
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)
call prim2conpolytypeM(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), Bcons_p(i,j,k,1),Bcons_p(i,j,k,2),Bcons_p(i,j,k,3),&
rho_p(i,j,k), vel_p(i,j,k,1), vel_p(i,j,k,2), &
vel_p(i,j,k,3), eps_p(i,j,k), press_p(i,j,k), &
Bvec_p(i,j,k,1),Bvec_p(i,j,k,2),Bvec_p(i,j,k,3),w_lorentz_p(i,j,k))
end if
end if
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
vel_p_p(i,j,k,1) = 0.0d0
vel_p_p(i,j,k,2) = 0.0d0
vel_p_p(i,j,k,3) = 0.0d0
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))
if(evolve_temper.ne.0) then
! ! set the temperature to be relatively low
temperature_p_p(i,j,k) = grhydro_hot_atmo_temp
y_e_p_p(i,j,k) = grhydro_hot_atmo_Y_e
keytemp = 1
call EOS_Omni_press(GRHydro_eos_handle,keytemp,GRHydro_eos_rf_prec,n,&
rho_p_p(i,j,k),eps_p_p(i,j,k),temperature_p_p(i,j,k),y_e_p_p(i,j,k),&
press_p_p(i,j,k),keyerr,anyerr)
call prim2conM_hot(GRHydro_eos_handle, GRHydro_reflevel,&
i,j,k,x(i,j,k),y(i,j,k),z(i,j,k),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),Bcons_p_p(i,j,k,1),Bcons_p_p(i,j,k,2),Bcons_p_p(i,j,k,3),&
rho_p_p(i,j,k),vel_p_p(i,j,k,1),vel_p_p(i,j,k,2),vel_p_p(i,j,k,3),&
eps_p_p(i,j,k),press_p_p(i,j,k),Bvec_p_p(i,j,k,1), &
Bvec_p_p(i,j,k,2), Bvec_p_p(i,j,k,3), w_lorentz_p_p(i,j,k),&
temperature_p_p(i,j,k),y_e_p_p(i,j,k))
y_e_con_p_p(i,j,k) = dens_p_p(i,j,k) * y_e_p_p(i,j,k)
else
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)
call prim2conpolytypeM(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), Bcons_p_p(i,j,k,1),Bcons_p_p(i,j,k,2),Bcons_p_p(i,j,k,3),&
rho_p_p(i,j,k), vel_p_p(i,j,k,1), vel_p_p(i,j,k,2), &
vel_p_p(i,j,k,3), eps_p_p(i,j,k), press_p_p(i,j,k), &
Bvec_p_p(i,j,k,1),Bvec_p_p(i,j,k,2),Bvec_p_p(i,j,k,3),w_lorentz_p_p(i,j,k))
end if
end if
end if
end do
end do
end do
write(*,*) " GRHydro_InitialAtmosphereReset"
!!$ call GRHydro_BoundariesM(CCTK_PASS_FTOF)
end subroutine GRHydro_InitialAtmosphereResetM
|