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
|
/*@@
@file GRHydro_ShockTubeM.F90
@date Sep 23, 2010
@author Joshua Faber, Scott Noble, Bruno Mundim, Ian Hawke
@desc
Initial data of the shock tube type - MHD version.
@enddesc
@@*/
#include "cctk.h"
#include "cctk_Parameters.h"
#include "cctk_Arguments.h"
#include "cctk_Functions.h"
#include "GRHydro_Macros.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 sx(i,j,k) scon(i,j,k,1)
#define sy(i,j,k) scon(i,j,k,2)
#define sz(i,j,k) scon(i,j,k,3)
#define Bvecx(i,j,k) Bvec(i,j,k,1)
#define Bvecy(i,j,k) Bvec(i,j,k,2)
#define Bvecz(i,j,k) Bvec(i,j,k,3)
/*@@
@routine GRHydro_shocktubeM
@date Sat Jan 26 02:53:49 2002
@author Joshua Faber, Scott Noble, Bruno Mundim, Ian Hawke
@desc
Initial data for shock tubes, parallel to
a coordinate axis. Either Sods problem or the standard shock tube.
@enddesc
@calls
@calledby
@history
Expansion and alteration of the test code from GRAstro_Hydro,
written by Mark Miller.
@endhistory
@@*/
subroutine GRHydro_shocktubeM(CCTK_ARGUMENTS)
implicit none
DECLARE_CCTK_ARGUMENTS
DECLARE_CCTK_PARAMETERS
DECLARE_CCTK_FUNCTIONS
CCTK_INT :: i, j, k, nx, ny, nz
CCTK_REAL :: direction, det
CCTK_REAL :: rhol, rhor, velxl, velxr, velyl, velyr, &
velzl, velzr, epsl, epsr
CCTK_REAL :: bvcxl,bvcyl,bvczl,bvcxr,bvcyr,bvczr
bvcxl = Bx_init
bvcyl = By_init
bvczl = Bz_init
bvcxr = Bx_init
bvcyr = By_init
bvczr = Bz_init
nx = cctk_lsh(1)
ny = cctk_lsh(2)
nz = cctk_lsh(3)
do i=1,nx
do j=1,ny
do k=1,nz
if (CCTK_EQUALS(shocktube_type,"diagshock")) then
direction = x(i,j,k) - shock_xpos + &
y(i,j,k) - shock_ypos + z(i,j,k) - shock_zpos
else if (CCTK_EQUALS(shocktube_type,"xshock")) then
direction = x(i,j,k) - shock_xpos
else if (CCTK_EQUALS(shocktube_type,"yshock")) then
direction = y(i,j,k) - shock_ypos
else if (CCTK_EQUALS(shocktube_type,"zshock")) then
direction = z(i,j,k) - shock_zpos
else if (CCTK_EQUALS(shocktube_type,"sphere")) then
direction = sqrt((x(i,j,k)-shock_xpos)**2+&
(y(i,j,k)-shock_ypos)**2+&
(z(i,j,k)-shock_zpos)**2)-shock_radius
end if
if (CCTK_EQUALS(shock_case,"Simple")) then
rhol = 10.d0
rhor = 1.d0
velxl = 0.d0
velxr = 0.d0
velyl = 0.d0
velyr = 0.d0
velzl = 0.d0
velzr = 0.d0
epsl = 2.d0
epsr = 1.d-6
else if (CCTK_EQUALS(shock_case,"Sod")) then
rhol = 1.d0
rhor = 0.125d0
velxl = 0.d0
velxr = 0.d0
velyl = 0.d0
velyr = 0.d0
velzl = 0.d0
velzr = 0.d0
epsl = 1.5d0
epsr = 1.2d0
!!$This line only for polytrope, k=1
!!$ epsr = 0.375d0
else if (CCTK_EQUALS(shock_case,"Blast")) then
rhol = 1.d0
rhor = 1.d0
velxl = 0.d0
velxr = 0.d0
velyl = 0.d0
velyr = 0.d0
velzl = 0.d0
velzr = 0.d0
epsl = 1500.d0
epsr = 1.5d-2
else if (CCTK_EQUALS(shock_case,"Balsara1")) then
rhol = 1.d0
rhor = 0.125d0
velxl = 0.d0
velxr = 0.d0
velyl = 0.d0
velyr = 0.d0
velzl = 0.d0
velzr = 0.d0
bvcxl=0.5d0
bvcxr=0.5d0
bvcyl=1.d0
bvcyr=-1.d0
bvczl=0.d0
bvczr=0.d0
epsl = 1.0/rhol
epsr = 0.1/rhor
else if (CCTK_EQUALS(shock_case,"Balsara2")) then
rhol = 1.d0
rhor = 1.d0
velxl = 0.d0
velxr = 0.d0
velyl = 0.d0
velyr = 0.d0
velzl = 0.d0
velzr = 0.d0
bvcxl=5.0d0
bvcxr=5.0d0
bvcyl=6.d0
bvcyr=0.7d0
bvczl=6.d0
bvczr=0.7d0
epsl = 1.5d0*30.0d0/rhol
epsr = 1.5d0*1.0d0/rhor
else if (CCTK_EQUALS(shock_case,"Balsara3")) then
rhol = 1.d0
rhor = 1.d0
velxl = 0.d0
velxr = 0.d0
velyl = 0.d0
velyr = 0.d0
velzl = 0.d0
velzr = 0.d0
bvcxl=10.0d0
bvcxr=10.0d0
bvcyl=7.d0
bvcyr=0.7d0
bvczl=7.d0
bvczr=0.7d0
epsl = 1.5d0*1000.0d0/rhol
epsr = 1.5d0*0.1d0/rhor
else if (CCTK_EQUALS(shock_case,"Balsara4")) then
rhol = 1.d0
rhor = 1.d0
velxl = 0.999d0
velxr = -0.999d0
velyl = 0.d0
velyr = 0.d0
velzl = 0.d0
velzr = 0.d0
bvcxl=10.0d0
bvcxr=10.0d0
bvcyl=7.d0
bvcyr=-7.d0
bvczl=7.d0
bvczr=-7.d0
epsl = 1.5d0*0.1d0/rhol
epsr = 1.5d0*0.1d0/rhor
else if (CCTK_EQUALS(shock_case,"Balsara5")) then
rhol = 1.08d0
rhor = 1.d0
velxl = 0.4d0
velxr = -0.45d0
velyl = 0.3d0
velyr = -0.2d0
velzl = 0.2d0
velzr = 0.2d0
bvcxl=2.0d0
bvcxr=2.0d0
bvcyl=0.3d0
bvcyr=-0.7d0
bvczl=0.3d0
bvczr=0.5d0
epsl = 1.5d0*0.95d0/rhol
epsr = 1.5d0*1.0d0/rhor
else
call CCTK_WARN(0,"Shock case not recognized")
end if
if ( ((change_shock_direction==0).and.(direction .lt. 0.0d0)).or.&
((change_shock_direction==1).and.(direction .gt. 0.0d0)) ) then
rho(i,j,k) = rhol
velx(i,j,k) = velxl
vely(i,j,k) = velyl
velz(i,j,k) = velzl
eps(i,j,k) = epsl
Bvecx(i,j,k)=bvcxl
Bvecy(i,j,k)=bvcyl
Bvecz(i,j,k)=bvczl
else
rho(i,j,k) = rhor
velx(i,j,k) = velxr
vely(i,j,k) = velyr
velz(i,j,k) = velzr
eps(i,j,k) = epsr
Bvecx(i,j,k)=bvcxr
Bvecy(i,j,k)=bvcyr
Bvecz(i,j,k)=bvczr
end if
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 (CCTK_EQUALS(GRHydro_eos_type,"Polytype")) then
call Prim2ConPolyM(GRHydro_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),sx(i,j,k),sy(i,j,k),sz(i,j,k),&
tau(i,j,k),Bvecx(i,j,k),Bvecy(i,j,k),Bvecz(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))
else
call Prim2ConGenM(GRHydro_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),sx(i,j,k),sy(i,j,k),sz(i,j,k),&
tau(i,j,k),Bvecx(i,j,k),Bvecy(i,j,k),Bvecz(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
enddo
enddo
enddo
densrhs = 0.d0
srhs = 0.d0
taurhs = 0.d0
Bvecrhs = 0.d0
return
end subroutine GRHydro_shocktubeM
|