-
Notifications
You must be signed in to change notification settings - Fork 10
/
Copy pathbcs.f90
384 lines (248 loc) · 10.6 KB
/
bcs.f90
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
371
372
373
374
375
376
377
378
379
380
381
382
383
384
module bcs_module
use datatypes_module
use grid_module
use variables_module
use params_module
use eos_module
use user_bc_module
implicit none
private
public :: fillBCs
contains
subroutine fillBCs(U, g)
type(gridvar_t), intent(inout) :: U, g
integer :: i, ir, ip, ic
real (kind=dp_t) :: p, e, p_above, e_above, p_below, e_below
real (kind=dp_t) :: econst, dens, vel
! sanity check
if ( (U%grid%xlboundary == "periodic" .and. &
U%grid%xrboundary /= "periodic") .or. &
(U%grid%xrboundary == "periodic" .and. &
U%grid%xlboundary /= "periodic") ) then
print *, "ERROR: both boundaries must be periodic"
stop
endif
!-------------------------------------------------------------------------
! lower boundary (-x)
!-------------------------------------------------------------------------
select case (U%grid%xlboundary)
case ("reflect")
! ir is the index that corresponds to the reflected ghostcell
ir = U%grid%lo
do i = U%grid%lo-1, U%grid%lo-U%grid%ng, -1
! even-reflect density and energy
U%data(i,iudens) = U%data(ir,iudens)
U%data(i,iuener) = U%data(ir,iuener)
! odd-reflect the velocities.
U%data(i,iumomx) = -U%data(ir,iumomx)
ir = ir + 1
enddo
case ("hse")
if (gravity_monopole == 1) then
print *, "HSE boundary conditions not supported for monopole gravity"
stop
endif
ir = U%grid%lo
ic = U%grid%lo
do i = U%grid%lo-1, U%grid%lo-U%grid%ng, -1
if (hse_bc_const == "density") then
! constant density in the ghost cells
! zero gradient to rho
U%data(i,iudens) = U%data(ic,iudens)
! velocity depends on hse_vel_type
select case (hse_vel_type)
case ("outflow")
vel = U%data(ic,iumomx)/U%data(ic,iudens)
case ("reflect")
vel = -U%data(ir,iumomx)/U%data(ir,iudens)
ir = ir + 1
case default
print *, "ERROR: invalid hse_vel_type"
stop
end select
U%data(i,iumomx) = U%data(i,iudens)*vel
! p via HSE
e_above = (U%data(i+1,iuener) - &
0.5_dp_t*U%data(i+1,iumomx)**2/U%data(i+1,iudens))/ &
U%data(i+1,iudens)
call eos(eos_input_e, p_above, e_above, U%data(i+1,iudens))
p = p_above - 0.5_dp_t*(U%data(i,iudens) + &
U%data(i+1,iudens))*grav*U%grid%dx
call eos(eos_input_p, p, e, U%data(i,iudens))
! compute E
U%data(i,iuener) = U%data(i,iudens)*e + &
0.5_dp_t*U%data(i,iumomx)**2/U%data(i,iudens)
else if (hse_bc_const == "temperature") then
! constant temperature (or internal energy) in the ghost
! cells
econst = (U%data(ic,iuener) - &
0.5_dp_t*U%data(ic,iumomx)**2/U%data(ic,iudens))/ &
U%data(ic,iudens)
! velocity depends on hse_vel_type
select case (hse_vel_type)
case ("outflow")
vel = U%data(ic,iumomx)/U%data(ic,iudens)
case ("reflect")
vel = -U%data(ir,iumomx)/U%data(ir,iudens)
ir = ir + 1
case default
print *, "ERROR: invalid hse_vel_type"
stop
end select
! get the pressure above (we already know that it has the
! internal energy econst)
call eos(eos_input_e, p_above, econst, U%data(i+1, iudens))
! p via HSE. Here we explicitly make use of the
! gamma-law nature of the EOS. Otherwise, we would need
! to iterate to simultaneously satisfy HSE + the EOS.
p = (p_above - 0.5_dp_t*U%data(i+1,iudens)*grav*U%grid%dx) / &
(1.0_dp_t + 0.5_dp_t*grav*U%grid%dx/ &
(econst*(gamma - 1.0_dp_t)))
! now find the density that corresponds to this p, e
call eos(eos_input_pe, p, econst, dens)
U%data(i,iudens) = dens
U%data(i,iumomx) = dens*vel
! compute E
U%data(i,iuener) = dens*econst + &
0.5_dp_t*U%data(i,iumomx)**2/dens
else
print *, "invalid hse_bc_const type"
stop
endif
enddo
case ("outflow")
do i = U%grid%lo-1, U%grid%lo-U%grid%ng, -1
! give all quantities a zero-gradient
U%data(i,:) = U%data(i+1,:)
enddo
case ("user")
call user_bc_xm(U)
case ("periodic")
ip = U%grid%hi
do i = U%grid%lo-1, U%grid%lo-U%grid%ng, -1
U%data(i,:) = U%data(ip,:)
ip = ip - 1
enddo
case default
print *, "ERROR: xlboundary not valid"
stop
end select
!-------------------------------------------------------------------------
! upper boundary (+x)
!-------------------------------------------------------------------------
select case (U%grid%xrboundary)
case ("reflect")
! ir is the index that corresponds to the reflected ghostcell
ir = U%grid%hi
do i = U%grid%hi+1, U%grid%hi+U%grid%ng
! even reflect density and energy
U%data(i,iudens) = U%data(ir,iudens)
U%data(i,iuener) = U%data(ir,iuener)
! reflect the velocities.
U%data(i,iumomx) = -U%data(ir,iumomx)
ir = ir - 1
enddo
case ("hse")
if (gravity_monopole == 1) then
print *, "HSE boundary conditions not supported for monopole gravity"
stop
endif
ir = U%grid%hi
ic = U%grid%hi
do i = U%grid%hi+1, U%grid%hi+U%grid%ng
if (hse_bc_const == "density") then
! constant density in the ghost cells
! zero gradient to rho
U%data(i,iudens) = U%data(ic,iudens)
! velocity depends on hse_vel_type
select case (hse_vel_type)
case ("outflow")
vel = U%data(ic,iumomx)/U%data(ic,iudens)
case ("reflect")
vel = -U%data(ir,iumomx)/U%data(ir,iudens)
ir = ir - 1
case default
print *, "ERROR: invalid hse_vel_type"
stop
end select
U%data(i,iumomx) = U%data(i,iudens)*vel
! p via HSE
e_below = (U%data(i-1,iuener) - &
0.5_dp_t*U%data(i-1,iumomx)**2/U%data(i-1,iudens))/ &
U%data(i-1,iudens)
call eos(eos_input_e, p_below, e_below, U%data(i-1,iudens))
p = p_below + 0.5_dp_t*(U%data(i,iudens) + &
U%data(i-1,iudens))*grav*U%grid%dx
call eos(eos_input_p, p, e, U%data(i,iudens))
! compute E
U%data(i,iuener) = U%data(i,iudens)*e + &
0.5_dp_t*U%data(i,iumomx)**2/U%data(i,iudens)
else if (hse_bc_const == "temperature") then
! constant temperature (or internal energy) in the ghost
! cells
econst = (U%data(ic,iuener) - &
0.5_dp_t*U%data(ic,iumomx)**2/U%data(ic,iudens))/ &
U%data(ic,iudens)
! velocity depends on hse_vel_type
select case(hse_vel_type)
case ("outflow")
vel = U%data(ic,iumomx)/U%data(ic,iudens)
case ("reflect")
vel = -U%data(ir,iumomx)/U%data(ir,iudens)
ir = ir - 1
case default
print *, "ERROR: invalid hse_vel_type"
stop
end select
! get the pressure below (we already know that it has the
! internal energy econst)
call eos(eos_input_e, p_below, econst, U%data(i-1, iudens))
! p via HSE. Here we explicitly make use of the
! gamma-law nature of the EOS. Otherwise, we would need
! to iterate to simultaneously satisfy HSE + the EOS.
p = (p_below + 0.5_dp_t*U%data(i-1,iudens)*grav*U%grid%dx) / &
(1.0_dp_t - 0.5_dp_t*grav*U%grid%dx/ &
(econst*(gamma - 1.0_dp_t)))
! now find the density that corresponds to this p, e
call eos(eos_input_pe, p, econst, dens)
U%data(i,iudens) = dens
U%data(i,iumomx) = dens*vel
! compute E
U%data(i,iuener) = dens*econst + &
0.5_dp_t*U%data(i,iumomx)**2/dens
else
print *, "invalid hse_bc_const type"
stop
endif
enddo
case ("outflow")
do i = U%grid%hi+1, U%grid%hi+U%grid%ng
! give all quantities a zero-gradient
U%data(i,:) = U%data(i-1,:)
enddo
case ("diode")
do i = U%grid%hi+1, U%grid%hi+U%grid%ng
! give all quantities a zero-gradient
U%data(i,:) = U%data(i-1,:)
! store internal energy
e = (U%data(i,iuener) - &
0.5_dp_t*U%data(i,iumomx)**2/U%data(i,iudens))/U%data(i,iudens)
U%data(i,iumomx) = max(0.0_dp_t, U%data(i,iumomx))
! recompute energy
U%data(i,iuener) = U%data(i,iudens)*e + &
0.5_dp_t*U%data(i,iumomx)**2/U%data(i,iudens)
enddo
case ("user")
call user_bc_xp(U)
case ("periodic")
ip = U%grid%lo
do i = U%grid%hi+1, U%grid%hi+U%grid%ng
U%data(i,:) = U%data(ip,:)
ip = ip + 1
enddo
case default
print *, "ERROR: xrboundary not valid"
stop
end select
end subroutine fillBCs
end module bcs_module