MPI-AMRVAC 3.2
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
Loading...
Searching...
No Matches
mod_finite_volume.t
Go to the documentation of this file.
1!> Module with finite volume methods for fluxes
3#include "amrvac.h"
4 implicit none
5 private
6
7 public :: finite_volume
8 public :: hancock
9 public :: reconstruct_lr
10
11contains
12
13 !> The non-conservative Hancock predictor for TVDLF
14 !>
15 !> on entry:
16 !> input available on ixI^L=ixG^L asks for output on ixO^L=ixG^L^LSUBnghostcells
17 !> one entry: (predictor): wCT -- w_n wnew -- w_n qdt=dt/2
18 !> on exit : (predictor): wCT -- w_n wnew -- w_n+1/2
19 subroutine hancock(qdt,dtfactor,ixI^L,ixO^L,idims^LIM,qtC,sCT,qt,snew,dxs,x)
20 use mod_physics
22 use mod_source, only: addsource2
23 use mod_comm_lib, only: mpistop
24
25 integer, intent(in) :: ixi^l, ixo^l, idims^lim
26 double precision, intent(in) :: qdt, dtfactor,qtc, qt, dxs(ndim), x(ixi^s,1:ndim)
27 type(state) :: sct, snew
28
29 double precision, dimension(ixI^S,1:nw) :: wprim, wlc, wrc
30 ! left and right constructed status in primitive form, needed for better performance
31 double precision, dimension(ixI^S,1:nw) :: wlp, wrp
32 double precision, dimension(ixO^S) :: inv_volume
33 double precision :: flc(ixi^s, nwflux), frc(ixi^s, nwflux)
34 double precision :: dxinv(1:ndim)
35 integer :: idims, iw, ix^l, hxo^l
36 logical :: active=.false.
37
38 associate(wct=>sct%w,wnew=>snew%w)
39 ! Expand limits in each idims direction in which fluxes are added
40 ix^l=ixo^l;
41 do idims= idims^lim
42 ix^l=ix^l^laddkr(idims,^d);
43 end do
44 if (ixi^l^ltix^l|.or.|.or.) &
45 call mpistop("Error in Hancock: Nonconforming input limits")
46
47 wrp=0.d0
48 wlp=0.d0
49 wprim=wct
50 call phys_to_primitive(ixi^l,ixi^l,wprim,x)
51
52 dxinv=-qdt/dxs
53 if(.not.slab_uniform) inv_volume = 1.d0/block%dvolume(ixo^s)
54 do idims= idims^lim
55 b0i=idims
56 ! Calculate w_j+g_j/2 and w_j-g_j/2
57 ! First copy all variables, then upwind wLC and wRC.
58 ! wLC is to the left of ixO, wRC is to the right of wCT.
59 hxo^l=ixo^l-kr(idims,^d);
60
61 wrp(hxo^s,1:nwflux)=wprim(ixo^s,1:nwflux)
62 wlp(ixo^s,1:nwflux)=wprim(ixo^s,1:nwflux)
63
64 ! apply limited reconstruction for left and right status at cell interfaces
65 call reconstruct_lr(ixi^l,ixo^l,hxo^l,idims,wprim,wlc,wrc,wlp,wrp,x,dxs(idims))
66
67 ! Calculate the fLC and fRC fluxes
68 call phys_get_flux(wrc,wrp,x,ixi^l,hxo^l,idims,frc)
69 call phys_get_flux(wlc,wlp,x,ixi^l,ixo^l,idims,flc)
70
71 ! Advect w(iw)
72 if (slab_uniform) then
73 if(local_timestep) then
74 do iw=1,nwflux
75 wnew(ixo^s,iw)=wnew(ixo^s,iw)-block%dt(ixo^s)*dtfactor/dxs(idims)* &
76 (flc(ixo^s, iw)-frc(hxo^s, iw))
77 end do
78 else
79 do iw=1,nwflux
80 wnew(ixo^s,iw)=wnew(ixo^s,iw)+dxinv(idims)* &
81 (flc(ixo^s, iw)-frc(hxo^s, iw))
82 end do
83 endif
84 else
85 if(local_timestep) then
86 do iw=1,nwflux
87 wnew(ixo^s,iw)=wnew(ixo^s,iw) - block%dt(ixo^s)*dtfactor * inv_volume &
88 *(block%surfaceC(ixo^s,idims)*flc(ixo^s, iw) &
89 -block%surfaceC(hxo^s,idims)*frc(hxo^s, iw))
90 end do
91 else
92 do iw=1,nwflux
93 wnew(ixo^s,iw)=wnew(ixo^s,iw) - qdt * inv_volume &
94 *(block%surfaceC(ixo^s,idims)*flc(ixo^s, iw) &
95 -block%surfaceC(hxo^s,idims)*frc(hxo^s, iw))
96 end do
97 end if
98 end if
99 end do ! next idims
100 b0i=0
101
102 if (.not.slab.and.idimsmin==1) call phys_add_source_geom(qdt,dtfactor,ixi^l,ixo^l,wct,wprim,wnew,x)
103
104 call addsource2(qdt*dble(idimsmax-idimsmin+1)/dble(ndim), &
105 dtfactor*dble(idimsmax-idimsmin+1)/dble(ndim),&
106 ixi^l,ixo^l,1,nw,qtc,wct,wprim,qt,wnew,x,.false.,active)
107
108 ! check and optionally correct unphysical values
109 if(fix_small_values) then
110 call phys_handle_small_values(.false.,wnew,x,ixi^l,ixo^l,'exit hancock finite_volume')
111 endif
112 end associate
113 end subroutine hancock
114
115 !> finite volume method
116 subroutine finite_volume(method,qdt,dtfactor,ixI^L,ixO^L,idims^LIM, &
117 qtC,sCT,qt,snew,fC,fE,dxs,x)
118 use mod_physics
119 use mod_variables
121 use mod_tvd, only:tvdlimit2
122 use mod_source, only: addsource2
124 use mod_comm_lib, only: mpistop
125
126 integer, intent(in) :: method
127 double precision, intent(in) :: qdt, dtfactor, qtc, qt, dxs(ndim)
128 integer, intent(in) :: ixi^l, ixo^l, idims^lim
129 double precision, dimension(ixI^S,1:ndim), intent(in) :: x
130 double precision, dimension(ixI^S,1:nwflux,1:ndim) :: fc
131 double precision, dimension(ixI^S,sdim:3) :: fe
132 type(state) :: sct, snew
133
134 ! primitive w at cell center
135 double precision, dimension(ixI^S,1:nw) :: wprim
136 ! left and right constructed status in conservative form
137 double precision, dimension(ixI^S,1:nw) :: wlc, wrc
138 ! left and right constructed status in primitive form, needed for better performance
139 double precision, dimension(ixI^S,1:nw) :: wlp, wrp
140 double precision, dimension(ixI^S,1:nwflux) :: flc, frc
141 double precision, dimension(ixI^S,1:number_species) :: cmaxc
142 double precision, dimension(ixI^S,1:number_species) :: cminc
143 double precision, dimension(ixI^S) :: hspeed
144 double precision, dimension(ixO^S) :: inv_volume
145 double precision, dimension(1:ndim) :: dxinv
146 integer :: idims, iw, ix^d, hx^d, ix^l, hxo^l, ixc^l, ixcr^l, kxc^l, kxr^l, ii
147 logical :: active=.false.
148 type(ct_velocity) :: vcts
149
150 associate(wct=>sct%w, wnew=>snew%w)
151
152 ! The flux calculation contracts by one in the idims direction it is applied.
153 ! The limiter contracts the same directions by one more, so expand ixO by 2.
154 ix^l=ixo^l;
155 do idims= idims^lim
156 ix^l=ix^l^ladd2*kr(idims,^d);
157 end do
158 if (ixi^l^ltix^l|.or.|.or.) &
159 call mpistop("Error in fv : Nonconforming input limits")
160
161 wprim=wct
162 call phys_to_primitive(ixi^l,ixi^l,wprim,x)
163 do idims= idims^lim
164 ! use interface value of w0 at idims
165 b0i=idims
166
167 kxcmin^d=iximin^d; kxcmax^d=iximax^d-kr(idims,^d);
168 kxr^l=kxc^l+kr(idims,^d);
169 ! wRp and wLp are defined at the same locations, and will correspond to
170 ! the left and right reconstructed values at a cell face. Their indexing
171 ! is similar to cell-centered values, but in direction idims they are
172 ! shifted half a cell towards the 'lower' direction.
173 do iw=1,nwflux
174 {do ix^db=iximin^db,iximax^db\}
175 ! fill all cells for averaging to fix small values
176 wrp(ix^d,iw)=wprim(ix^d,iw)
177 wlp(ix^d,iw)=wprim(ix^d,iw)
178 {end do\}
179 wrp(kxc^s,iw)=wprim(kxr^s,iw)
180 end do
181
182 hxo^l=ixo^l-kr(idims,^d);
183 if(stagger_grid) then
184 ! ct needs 1 or 2 (hll) ghost cells in the transverse dimensions
185 ixcmax^d=ixomax^d+transverse_ghost_cells-transverse_ghost_cells*kr(idims,^d);
186 ixcmin^d=hxomin^d-transverse_ghost_cells+transverse_ghost_cells*kr(idims,^d);
187 else
188 ! ixC is centered index in the idims direction from ixOmin-1/2 to ixOmax+1/2
189 ixcmax^d=ixomax^d; ixcmin^d=hxomin^d;
190 end if
191
192
193 ! Determine stencil size
194 {ixcrmin^d = max(ixcmin^d - phys_wider_stencil,ixglo^d)\}
195 {ixcrmax^d = min(ixcmax^d + phys_wider_stencil,ixghi^d)\}
196
197 ! apply limited reconstruction for left and right status at cell interfaces
198 call reconstruct_lr(ixi^l,ixcr^l,ixcr^l,idims,wprim,wlc,wrc,wlp,wrp,x,dxs(idims))
199
200 ! special modification of left and right status before flux evaluation
201 call phys_modify_wlr(ixi^l,ixcr^l,qt,wlc,wrc,wlp,wrp,sct,idims)
202
203 ! evaluate physical fluxes according to reconstructed status
204 call phys_get_flux(wlc,wlp,x,ixi^l,ixc^l,idims,flc)
205 call phys_get_flux(wrc,wrp,x,ixi^l,ixc^l,idims,frc)
206 if(h_correction) then
207 call phys_get_h_speed(wprim,x,ixi^l,ixo^l,idims,hspeed)
208 end if
209 ! estimating bounds for the minimum and maximum signal velocities
210 if(method==fs_tvdlf.or.method==fs_tvdmu) then
211 call phys_get_cbounds(wlc,wrc,wlp,wrp,x,ixi^l,ixc^l,idims,hspeed,cmaxc)
212 ! index of var velocity appears in the induction eq.
213 if(stagger_grid) call phys_get_ct_velocity(vcts,wlp,wrp,ixi^l,ixc^l,idims,cmaxc(ixi^s,index_v_mag))
214 else
215 call phys_get_cbounds(wlc,wrc,wlp,wrp,x,ixi^l,ixc^l,idims,hspeed,cmaxc,cminc)
216 if(stagger_grid) call phys_get_ct_velocity(vcts,wlp,wrp,ixi^l,ixc^l,idims,cmaxc(ixi^s,index_v_mag),cminc(ixi^s,index_v_mag))
217 end if
218
219 ! use approximate Riemann solver to get flux at interfaces
220 select case(method)
221 case(fs_hll)
222 do ii=1,number_species
223 call get_riemann_flux_hll(start_indices(ii),stop_indices(ii))
224 end do
225 case(fs_hllc,fs_hllcd)
226 do ii=1,number_species
227 call get_riemann_flux_hllc(start_indices(ii),stop_indices(ii))
228 end do
229 case(fs_hlld)
230 do ii=1,number_species
231 if(ii==index_v_mag) then
232 call get_riemann_flux_hlld(start_indices(ii),stop_indices(ii))
233 else
234 call get_riemann_flux_hll(start_indices(ii),stop_indices(ii))
235 endif
236 end do
237 case(fs_tvdlf)
238 do ii=1,number_species
239 call get_riemann_flux_tvdlf(start_indices(ii),stop_indices(ii))
240 end do
241 case(fs_tvdmu)
243 case default
244 call mpistop('unkown Riemann flux in finite volume')
245 end select
246
247 end do ! Next idims
248 b0i=0
249 if(stagger_grid) call phys_update_faces(ixi^l,ixo^l,qt,qdt,wprim,fc,fe,sct,snew,vcts)
250 if(slab_uniform) then
251 if(local_timestep) then
252 dxinv(1:ndim)=-dtfactor/dxs(1:ndim)
253 do idims= idims^lim
254 hx^d=kr(idims,^d)\
255 hxomin^d=ixomin^d-hx^d\
256 do iw=iwstart,nwflux
257 {do ix^db=hxomin^db,ixomax^db\}
258 fc(ix^d,iw,idims)=block%dt(ix^d)*dxinv(idims)*fc(ix^d,iw,idims)
259 {end do\}
260 {do ix^db=ixomin^db,ixomax^db\}
261 wnew(ix^d,iw)=wnew(ix^d,iw)+fc(ix^d,iw,idims)-fc(ix^d-hx^d,iw,idims)
262 {end do\}
263 end do
264 ! For the MUSCL scheme apply the characteristic based limiter
265 if(method==fs_tvdmu) &
266 call tvdlimit2(method,qdt,ixi^l,ixc^l,ixo^l,idims,wlc,wrc,wnew,x,fc,dxs)
267 end do
268 else
269 dxinv(1:ndim)=-qdt/dxs(1:ndim)
270 do idims= idims^lim
271 hx^d=kr(idims,^d)\
272 hxomin^d=ixomin^d-hx^d\
273 do iw=iwstart,nwflux
274 {do ix^db=hxomin^db,ixomax^db\}
275 fc(ix^d,iw,idims)=dxinv(idims)*fc(ix^d,iw,idims)
276 {end do\}
277 {do ix^db=ixomin^db,ixomax^db\}
278 wnew(ix^d,iw)=wnew(ix^d,iw)+fc(ix^d,iw,idims)-fc(ix^d-hx^d,iw,idims)
279 {end do\}
280 end do
281 ! For the MUSCL scheme apply the characteristic based limiter
282 if(method==fs_tvdmu) &
283 call tvdlimit2(method,qdt,ixi^l,ixc^l,ixo^l,idims,wlc,wrc,wnew,x,fc,dxs)
284 end do
285 end if
286 else
287 inv_volume(ixo^s) = 1.d0/block%dvolume(ixo^s)
288 if(local_timestep) then
289 do idims= idims^lim
290 hx^d=kr(idims,^d)\
291 hxomin^d=ixomin^d-hx^d\
292 do iw=iwstart,nwflux
293 {do ix^db=hxomin^db,ixomax^db\}
294 fc(ix^d,iw,idims)=-block%dt(ix^d)*dtfactor*fc(ix^d,iw,idims)*block%surfaceC(ix^d,idims)
295 {end do\}
296 {do ix^db=ixomin^db,ixomax^db\}
297 wnew(ix^d,iw)=wnew(ix^d,iw)+(fc(ix^d,iw,idims)-fc(ix^d-hx^d,iw,idims))*inv_volume(ix^d)
298 {end do\}
299 end do
300 ! For the MUSCL scheme apply the characteristic based limiter
301 if (method==fs_tvdmu) &
302 call tvdlimit2(method,qdt,ixi^l,ixc^l,ixo^l,idims,wlc,wrc,wnew,x,fc,dxs)
303 end do
304 else
305 do idims= idims^lim
306 hx^d=kr(idims,^d)\
307 hxomin^d=ixomin^d-hx^d\
308 do iw=iwstart,nwflux
309 {do ix^db=hxomin^db,ixomax^db\}
310 fc(ix^d,iw,idims)=-qdt*fc(ix^d,iw,idims)*block%surfaceC(ix^d,idims)
311 {end do\}
312 {do ix^db=ixomin^db,ixomax^db\}
313 wnew(ix^d,iw)=wnew(ix^d,iw)+(fc(ix^d,iw,idims)-fc(ix^d-hx^d,iw,idims))*inv_volume(ix^d)
314 {end do\}
315 end do
316 end do
317 end if
318 end if
319
320 if (.not.slab.and.idimsmin==1) &
321 call phys_add_source_geom(qdt,dtfactor,ixi^l,ixo^l,wct,wprim,wnew,x)
322
323 if(stagger_grid) call phys_face_to_center(ixo^l,snew)
324
325 ! check and optionally correct unphysical values
326 if(fix_small_values) then
327 call phys_handle_small_values(.false.,wnew,x,ixi^l,ixo^l,'multi-D finite_volume')
328 if(crash) then
329 ! replace erroneous values with values at previous step
330 wnew=wct
331 if(stagger_grid) snew%ws=sct%ws
332 end if
333 end if
334
335 call addsource2(qdt*dble(idimsmax-idimsmin+1)/dble(ndim),&
336 dtfactor*dble(idimsmax-idimsmin+1)/dble(ndim),&
337 ixi^l,ixo^l,1,nw,qtc,wct,wprim,qt,wnew,x,.false.,active)
338
339 end associate
340 contains
341
343 do iw=iwstart,nwflux
344 fc(ixc^s,iw,idims)=half*(flc(ixc^s,iw)+frc(ixc^s,iw))
345 end do
346 end subroutine get_riemann_flux_tvdmu
347
348 subroutine get_riemann_flux_tvdlf(iws,iwe)
349 integer, intent(in) :: iws,iwe
350
351 integer :: ix^D,jx^D
352 double precision :: fac(ixC^S),phi
353
354 fac(ixc^s) = -0.5d0*tvdlfeps*cmaxc(ixc^s,ii)
355 do iw=iws,iwe
356 if(flux_energy_only .and. iw /= iw_e) then
357 fc(ixc^s,iw,idims)=zero
358 else
359 fc(ixc^s,iw,idims)=0.5d0*(flc(ixc^s, iw)+frc(ixc^s, iw))
360 ! Add TVDLF dissipation to the flux
361 if(flux_type(idims, iw) /= flux_no_dissipation) then
362 if(flux_adaptive_diffusion) then
363 {do ix^db=ixcmin^db,ixcmax^db\}
364 jx^d=ix^d+kr(idims,^d)\
365 !> adaptive diffusion from Rempel et al. 2009, see also Rempel et al. 2014
366 !> the previous version is adopt
367 !if(((wRC(ix^D,iw)-wLC(ix^D,iw))*(sCT%w(jx^D,iw)-sCT%w(ix^D,iw))) .gt. 1.e-18) then
368 ! phi=min((wRC(ix^D,iw)-wLC(ix^D,iw))**2/((sCT%w(jx^D,iw)-sCT%w(ix^D,iw))**2+1.e-18),one)
369 !else
370 ! phi=1.d0
371 !end if
372 phi = flux_adaptive_diffusion_min
373 if(((wrc(ix^d,iw)-wlc(ix^d,iw))*(sct%w(jx^d,iw)-sct%w(ix^d,iw))) .gt. 1.d-18) then
374 phi = max(flux_adaptive_diffusion_min, &
375 min(flux_adaptive_diffusion_scale * &
376 (wrc(ix^d,iw)-wlc(ix^d,iw))**2 / &
377 ((sct%w(jx^d,iw)-sct%w(ix^d,iw))**2 + 1.d-18), one))
378 else
379 phi = one
380 end if
381 fc(ix^d,iw,idims)=fc(ix^d,iw,idims)+fac(ix^d)*(wrc(ix^d,iw)-wlc(ix^d,iw))*phi
382 {end do\}
383 else
384 {do ix^db=ixcmin^db,ixcmax^db\}
385 fc(ix^d,iw,idims)=fc(ix^d,iw,idims)+fac(ix^d)*(wrc(ix^d,iw)-wlc(ix^d,iw))
386 {end do\}
387 end if
388 end if
389 end if
390 end do
391 end subroutine get_riemann_flux_tvdlf
392
393 subroutine get_riemann_flux_hll(iws,iwe)
394 integer, intent(in) :: iws,iwe
395 integer :: ix^D
396 double precision :: phi
397
398 if(flux_adaptive_diffusion) then
399 do iw=iws,iwe
400 if(flux_type(idims, iw) == flux_tvdlf) then
401 if(stagger_grid) then
402 ! CT MHD set zero normal B flux
403 fc(ixc^s,iw,idims)=0.d0
404 else
405 fc(ixc^s,iw,idims)=-tvdlfeps*half*max(cmaxc(ixc^s,ii),dabs(cminc(ixc^s,ii)))*&
406 (wrc(ixc^s,iw)-wlc(ixc^s,iw))
407 end if
408 else
409 {do ix^db=ixcmin^db,ixcmax^db\}
410 if(cminc(ix^d,ii) >= zero) then
411 fc(ix^d,iw,idims)=flc(ix^d,iw)
412 else if(cmaxc(ix^d,ii) <= zero) then
413 fc(ix^d,iw,idims)=frc(ix^d,iw)
414 else
415 !> reduced diffusion is from Wang et al. 2024
416 phi=max(abs(cmaxc(ix^d,ii)),abs(cminc(ix^d,ii)))/(cmaxc(ix^d,ii)-cminc(ix^d,ii))
417 fc(ix^d,iw,idims)=(cmaxc(ix^d,ii)*flc(ix^d, iw)-cminc(ix^d,ii)*frc(ix^d,iw)&
418 +phi*cminc(ix^d,ii)*cmaxc(ix^d,ii)*(wrc(ix^d,iw)-wlc(ix^d,iw)))&
419 /(cmaxc(ix^d,ii)-cminc(ix^d,ii))
420 end if
421 {end do\}
422 end if
423 end do
424 else
425 do iw=iws,iwe
426 if(flux_type(idims, iw) == flux_tvdlf) then
427 if(stagger_grid) then
428 ! CT MHD set zero normal B flux
429 fc(ixc^s,iw,idims)=0.d0
430 else
431 fc(ixc^s,iw,idims)=-tvdlfeps*half*max(cmaxc(ixc^s,ii),dabs(cminc(ixc^s,ii)))*&
432 (wrc(ixc^s,iw)-wlc(ixc^s,iw))
433 end if
434 else
435 {do ix^db=ixcmin^db,ixcmax^db\}
436 if(cminc(ix^d,ii) >= zero) then
437 fc(ix^d,iw,idims)=flc(ix^d,iw)
438 else if(cmaxc(ix^d,ii) <= zero) then
439 fc(ix^d,iw,idims)=frc(ix^d,iw)
440 else
441 fc(ix^d,iw,idims)=(cmaxc(ix^d,ii)*flc(ix^d, iw)-cminc(ix^d,ii)*frc(ix^d,iw)&
442 +cminc(ix^d,ii)*cmaxc(ix^d,ii)*(wrc(ix^d,iw)-wlc(ix^d,iw)))&
443 /(cmaxc(ix^d,ii)-cminc(ix^d,ii))
444 end if
445 {end do\}
446 end if
447 end do
448 end if
449 end subroutine get_riemann_flux_hll
450
451 subroutine get_riemann_flux_hllc(iws,iwe)
452 integer, intent(in) :: iws, iwe
453 double precision, dimension(ixI^S,1:nwflux) :: whll, Fhll, fCD
454 double precision, dimension(ixI^S) :: lambdaCD
455
456 integer, dimension(ixI^S) :: patchf
457 integer :: rho_, p_, e_, mom(1:ndir)
458
459 rho_ = iw_rho
460 if (allocated(iw_mom)) mom(:) = iw_mom(:)
461 e_ = iw_e
462
463 if(associated(phys_hllc_init_species)) then
464 call phys_hllc_init_species(ii, rho_, mom(:), e_)
465 endif
466
467 p_ = e_
468
469 patchf(ixc^s) = 1
470 where(cminc(ixc^s,1) >= zero)
471 patchf(ixc^s) = -2
472 elsewhere(cmaxc(ixc^s,1) <= zero)
473 patchf(ixc^s) = 2
474 endwhere
475 ! Use more diffusive scheme, is actually TVDLF and selected by patchf=4
476 if(method==fs_hllcd) &
477 call phys_diffuse_hllcd(ixi^l,ixc^l,idims,wlc,wrc,flc,frc,patchf)
478
479 !---- calculate speed lambda at CD ----!
480 if(any(patchf(ixc^s)==1)) &
481 call phys_get_lcd(wlc,wrc,flc,frc,cminc(ixi^s,ii),cmaxc(ixi^s,ii),idims,ixi^l,ixc^l, &
482 whll,fhll,lambdacd,patchf)
483
484 ! now patchf may be -1 or 1 due to phys_get_lCD
485 if(any(abs(patchf(ixc^s))== 1))then
486 !======== flux at intermediate state ========!
487 call phys_get_wcd(wlc,wrc,whll,frc,flc,fhll,patchf,lambdacd,&
488 cminc(ixi^s,ii),cmaxc(ixi^s,ii),ixi^l,ixc^l,idims,fcd)
489 endif ! Calculate the CD flux
490
491 do iw=iws,iwe
492 if (flux_type(idims, iw) == flux_tvdlf) then
493 flc(ixc^s,iw)=-tvdlfeps*half*max(cmaxc(ixc^s,ii),abs(cminc(ixc^s,ii))) * &
494 (wrc(ixc^s,iw) - wlc(ixc^s,iw))
495 else
496 where(patchf(ixc^s)==-2)
497 flc(ixc^s,iw)=flc(ixc^s,iw)
498 elsewhere(abs(patchf(ixc^s))==1)
499 flc(ixc^s,iw)=fcd(ixc^s,iw)
500 elsewhere(patchf(ixc^s)==2)
501 flc(ixc^s,iw)=frc(ixc^s,iw)
502 elsewhere(patchf(ixc^s)==3)
503 ! fallback option, reducing to HLL flux
504 flc(ixc^s,iw)=fhll(ixc^s,iw)
505 elsewhere(patchf(ixc^s)==4)
506 ! fallback option, reducing to TVDLF flux
507 flc(ixc^s,iw) = half*((flc(ixc^s,iw)+frc(ixc^s,iw)) &
508 -tvdlfeps * max(cmaxc(ixc^s,ii), dabs(cminc(ixc^s,ii))) * &
509 (wrc(ixc^s,iw)-wlc(ixc^s,iw)))
510 endwhere
511 end if
512
513 fc(ixc^s,iw,idims)=flc(ixc^s,iw)
514
515 end do ! Next iw
516 end subroutine get_riemann_flux_hllc
517
518 !> HLLD Riemann flux from Miyoshi 2005 JCP, 208, 315 and Guo 2016 JCP, 327, 543
519 subroutine get_riemann_flux_hlld(iws,iwe)
520 integer, intent(in) :: iws, iwe
521 double precision, dimension(ixI^S,1:nwflux) :: w1R,w1L,w2R,w2L
522 double precision, dimension(ixI^S) :: sm,s1R,s1L,suR,suL,Bx
523 double precision, dimension(ixI^S) :: pts,ptR,ptL,signBx,r1L,r1R,tmp
524 ! magnetic field from the right and the left reconstruction
525 double precision, dimension(ixI^S,ndir) :: BR, BL
526 integer :: ip1,ip2,ip3,idir,ix^D,^C&b^C_,^C&m^C_
527 integer :: rho_, p_, e_, mom(1:ndir), mag(1:ndir)
528
529 associate(sr=>cmaxc,sl=>cminc)
530
531 rho_=iw_rho
532 ^c&mom(^c)=iw_mom(^c)\
533 m^c_=mom(^c);
534 ^c&mag(^c)=iw_mag(^c)\
535 b^c_=mag(^c);
536 e_ = iw_e
537 p_ = e_
538
539 ip1=idims
540 ip3=3
541 if(b0field) then
542 br(ixc^s,:)=wrc(ixc^s,mag(:))+block%B0(ixc^s,:,ip1)
543 bl(ixc^s,:)=wlc(ixc^s,mag(:))+block%B0(ixc^s,:,ip1)
544 else
545 br(ixc^s,:)=wrc(ixc^s,mag(:))
546 bl(ixc^s,:)=wlc(ixc^s,mag(:))
547 end if
548 if(stagger_grid) then
549 bx(ixc^s)=block%ws(ixc^s,ip1)
550 else
551 ! HLL estimation of normal magnetic field at cell interfaces
552 ! Li, Shenghai, 2005 JCP, 203, 344, equation (33)
553 bx(ixc^s)=(sr(ixc^s,ii)*br(ixc^s,ip1)-sl(ixc^s,ii)*bl(ixc^s,ip1))/(sr(ixc^s,ii)-sl(ixc^s,ii))
554 end if
555 {!DEC$ VECTOR ALWAYS
556 do ix^db=ixcmin^db,ixcmax^db\}
557 ptr(ix^d)=wrp(ix^d,p_)+0.5d0*(^c&br(ix^d,^c)**2+)
558 ptl(ix^d)=wlp(ix^d,p_)+0.5d0*(^c&bl(ix^d,^c)**2+)
559 if(iw_equi_rho>0) then
560 sur(ix^d)=(sr(ix^d,ii)-wrp(ix^d,mom(ip1)))*(wrc(ix^d,rho_)+block%equi_vars(ix^d,iw_equi_rho,ip1))
561 sul(ix^d)=(sl(ix^d,ii)-wlp(ix^d,mom(ip1)))*(wlc(ix^d,rho_)+block%equi_vars(ix^d,iw_equi_rho,ip1))
562 else
563 sur(ix^d)=(sr(ix^d,ii)-wrp(ix^d,mom(ip1)))*wrc(ix^d,rho_)
564 sul(ix^d)=(sl(ix^d,ii)-wlp(ix^d,mom(ip1)))*wlc(ix^d,rho_)
565 end if
566 ! Miyoshi equation (38) and Guo euqation (20)
567 sm(ix^d)=(sur(ix^d)*wrp(ix^d,mom(ip1))-sul(ix^d)*wlp(ix^d,mom(ip1))-&
568 ptr(ix^d)+ptl(ix^d))/(sur(ix^d)-sul(ix^d))
569 ! Miyoshi equation (39) and Guo euqation (28)
570 w1r(ix^d,mom(ip1))=sm(ix^d)
571 w1l(ix^d,mom(ip1))=sm(ix^d)
572 w2r(ix^d,mom(ip1))=sm(ix^d)
573 w2l(ix^d,mom(ip1))=sm(ix^d)
574 ! Guo equation (22)
575 w1r(ix^d,mag(ip1))=bx(ix^d)
576 w1l(ix^d,mag(ip1))=bx(ix^d)
577 if(b0field) then
578 ptr(ix^d)=wrp(ix^d,p_)+0.5d0*(^c&wrc(ix^d,b^c_)**2+)
579 ptl(ix^d)=wlp(ix^d,p_)+0.5d0*(^c&wlc(ix^d,b^c_)**2+)
580 end if
581 ! Miyoshi equation (43) and Guo equation (27)
582 w1r(ix^d,rho_)=sur(ix^d)/(sr(ix^d,ii)-sm(ix^d))
583 w1l(ix^d,rho_)=sul(ix^d)/(sl(ix^d,ii)-sm(ix^d))
584 ip2=mod(ip1+1,ndir)
585 if(ip2==0) ip2=ndir
586 r1r(ix^d)=sur(ix^d)*(sr(ix^d,ii)-sm(ix^d))-bx(ix^d)**2
587 if(r1r(ix^d)/=0.d0) r1r(ix^d)=1.d0/r1r(ix^d)
588 r1l(ix^d)=sul(ix^d)*(sl(ix^d,ii)-sm(ix^d))-bx(ix^d)**2
589 if(r1l(ix^d)/=0.d0) r1l(ix^d)=1.d0/r1l(ix^d)
590 ! Miyoshi equation (44)
591 w1r(ix^d,mom(ip2))=wrp(ix^d,mom(ip2))-bx(ix^d)*br(ix^d,ip2)*&
592 (sm(ix^d)-wrp(ix^d,mom(ip1)))*r1r(ix^d)
593 w1l(ix^d,mom(ip2))=wlp(ix^d,mom(ip2))-bx(ix^d)*bl(ix^d,ip2)*&
594 (sm(ix^d)-wlp(ix^d,mom(ip1)))*r1l(ix^d)
595 ! partial solution for later usage
596 w1r(ix^d,mag(ip2))=(sur(ix^d)*(sr(ix^d,ii)-wrp(ix^d,mom(ip1)))-bx(ix^d)**2)*r1r(ix^d)
597 w1l(ix^d,mag(ip2))=(sul(ix^d)*(sl(ix^d,ii)-wlp(ix^d,mom(ip1)))-bx(ix^d)**2)*r1l(ix^d)
598 {^ifthreec
599 ip3=mod(ip1+2,ndir)
600 if(ip3==0) ip3=ndir
601 ! Miyoshi equation (46)
602 w1r(ix^d,mom(ip3))=wrp(ix^d,mom(ip3))-bx(ix^d)*br(ix^d,ip3)*&
603 (sm(ix^d)-wrp(ix^d,mom(ip1)))*r1r(ix^d)
604 w1l(ix^d,mom(ip3))=wlp(ix^d,mom(ip3))-bx(ix^d)*bl(ix^d,ip3)*&
605 (sm(ix^d)-wlp(ix^d,mom(ip1)))*r1l(ix^d)
606 ! Miyoshi equation (47)
607 w1r(ix^d,mag(ip3))=br(ix^d,ip3)*w1r(ix^d,mag(ip2))
608 w1l(ix^d,mag(ip3))=bl(ix^d,ip3)*w1l(ix^d,mag(ip2))
609 }
610 ! Miyoshi equation (45)
611 w1r(ix^d,mag(ip2))=br(ix^d,ip2)*w1r(ix^d,mag(ip2))
612 w1l(ix^d,mag(ip2))=bl(ix^d,ip2)*w1l(ix^d,mag(ip2))
613 if(b0field) then
614 ! Guo equation (26)
615 ^c&w1r(ix^d,b^c_)=w1r(ix^d,b^c_)-block%B0(ix^d,^c,ip1)\
616 ^c&w1l(ix^d,b^c_)=w1l(ix^d,b^c_)-block%B0(ix^d,^c,ip1)\
617 end if
618 ! equation (48)
619 if(phys_energy) then
620 ! Guo equation (25) equivalent to Miyoshi equation (41)
621 w1r(ix^d,p_)=sur(ix^d)*(sm(ix^d)-wrp(ix^d,mom(ip1)))+ptr(ix^d)
622 w1l(ix^d,p_)=w1r(ix^d,p_)
623 if(b0field) then
624 ! Guo equation (32)
625 w1r(ix^d,p_)=w1r(ix^d,p_)+(^c&block%B0(ix^d,^c,ip1)*(wrc(ix^d,b^c_)-w1r(ix^d,b^c_))+)
626 w1l(ix^d,p_)=w1l(ix^d,p_)+(^c&block%B0(ix^d,^c,ip1)*(wlc(ix^d,b^c_)-w1l(ix^d,b^c_))+)
627 end if
628 ! Miyoshi equation (48) and main part of Guo euqation (31)
629 w1r(ix^d,e_)=((sr(ix^d,ii)-wrp(ix^d,mom(ip1)))*wrc(ix^d,e_)-ptr(ix^d)*wrp(ix^d,mom(ip1))+&
630 w1r(ix^d,p_)*sm(ix^d)+bx(ix^d)*((^c&wrp(ix^d,m^c_)*wrc(ix^d,b^c_)+)-&
631 (^c&w1r(ix^d,m^c_)*w1r(ix^d,b^c_)+)))/(sr(ix^d,ii)-sm(ix^d))
632 w1l(ix^d,e_)=((sl(ix^d,ii)-wlp(ix^d,mom(ip1)))*wlc(ix^d,e_)-ptl(ix^d)*wlp(ix^d,mom(ip1))+&
633 w1l(ix^d,p_)*sm(ix^d)+bx(ix^d)*((^c&wlp(ix^d,m^c_)*wlc(ix^d,b^c_)+)-&
634 (^c&w1l(ix^d,m^c_)*w1l(ix^d,b^c_)+)))/(sl(ix^d,ii)-sm(ix^d))
635 if(b0field) then
636 ! Guo equation (31)
637 w1r(ix^d,e_)=w1r(ix^d,e_)+((^c&w1r(ix^d,b^c_)*block%B0(ix^d,^c,ip1)+)*sm(ix^d)-&
638 (^c&wrc(ix^d,b^c_)*block%B0(ix^d,^c,ip1)+)*wrp(ix^d,mom(ip1)))/(sr(ix^d,ii)-sm(ix^d))
639 w1l(ix^d,e_)=w1l(ix^d,e_)+((^c&w1l(ix^d,b^c_)*block%B0(ix^d,^c,ip1)+)*sm(ix^d)-&
640 (^c&wlc(ix^d,b^c_)*block%B0(ix^d,^c,ip1)+)*wlp(ix^d,mom(ip1)))/(sl(ix^d,ii)-sm(ix^d))
641 end if
642 if(iw_equi_p>0) then
643 w1r(ix^d,e_)=w1r(ix^d,e_)+1d0/(phys_gamma-1)*block%equi_vars(ix^d,iw_equi_p,ip1)*&
644 (sm(ix^d)-wrp(ix^d,mom(ip1)))/(sr(ix^d,ii)-sm(ix^d))
645 w1l(ix^d,e_)=w1l(ix^d,e_)+1d0/(phys_gamma-1)*block%equi_vars(ix^d,iw_equi_p,ip1)*&
646 (sm(ix^d)-wlp(ix^d,mom(ip1)))/(sl(ix^d,ii)-sm(ix^d))
647 end if
648 end if
649
650 ! Miyoshi equation (49) and Guo equation (35)
651 w2r(ix^d,rho_)=w1r(ix^d,rho_)
652 w2l(ix^d,rho_)=w1l(ix^d,rho_)
653 w2r(ix^d,mag(ip1))=w1r(ix^d,mag(ip1))
654 w2l(ix^d,mag(ip1))=w1l(ix^d,mag(ip1))
655 r1r(ix^d)=sqrt(w1r(ix^d,rho_))
656 r1l(ix^d)=sqrt(w1l(ix^d,rho_))
657 tmp(ix^d)=1.d0/(r1r(ix^d)+r1l(ix^d))
658 signbx(ix^d)=sign(1.d0,bx(ix^d))
659 ! Miyoshi equation (51) and Guo equation (33)
660 s1r(ix^d)=sm(ix^d)+abs(bx(ix^d))/r1r(ix^d)
661 s1l(ix^d)=sm(ix^d)-abs(bx(ix^d))/r1l(ix^d)
662 ! Miyoshi equation (59) and Guo equation (41)
663 w2r(ix^d,mom(ip2))=(r1l(ix^d)*w1l(ix^d,mom(ip2))+r1r(ix^d)*w1r(ix^d,mom(ip2))+&
664 (w1r(ix^d,mag(ip2))-w1l(ix^d,mag(ip2)))*signbx(ix^d))*tmp(ix^d)
665 w2l(ix^d,mom(ip2))=w2r(ix^d,mom(ip2))
666 ! Miyoshi equation (61) and Guo equation (43)
667 w2r(ix^d,mag(ip2))=(r1l(ix^d)*w1r(ix^d,mag(ip2))+r1r(ix^d)*w1l(ix^d,mag(ip2))+&
668 r1l(ix^d)*r1r(ix^d)*(w1r(ix^d,mom(ip2))-w1l(ix^d,mom(ip2)))*signbx(ix^d))*tmp(ix^d)
669 w2l(ix^d,mag(ip2))=w2r(ix^d,mag(ip2))
670 {^ifthreec
671 ! Miyoshi equation (60) and Guo equation (42)
672 w2r(ix^d,mom(ip3))=(r1l(ix^d)*w1l(ix^d,mom(ip3))+r1r(ix^d)*w1r(ix^d,mom(ip3))+&
673 (w1r(ix^d,mag(ip3))-w1l(ix^d,mag(ip3)))*signbx(ix^d))*tmp(ix^d)
674 w2l(ix^d,mom(ip3))=w2r(ix^d,mom(ip3))
675 ! Miyoshi equation (62) and Guo equation (44)
676 w2r(ix^d,mag(ip3))=(r1l(ix^d)*w1r(ix^d,mag(ip3))+r1r(ix^d)*w1l(ix^d,mag(ip3))+&
677 r1l(ix^d)*r1r(ix^d)*(w1r(ix^d,mom(ip3))-w1l(ix^d,mom(ip3)))*signbx(ix^d))*tmp(ix^d)
678 w2l(ix^d,mag(ip3))=w2r(ix^d,mag(ip3))
679 }
680 ! Miyoshi equation (63) and Guo equation (45)
681 if(phys_energy) then
682 w2r(ix^d,e_)=w1r(ix^d,e_)+r1r(ix^d)*((^c&w1r(ix^d,m^c_)*w1r(ix^d,b^c_)+)-&
683 (^c&w2r(ix^d,m^c_)*w2r(ix^d,b^c_)+))*signbx(ix^d)
684 w2l(ix^d,e_)=w1l(ix^d,e_)-r1l(ix^d)*((^c&w1l(ix^d,m^c_)*w1l(ix^d,b^c_)+)-&
685 (^c&w2l(ix^d,m^c_)*w2l(ix^d,b^c_)+))*signbx(ix^d)
686 end if
687
688 ! convert velocity to momentum
689 ^c&w1r(ix^d,m^c_)=w1r(ix^d,m^c_)*w1r(ix^d,rho_)\
690 ^c&w1l(ix^d,m^c_)=w1l(ix^d,m^c_)*w1l(ix^d,rho_)\
691 ^c&w2r(ix^d,m^c_)=w2r(ix^d,m^c_)*w2r(ix^d,rho_)\
692 ^c&w2l(ix^d,m^c_)=w2l(ix^d,m^c_)*w2l(ix^d,rho_)\
693 if(iw_equi_rho>0) then
694 w1r(ix^d,rho_)=w1r(ix^d,rho_)-block%equi_vars(ix^d,iw_equi_rho,ip1)
695 w1l(ix^d,rho_)=w1l(ix^d,rho_)-block%equi_vars(ix^d,iw_equi_rho,ip1)
696 w2r(ix^d,rho_)=w2r(ix^d,rho_)-block%equi_vars(ix^d,iw_equi_rho,ip1)
697 w2l(ix^d,rho_)=w2l(ix^d,rho_)-block%equi_vars(ix^d,iw_equi_rho,ip1)
698 end if
699 {end do\}
700
701 do iw=iws,iwe
702 if(flux_type(idims, iw)==flux_special) then
703 ! known flux (fLC=fRC) for normal B and psi_ in GLM method
704 {!dir$ ivdep
705 do ix^db=ixcmin^db,ixcmax^db\}
706 fc(ix^d,iw,ip1)=flc(ix^d,iw)
707 {end do\}
708 else if(flux_type(idims, iw)==flux_hll) then
709 ! using hll flux for tracers
710 {!dir$ ivdep
711 do ix^db=ixcmin^db,ixcmax^db\}
712 fc(ix^d,iw,ip1)=(sr(ix^d,ii)*flc(ix^d,iw)-sl(ix^d,ii)*frc(ix^d,iw) &
713 +sr(ix^d,ii)*sl(ix^d,ii)*(wrc(ix^d,iw)-wlc(ix^d,iw)))/(sr(ix^d,ii)-sl(ix^d,ii))
714 {end do\}
715 else
716 ! construct hlld flux
717 ! Miyoshi equation (66) and Guo equation (46)
718 {!dir$ ivdep
719 !DEC$ VECTOR ALWAYS
720 do ix^db=ixcmin^db,ixcmax^db\}
721 if(sl(ix^d,ii)>0.d0) then
722 fc(ix^d,iw,ip1)=flc(ix^d,iw)
723 else if(s1l(ix^d)>=0.d0) then
724 fc(ix^d,iw,ip1)=flc(ix^d,iw)+sl(ix^d,ii)*(w1l(ix^d,iw)-wlc(ix^d,iw))
725 else if(sm(ix^d)>=0.d0) then
726 fc(ix^d,iw,ip1)=flc(ix^d,iw)+sl(ix^d,ii)*(w1l(ix^d,iw)-wlc(ix^d,iw))+&
727 s1l(ix^d)*(w2l(ix^d,iw)-w1l(ix^d,iw))
728 else if(s1r(ix^d)>=0.d0) then
729 fc(ix^d,iw,ip1)=frc(ix^d,iw)+sr(ix^d,ii)*(w1r(ix^d,iw)-wrc(ix^d,iw))+&
730 s1r(ix^d)*(w2r(ix^d,iw)-w1r(ix^d,iw))
731 else if(sr(ix^d,ii)>=0.d0) then
732 fc(ix^d,iw,ip1)=frc(ix^d,iw)+sr(ix^d,ii)*(w1r(ix^d,iw)-wrc(ix^d,iw))
733 else if(sr(ix^d,ii)<0.d0) then
734 fc(ix^d,iw,ip1)=frc(ix^d,iw)
735 end if
736 {end do\}
737 end if
738 end do
739
740 end associate
741 end subroutine get_riemann_flux_hlld
742
743 !> HLLD Riemann flux from Miyoshi 2005 JCP, 208, 315 and Guo 2016 JCP, 327, 543
744 !> https://arxiv.org/pdf/2108.04991.pdf
745 subroutine get_riemann_flux_hlld_mag2(iws,iwe)
746 implicit none
747 integer, intent(in) :: iws, iwe
748
749 double precision, dimension(ixI^S,1:nwflux) :: w1R,w1L,f1R,f1L,f2R,f2L
750 double precision, dimension(ixI^S,1:nwflux) :: w2R,w2L
751 double precision, dimension(ixI^S) :: sm,s1R,s1L,suR,suL,Bx
752 double precision, dimension(ixI^S) :: pts,ptR,ptL,signBx,r1L,r1R,tmp
753 ! velocity from the right and the left reconstruction
754 double precision, dimension(ixI^S,ndir) :: vRC, vLC
755 ! magnetic field from the right and the left reconstruction
756 double precision, dimension(ixI^S,ndir) :: BR, BL
757 integer :: ip1,ip2,ip3,idir,ix^D
758 double precision :: phiPres, thetaSM, du, dv, dw
759 integer :: ixV^L, ixVb^L, ixVc^L, ixVd^L, ixVe^L, ixVf^L
760 integer :: rho_, p_, e_, mom(1:ndir), mag(1:ndir)
761 double precision, parameter :: aParam = 4d0
762
763 rho_ = iw_rho
764 mom(:) = iw_mom(:)
765 mag(:) = iw_mag(:)
766 p_ = iw_e
767 e_ = iw_e
768
769 associate(sr=>cmaxc,sl=>cminc)
770
771 f1r=0.d0
772 f1l=0.d0
773 ip1=idims
774 ip3=3
775
776 vrc(ixc^s,:)=wrp(ixc^s,mom(:))
777 vlc(ixc^s,:)=wlp(ixc^s,mom(:))
778
779 ! reuse s1L s1R
780 call get_hlld2_modif_c(wlp,x,ixi^l,ixo^l,s1l)
781 call get_hlld2_modif_c(wrp,x,ixi^l,ixo^l,s1r)
782 !phiPres = min(1, maxval(max(s1L(ixO^S),s1R(ixO^S))/cmaxC(ixO^S,1)))
783 phipres = min(1d0, maxval(max(s1l(ixo^s),s1r(ixo^s)))/maxval(cmaxc(ixo^s,1)))
784 phipres = phipres*(2d0 - phipres)
785
786 !we use here not reconstructed velocity: wprim?
787 ixv^l=ixo^l;
788 !first dim
789 ixvmin1=ixomin1+1
790 ixvmax1=ixomax1+1
791 du = minval(wprim(ixv^s,mom(1))-wprim(ixo^s,mom(1)))
792 if(du>0d0) du=0d0
793 dv = 0d0
794 dw = 0d0
795
796 {^nooned
797 !second dim
798 !i,j-1,k
799 ixv^l=ixo^l;
800 ixvmin2=ixomin2-1
801 ixvmax2=ixomax2-1
802
803 !i,j+1,k
804 ixvb^l=ixo^l;
805 ixvbmin2=ixomin2+1
806 ixvbmax2=ixomax2+1
807
808 !i+1,j,k
809 ixvc^l=ixo^l;
810 ixvcmin1=ixomin1+1
811 ixvcmax1=ixomax1+1
812
813 !i+1,j-1,k
814 ixvd^l=ixo^l;
815 ixvdmin1=ixomin1+1
816 ixvdmax1=ixomax1+1
817 ixvdmin2=ixomin2-1
818 ixvdmax2=ixomax2-1
819
820 !i+1,j+1,k
821 ixve^l=ixo^l;
822 ixvemin1=ixomin1+1
823 ixvemax1=ixomax1+1
824 ixvemin2=ixomin2+1
825 ixvemax2=ixomax2+1
826
827 dv = minval(min(wprim(ixo^s,mom(2))-wprim(ixv^s,mom(2)),&
828 wprim(ixvb^s,mom(2))-wprim(ixo^s,mom(2)),&
829 wprim(ixvc^s,mom(2))-wprim(ixvd^s,mom(2)),&
830 wprim(ixve^s,mom(2))-wprim(ixvc^s,mom(2))&
831 ))
832 if(dv>0d0) dv=0d0}
833
834 {^ifthreed
835 !third dim
836 !i,j,k-1
837 ixv^l=ixo^l;
838 ixvmin3=ixomin3-1
839 ixvmax3=ixomax3-1
840
841 !i,j,k+1
842 ixvb^l=ixo^l;
843 ixvbmin3=ixomin3+1
844 ixvbmax3=ixomax3+1
845
846 !i+1,j,k
847 ixvc^l=ixo^l;
848 ixvcmin1=ixomin1+1
849 ixvcmax1=ixomax1+1
850
851 !i+1,j,k-1
852 ixvd^l=ixo^l;
853 ixvdmin1=ixomin1+1
854 ixvdmax1=ixomax1+1
855 ixvdmin3=ixomin3-1
856 ixvdmax3=ixomax3-1
857
858 !i+1,j,k+1
859 ixve^l=ixo^l;
860 ixvemin1=ixomin1+1
861 ixvemax1=ixomax1+1
862 ixvemin3=ixomin3+1
863 ixvemax3=ixomax3+1
864 dw = minval(min(wprim(ixo^s,mom(3))-wprim(ixv^s,mom(3)),&
865 wprim(ixvb^s,mom(3))-wprim(ixo^s,mom(3)),&
866 wprim(ixvc^s,mom(3))-wprim(ixvd^s,mom(3)),&
867 wprim(ixve^s,mom(3))-wprim(ixvc^s,mom(3))&
868 ))
869 if(dw>0d0) dw=0d0}
870 thetasm = maxval(cmaxc(ixo^s,1))
871
872 thetasm = (min(1d0, (thetasm-du)/(thetasm-min(dv,dw))))**aparam
873
874 if(b0field) then
875 br(ixc^s,:)=wrc(ixc^s,mag(:))+block%B0(ixc^s,:,ip1)
876 bl(ixc^s,:)=wlc(ixc^s,mag(:))+block%B0(ixc^s,:,ip1)
877 else
878 br(ixc^s,:)=wrc(ixc^s,mag(:))
879 bl(ixc^s,:)=wlc(ixc^s,mag(:))
880 end if
881 ! HLL estimation of normal magnetic field at cell interfaces
882 bx(ixc^s)=(sr(ixc^s,index_v_mag)*br(ixc^s,ip1)-sl(ixc^s,index_v_mag)*bl(ixc^s,ip1)-&
883 flc(ixc^s,mag(ip1))-frc(ixc^s,mag(ip1)))/(sr(ixc^s,index_v_mag)-sl(ixc^s,index_v_mag))
884 ptr(ixc^s)=wrp(ixc^s,p_)+0.5d0*sum(br(ixc^s,:)**2,dim=ndim+1)
885 ptl(ixc^s)=wlp(ixc^s,p_)+0.5d0*sum(bl(ixc^s,:)**2,dim=ndim+1)
886 sur(ixc^s)=(sr(ixc^s,index_v_mag)-vrc(ixc^s,ip1))*wrc(ixc^s,rho_)
887 sul(ixc^s)=(sl(ixc^s,index_v_mag)-vlc(ixc^s,ip1))*wlc(ixc^s,rho_)
888 ! Miyoshi equation (38) and Guo euqation (20)
889 sm(ixc^s)=(sur(ixc^s)*vrc(ixc^s,ip1)-sul(ixc^s)*vlc(ixc^s,ip1)-&
890 thetasm*(ptr(ixc^s)-ptl(ixc^s)) )/(sur(ixc^s)-sul(ixc^s))
891 ! Miyoshi equation (39) and Guo euqation (28)
892 w1r(ixc^s,mom(ip1))=sm(ixc^s)
893 w1l(ixc^s,mom(ip1))=sm(ixc^s)
894 w2r(ixc^s,mom(ip1))=sm(ixc^s)
895 w2l(ixc^s,mom(ip1))=sm(ixc^s)
896 ! Guo equation (22)
897 w1r(ixc^s,mag(ip1))=bx(ixc^s)
898 w1l(ixc^s,mag(ip1))=bx(ixc^s)
899 if(b0field) then
900 ptr(ixc^s)=wrp(ixc^s,p_)+0.5d0*sum(wrc(ixc^s,mag(:))**2,dim=ndim+1)
901 ptl(ixc^s)=wlp(ixc^s,p_)+0.5d0*sum(wlc(ixc^s,mag(:))**2,dim=ndim+1)
902 end if
903
904 ! Miyoshi equation (43) and Guo equation (27)
905 w1r(ixc^s,rho_)=sur(ixc^s)/(sr(ixc^s,index_v_mag)-sm(ixc^s))
906 w1l(ixc^s,rho_)=sul(ixc^s)/(sl(ixc^s,index_v_mag)-sm(ixc^s))
907
908 ip2=mod(ip1+1,ndir)
909 if(ip2==0) ip2=ndir
910 r1r(ixc^s)=sur(ixc^s)*(sr(ixc^s,index_v_mag)-sm(ixc^s))-bx(ixc^s)**2
911 where(r1r(ixc^s)/=0.d0)
912 r1r(ixc^s)=1.d0/r1r(ixc^s)
913 endwhere
914 r1l(ixc^s)=sul(ixc^s)*(sl(ixc^s,index_v_mag)-sm(ixc^s))-bx(ixc^s)**2
915 where(r1l(ixc^s)/=0.d0)
916 r1l(ixc^s)=1.d0/r1l(ixc^s)
917 endwhere
918 ! Miyoshi equation (44)
919 w1r(ixc^s,mom(ip2))=vrc(ixc^s,ip2)-bx(ixc^s)*br(ixc^s,ip2)*&
920 (sm(ixc^s)-vrc(ixc^s,ip1))*r1r(ixc^s)
921 w1l(ixc^s,mom(ip2))=vlc(ixc^s,ip2)-bx(ixc^s)*bl(ixc^s,ip2)*&
922 (sm(ixc^s)-vlc(ixc^s,ip1))*r1l(ixc^s)
923 ! partial solution for later usage
924 w1r(ixc^s,mag(ip2))=(sur(ixc^s)*(sr(ixc^s,index_v_mag)-vrc(ixc^s,ip1))-bx(ixc^s)**2)*r1r(ixc^s)
925 w1l(ixc^s,mag(ip2))=(sul(ixc^s)*(sl(ixc^s,index_v_mag)-vlc(ixc^s,ip1))-bx(ixc^s)**2)*r1l(ixc^s)
926 if(ndir==3) then
927 ip3=mod(ip1+2,ndir)
928 if(ip3==0) ip3=ndir
929 ! Miyoshi equation (46)
930 w1r(ixc^s,mom(ip3))=vrc(ixc^s,ip3)-bx(ixc^s)*br(ixc^s,ip3)*&
931 (sm(ixc^s)-vrc(ixc^s,ip1))*r1r(ixc^s)
932 w1l(ixc^s,mom(ip3))=vlc(ixc^s,ip3)-bx(ixc^s)*bl(ixc^s,ip3)*&
933 (sm(ixc^s)-vlc(ixc^s,ip1))*r1l(ixc^s)
934 ! Miyoshi equation (47)
935 w1r(ixc^s,mag(ip3))=br(ixc^s,ip3)*w1r(ixc^s,mag(ip2))
936 w1l(ixc^s,mag(ip3))=bl(ixc^s,ip3)*w1l(ixc^s,mag(ip2))
937 end if
938 ! Miyoshi equation (45)
939 w1r(ixc^s,mag(ip2))=br(ixc^s,ip2)*w1r(ixc^s,mag(ip2))
940 w1l(ixc^s,mag(ip2))=bl(ixc^s,ip2)*w1l(ixc^s,mag(ip2))
941 if(b0field) then
942 ! Guo equation (26)
943 w1r(ixc^s,mag(:))=w1r(ixc^s,mag(:))-block%B0(ixc^s,:,ip1)
944 w1l(ixc^s,mag(:))=w1l(ixc^s,mag(:))-block%B0(ixc^s,:,ip1)
945 end if
946 ! equation (48)
947 if(phys_energy) then
948 ! Guo equation (25) equivalent to Miyoshi equation (41)
949 w1r(ixc^s,p_)=(sur(ixc^s)*ptl(ixc^s) - sul(ixc^s)*ptr(ixc^s) +&
950 phipres * sur(ixc^s)*sul(ixc^s)*(vrc(ixc^s,ip1)-vlc(ixc^s,ip1)))/&
951 (sur(ixc^s)-sul(ixc^s))
952 w1l(ixc^s,p_)=w1r(ixc^s,p_)
953 if(b0field) then
954 ! Guo equation (32)
955 w1r(ixc^s,p_)=w1r(ixc^s,p_)+sum(block%B0(ixc^s,:,ip1)*(wrc(ixc^s,mag(:))-w1r(ixc^s,mag(:))),dim=ndim+1)
956 w1l(ixc^s,p_)=w1l(ixc^s,p_)+sum(block%B0(ixc^s,:,ip1)*(wlc(ixc^s,mag(:))-w1l(ixc^s,mag(:))),dim=ndim+1)
957 end if
958 ! Miyoshi equation (48) and main part of Guo euqation (31)
959 w1r(ixc^s,e_)=((sr(ixc^s,index_v_mag)-vrc(ixc^s,ip1))*wrc(ixc^s,e_)-ptr(ixc^s)*vrc(ixc^s,ip1)+&
960 w1r(ixc^s,p_)*sm(ixc^s)+bx(ixc^s)*(sum(vrc(ixc^s,:)*wrc(ixc^s,mag(:)),dim=ndim+1)-&
961 sum(w1r(ixc^s,mom(:))*w1r(ixc^s,mag(:)),dim=ndim+1)))/(sr(ixc^s,index_v_mag)-sm(ixc^s))
962 w1l(ixc^s,e_)=((sl(ixc^s,index_v_mag)-vlc(ixc^s,ip1))*wlc(ixc^s,e_)-ptl(ixc^s)*vlc(ixc^s,ip1)+&
963 w1l(ixc^s,p_)*sm(ixc^s)+bx(ixc^s)*(sum(vlc(ixc^s,:)*wlc(ixc^s,mag(:)),dim=ndim+1)-&
964 sum(w1l(ixc^s,mom(:))*w1l(ixc^s,mag(:)),dim=ndim+1)))/(sl(ixc^s,index_v_mag)-sm(ixc^s))
965 if(b0field) then
966 ! Guo equation (31)
967 w1r(ixc^s,e_)=w1r(ixc^s,e_)+(sum(w1r(ixc^s,mag(:))*block%B0(ixc^s,:,ip1),dim=ndim+1)*sm(ixc^s)-&
968 sum(wrc(ixc^s,mag(:))*block%B0(ixc^s,:,ip1),dim=ndim+1)*vrc(ixc^s,ip1))/(sr(ixc^s,index_v_mag)-sm(ixc^s))
969 w1l(ixc^s,e_)=w1l(ixc^s,e_)+(sum(w1l(ixc^s,mag(:))*block%B0(ixc^s,:,ip1),dim=ndim+1)*sm(ixc^s)-&
970 sum(wlc(ixc^s,mag(:))*block%B0(ixc^s,:,ip1),dim=ndim+1)*vlc(ixc^s,ip1))/(sl(ixc^s,index_v_mag)-sm(ixc^s))
971 end if
972 end if
973
974 ! Miyoshi equation (49) and Guo equation (35)
975 w2r(ixc^s,rho_)=w1r(ixc^s,rho_)
976 w2l(ixc^s,rho_)=w1l(ixc^s,rho_)
977 w2r(ixc^s,mag(ip1))=w1r(ixc^s,mag(ip1))
978 w2l(ixc^s,mag(ip1))=w1l(ixc^s,mag(ip1))
979
980 r1r(ixc^s)=sqrt(w1r(ixc^s,rho_))
981 r1l(ixc^s)=sqrt(w1l(ixc^s,rho_))
982 tmp(ixc^s)=1.d0/(r1r(ixc^s)+r1l(ixc^s))
983 signbx(ixc^s)=sign(1.d0,bx(ixc^s))
984 ! Miyoshi equation (51) and Guo equation (33)
985 s1r(ixc^s)=sm(ixc^s)+abs(bx(ixc^s))/r1r(ixc^s)
986 s1l(ixc^s)=sm(ixc^s)-abs(bx(ixc^s))/r1l(ixc^s)
987 ! Miyoshi equation (59) and Guo equation (41)
988 w2r(ixc^s,mom(ip2))=(r1l(ixc^s)*w1l(ixc^s,mom(ip2))+r1r(ixc^s)*w1r(ixc^s,mom(ip2))+&
989 (w1r(ixc^s,mag(ip2))-w1l(ixc^s,mag(ip2)))*signbx(ixc^s))*tmp(ixc^s)
990 w2l(ixc^s,mom(ip2))=w2r(ixc^s,mom(ip2))
991 ! Miyoshi equation (61) and Guo equation (43)
992 w2r(ixc^s,mag(ip2))=(r1l(ixc^s)*w1r(ixc^s,mag(ip2))+r1r(ixc^s)*w1l(ixc^s,mag(ip2))+&
993 r1l(ixc^s)*r1r(ixc^s)*(w1r(ixc^s,mom(ip2))-w1l(ixc^s,mom(ip2)))*signbx(ixc^s))*tmp(ixc^s)
994 w2l(ixc^s,mag(ip2))=w2r(ixc^s,mag(ip2))
995 if(ndir==3) then
996 ! Miyoshi equation (60) and Guo equation (42)
997 w2r(ixc^s,mom(ip3))=(r1l(ixc^s)*w1l(ixc^s,mom(ip3))+r1r(ixc^s)*w1r(ixc^s,mom(ip3))+&
998 (w1r(ixc^s,mag(ip3))-w1l(ixc^s,mag(ip3)))*signbx(ixc^s))*tmp(ixc^s)
999 w2l(ixc^s,mom(ip3))=w2r(ixc^s,mom(ip3))
1000 ! Miyoshi equation (62) and Guo equation (44)
1001 w2r(ixc^s,mag(ip3))=(r1l(ixc^s)*w1r(ixc^s,mag(ip3))+r1r(ixc^s)*w1l(ixc^s,mag(ip3))+&
1002 r1l(ixc^s)*r1r(ixc^s)*(w1r(ixc^s,mom(ip3))-w1l(ixc^s,mom(ip3)))*signbx(ixc^s))*tmp(ixc^s)
1003 w2l(ixc^s,mag(ip3))=w2r(ixc^s,mag(ip3))
1004 end if
1005 ! Miyoshi equation (63) and Guo equation (45)
1006 if(phys_energy) then
1007 w2r(ixc^s,e_)=w1r(ixc^s,e_)+r1r(ixc^s)*(sum(w1r(ixc^s,mom(:))*w1r(ixc^s,mag(:)),dim=ndim+1)-&
1008 sum(w2r(ixc^s,mom(:))*w2r(ixc^s,mag(:)),dim=ndim+1))*signbx(ixc^s)
1009 w2l(ixc^s,e_)=w1l(ixc^s,e_)-r1l(ixc^s)*(sum(w1l(ixc^s,mom(:))*w1l(ixc^s,mag(:)),dim=ndim+1)-&
1010 sum(w2l(ixc^s,mom(:))*w2l(ixc^s,mag(:)),dim=ndim+1))*signbx(ixc^s)
1011 end if
1012
1013 ! convert velocity to momentum
1014 do idir=1,ndir
1015 w1r(ixc^s,mom(idir))=w1r(ixc^s,mom(idir))*w1r(ixc^s,rho_)
1016 w1l(ixc^s,mom(idir))=w1l(ixc^s,mom(idir))*w1l(ixc^s,rho_)
1017 w2r(ixc^s,mom(idir))=w2r(ixc^s,mom(idir))*w2r(ixc^s,rho_)
1018 w2l(ixc^s,mom(idir))=w2l(ixc^s,mom(idir))*w2l(ixc^s,rho_)
1019 end do
1020
1021 ! get fluxes of intermedate states
1022 do iw=iws,iwe
1023 ! CT MHD does not need normal B flux
1024 if(stagger_grid .and. flux_type(idims, iw) == flux_tvdlf) cycle
1025 if(flux_type(idims, iw) == flux_special) then
1026 ! known flux (fLC=fRC) for normal B and psi_ in GLM method
1027 f1l(ixc^s,iw)=flc(ixc^s,iw)
1028 f1r(ixc^s,iw)=f1l(ixc^s,iw)
1029 f2l(ixc^s,iw)=f1l(ixc^s,iw)
1030 f2r(ixc^s,iw)=f1l(ixc^s,iw)
1031 else if(flux_type(idims, iw) == flux_hll) then
1032 ! using hll flux for tracers
1033 f1l(ixc^s,iw)=(sr(ixc^s,index_v_mag)*flc(ixc^s, iw)-sl(ixc^s,index_v_mag)*frc(ixc^s, iw) &
1034 +sr(ixc^s,index_v_mag)*sl(ixc^s,index_v_mag)*(wrc(ixc^s,iw)-wlc(ixc^s,iw)))/(sr(ixc^s,index_v_mag)-sl(ixc^s,index_v_mag))
1035 f1r(ixc^s,iw)=f1l(ixc^s,iw)
1036 f2l(ixc^s,iw)=f1l(ixc^s,iw)
1037 f2r(ixc^s,iw)=f1l(ixc^s,iw)
1038 else
1039 f1l(ixc^s,iw)=flc(ixc^s,iw)+sl(ixc^s,index_v_mag)*(w1l(ixc^s,iw)-wlc(ixc^s,iw))
1040 f1r(ixc^s,iw)=frc(ixc^s,iw)+sr(ixc^s,index_v_mag)*(w1r(ixc^s,iw)-wrc(ixc^s,iw))
1041 f2l(ixc^s,iw)=f1l(ixc^s,iw)+s1l(ixc^s)*(w2l(ixc^s,iw)-w1l(ixc^s,iw))
1042 f2r(ixc^s,iw)=f1r(ixc^s,iw)+s1r(ixc^s)*(w2r(ixc^s,iw)-w1r(ixc^s,iw))
1043 end if
1044 end do
1045
1046 ! Miyoshi equation (66) and Guo equation (46)
1047 {do ix^db=ixcmin^db,ixcmax^db\}
1048 if(sl(ix^d,index_v_mag)>0.d0) then
1049 fc(ix^d,iws:iwe,ip1)=flc(ix^d,iws:iwe)
1050 else if(s1l(ix^d)>=0.d0) then
1051 fc(ix^d,iws:iwe,ip1)=f1l(ix^d,iws:iwe)
1052 else if(sm(ix^d)>=0.d0) then
1053 fc(ix^d,iws:iwe,ip1)=f2l(ix^d,iws:iwe)
1054 else if(s1r(ix^d)>=0.d0) then
1055 fc(ix^d,iws:iwe,ip1)=f2r(ix^d,iws:iwe)
1056 else if(sr(ix^d,index_v_mag)>=0.d0) then
1057 fc(ix^d,iws:iwe,ip1)=f1r(ix^d,iws:iwe)
1058 else if(sr(ix^d,index_v_mag)<0.d0) then
1059 fc(ix^d,iws:iwe,ip1)=frc(ix^d,iws:iwe)
1060 end if
1061 {end do\}
1062
1063 end associate
1064 end subroutine get_riemann_flux_hlld_mag2
1065
1066 !> Calculate fast magnetosonic wave speed
1067 subroutine get_hlld2_modif_c(w,x,ixI^L,ixO^L,csound)
1069
1070 integer, intent(in) :: ixI^L, ixO^L
1071 double precision, intent(in) :: w(ixI^S, nw), x(ixI^S,1:ndim)
1072 double precision, intent(out):: csound(ixI^S)
1073 double precision :: cfast2(ixI^S), AvMinCs2(ixI^S), b2(ixI^S), kmax
1074 double precision :: inv_rho(ixO^S), gamma_A2(ixO^S)
1075 integer :: rho_, p_, e_, mom(1:ndir), mag(1:ndir)
1076
1077 rho_ = iw_rho
1078 mom(:) = iw_mom(:)
1079 mag(:) = iw_mag(:)
1080 p_ = iw_e
1081 e_ = iw_e
1082
1083 inv_rho=1.d0/w(ixo^s,rho_)
1084
1085 ! store |B|^2 in v
1086
1087 if (b0field) then
1088 b2(ixo^s) = sum((w(ixo^s, mag(:))+block%B0(ixo^s,:,b0i))**2, dim=ndim+1)
1089 else
1090 b2(ixo^s) = sum(w(ixo^s, mag(:))**2, dim=ndim+1)
1091 end if
1092
1093
1094 if (b0field) then
1095 avmincs2= w(ixo^s, mag(idims))+block%B0(ixo^s,idims,b0i)
1096 else
1097 avmincs2= w(ixo^s, mag(idims))
1098 end if
1099
1100
1101 csound(ixo^s) = sum(w(ixo^s, mom(:))**2, dim=ndim+1)
1102
1103 cfast2(ixo^s) = b2(ixo^s) * inv_rho+csound(ixo^s)
1104 avmincs2(ixo^s) = cfast2(ixo^s)**2-4.0d0*csound(ixo^s) &
1105 * avmincs2(ixo^s)**2 &
1106 * inv_rho
1107
1108 where(avmincs2(ixo^s)<zero)
1109 avmincs2(ixo^s)=zero
1110 end where
1111
1112 avmincs2(ixo^s)=sqrt(avmincs2(ixo^s))
1113
1114 csound(ixo^s) = sqrt(half*(cfast2(ixo^s)+avmincs2(ixo^s)))
1115
1116 end subroutine get_hlld2_modif_c
1117
1118 end subroutine finite_volume
1119
1120 !> Determine the upwinded wLC(ixL) and wRC(ixR) from w.
1121 !> the wCT is only used when PPM is exploited.
1122 subroutine reconstruct_lr(ixI^L,ixL^L,ixR^L,idims,w,wLC,wRC,wLp,wRp,x,dxdim)
1123 use mod_physics
1125 use mod_limiter
1126 use mod_comm_lib, only: mpistop
1127
1128 integer, intent(in) :: ixi^l, ixl^l, ixr^l, idims
1129 double precision, intent(in) :: dxdim
1130 ! cell center w in primitive form
1131 double precision, dimension(ixI^S,1:nw) :: w
1132 ! left and right constructed status in conservative form
1133 double precision, dimension(ixI^S,1:nw) :: wlc, wrc
1134 ! left and right constructed status in primitive form
1135 double precision, dimension(ixI^S,1:nw) :: wlp, wrp
1136 double precision, dimension(ixI^S,1:ndim) :: x
1137
1138 integer :: jxr^l, ixc^l, jxc^l, iw
1139 double precision :: ldw(ixi^s), rdw(ixi^s), dwc(ixi^s)
1140
1141 select case (type_limiter(block%level))
1142 case (limiter_mp5)
1143 call mp5limiter(ixi^l,ixl^l,idims,w,wlp,wrp)
1144 case (limiter_weno3)
1145 call weno3limiter(ixi^l,ixl^l,idims,dxdim,w,wlp,wrp,1)
1146 case (limiter_wenoyc3)
1147 call weno3limiter(ixi^l,ixl^l,idims,dxdim,w,wlp,wrp,2)
1148 case (limiter_weno5)
1149 call weno5limiter(ixi^l,ixl^l,idims,dxdim,w,wlp,wrp,1)
1150 case (limiter_weno5nm)
1151 call weno5nmlimiter(ixi^l,ixl^l,idims,dxdim,w,wlp,wrp,1)
1152 case (limiter_wenoz5)
1153 call weno5limiter(ixi^l,ixl^l,idims,dxdim,w,wlp,wrp,2)
1154 case (limiter_wenoz5nm)
1155 call weno5nmlimiter(ixi^l,ixl^l,idims,dxdim,w,wlp,wrp,2)
1156 case (limiter_wenozp5)
1157 call weno5limiter(ixi^l,ixl^l,idims,dxdim,w,wlp,wrp,3)
1158 case (limiter_wenozp5nm)
1159 call weno5nmlimiter(ixi^l,ixl^l,idims,dxdim,w,wlp,wrp,3)
1160 case (limiter_weno5cu6)
1161 call weno5cu6limiter(ixi^l,ixl^l,idims,w,wlp,wrp)
1162 case (limiter_teno5ad)
1163 call teno5adlimiter(ixi^l,ixl^l,idims,dxdim,w,wlp,wrp)
1164 case (limiter_weno7)
1165 call weno7limiter(ixi^l,ixl^l,idims,w,wlp,wrp,1)
1166 case (limiter_mpweno7)
1167 call weno7limiter(ixi^l,ixl^l,idims,w,wlp,wrp,2)
1168 case (limiter_venk)
1169 call venklimiter(ixi^l,ixl^l,idims,dxdim,w,wlp,wrp)
1170 if(fix_small_values) then
1171 call phys_handle_small_values(.true.,wlp,x,ixi^l,ixl^l,'reconstruct left')
1172 call phys_handle_small_values(.true.,wrp,x,ixi^l,ixr^l,'reconstruct right')
1173 end if
1174 case (limiter_ppm)
1175 ixcmin^d=ixlmin^d+kr(idims,^d);
1176 ixcmax^d=ixlmax^d;
1177 call ppmlimiter(ixi^l,ixc^l,idims,w,w,wlp,wrp)
1178 if(fix_small_values) then
1179 call phys_handle_small_values(.true.,wlp,x,ixi^l,ixc^l,'reconstruct left')
1180 call phys_handle_small_values(.true.,wrp,x,ixi^l,ixc^l,'reconstruct right')
1181 end if
1182 case default
1183 jxr^l=ixr^l+kr(idims,^d);
1184 ixcmax^d=jxrmax^d; ixcmin^d=ixlmin^d-kr(idims,^d);
1185 jxc^l=ixc^l+kr(idims,^d);
1186 do iw=1,nwflux
1187 if (loglimit(iw)) then
1188 w(ixcmin^d:jxcmax^d,iw)=dlog10(w(ixcmin^d:jxcmax^d,iw))
1189 wlp(ixl^s,iw)=dlog10(wlp(ixl^s,iw))
1190 wrp(ixr^s,iw)=dlog10(wrp(ixr^s,iw))
1191 end if
1192
1193 dwc(ixc^s)=w(jxc^s,iw)-w(ixc^s,iw)
1194
1195 ! limit flux from left and/or right
1196 call dwlimiter2(dwc,ixi^l,ixc^l,idims,type_limiter(block%level),ldw,rdw)
1197 wlp(ixl^s,iw)=wlp(ixl^s,iw)+half*ldw(ixl^s)
1198 wrp(ixr^s,iw)=wrp(ixr^s,iw)-half*rdw(jxr^s)
1199
1200 if (loglimit(iw)) then
1201 w(ixcmin^d:jxcmax^d,iw)=10.0d0**w(ixcmin^d:jxcmax^d,iw)
1202 wlp(ixl^s,iw)=10.0d0**wlp(ixl^s,iw)
1203 wrp(ixr^s,iw)=10.0d0**wrp(ixr^s,iw)
1204 end if
1205 end do
1206 if(fix_small_values) then
1207 call phys_handle_small_values(.true.,wlp,x,ixi^l,ixl^l,'reconstruct left')
1208 call phys_handle_small_values(.true.,wrp,x,ixi^l,ixr^l,'reconstruct right')
1209 end if
1210 end select
1211
1212 wlc(ixl^s,1:nwflux)=wlp(ixl^s,1:nwflux)
1213 wrc(ixr^s,1:nwflux)=wrp(ixr^s,1:nwflux)
1214 call phys_to_conserved(ixi^l,ixl^l,wlc,x)
1215 call phys_to_conserved(ixi^l,ixr^l,wrc,x)
1216 if(nwaux>0)then
1217 wlp(ixl^s,nwflux+1:nwflux+nwaux)=wlc(ixl^s,nwflux+1:nwflux+nwaux)
1218 wrp(ixr^s,nwflux+1:nwflux+nwaux)=wrc(ixr^s,nwflux+1:nwflux+nwaux)
1219 endif
1220
1221 end subroutine reconstruct_lr
1222
1223end module mod_finite_volume
subroutine get_riemann_flux_tvdmu()
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Module with finite volume methods for fluxes.
subroutine, public finite_volume(method, qdt, dtfactor, ixil, ixol, idimslim, qtc, sct, qt, snew, fc, fe, dxs, x)
finite volume method
subroutine, public reconstruct_lr(ixil, ixll, ixrl, idims, w, wlc, wrc, wlp, wrp, x, dxdim)
Determine the upwinded wLC(ixL) and wRC(ixR) from w. the wCT is only used when PPM is exploited.
subroutine, public hancock(qdt, dtfactor, ixil, ixol, idimslim, qtc, sct, qt, snew, dxs, x)
The non-conservative Hancock predictor for TVDLF.
This module contains definitions of global parameters and variables and some generic functions/subrou...
type(state), pointer block
Block pointer for using one block and its previous state.
integer, dimension(3, 3) kr
Kronecker delta tensor.
logical, dimension(:), allocatable loglimit
integer, parameter ndim
Number of spatial dimensions for grid variables.
integer b0i
background magnetic field location indicator
double precision, dimension(:), allocatable, parameter d
logical local_timestep
each cell has its own timestep or not
logical slab
Cartesian geometry or not.
logical b0field
split magnetic field as background B0 field
integer, dimension(:), allocatable type_limiter
Type of slope limiter used for reconstructing variables on cell edges.
logical fix_small_values
fix small values with average or replace methods
logical slab_uniform
uniform Cartesian geometry or not (stretched Cartesian)
Module with slope/flux limiters.
Definition mod_limiter.t:2
integer, parameter limiter_weno5cu6
Definition mod_limiter.t:35
integer, parameter limiter_mpweno7
Definition mod_limiter.t:38
integer, parameter limiter_teno5ad
Definition mod_limiter.t:36
integer, parameter limiter_weno3
Definition mod_limiter.t:27
integer, parameter limiter_ppm
Definition mod_limiter.t:25
subroutine dwlimiter2(dwc, ixil, ixcl, idims, typelim, ldw, rdw)
Limit the centered dwC differences within ixC for iw in direction idim. The limiter is chosen accordi...
integer, parameter limiter_wenozp5
Definition mod_limiter.t:33
integer, parameter limiter_weno5
Definition mod_limiter.t:29
integer, parameter limiter_wenoz5
Definition mod_limiter.t:31
integer, parameter limiter_wenoz5nm
Definition mod_limiter.t:32
integer, parameter limiter_weno7
Definition mod_limiter.t:37
integer, parameter limiter_wenozp5nm
Definition mod_limiter.t:34
integer, parameter limiter_mp5
Definition mod_limiter.t:26
integer, parameter limiter_venk
Definition mod_limiter.t:23
integer, parameter limiter_weno5nm
Definition mod_limiter.t:30
integer, parameter limiter_wenoyc3
Definition mod_limiter.t:28
This module defines the procedures of a physics module. It contains function pointers for the various...
Definition mod_physics.t:4
procedure(sub_convert), pointer phys_to_primitive
Definition mod_physics.t:54
procedure(sub_small_values), pointer phys_handle_small_values
Definition mod_physics.t:76
procedure(sub_get_flux), pointer phys_get_flux
Definition mod_physics.t:61
procedure(sub_add_source_geom), pointer phys_add_source_geom
Definition mod_physics.t:65
procedure(sub_convert), pointer phys_to_conserved
Definition mod_physics.t:53
Module for handling split source terms (split from the fluxes)
Definition mod_source.t:2
subroutine, public addsource2(qdt, dtfactor, ixil, ixol, iwlim, qtc, wct, wctprim, qt, w, x, qsourcesplit, src_active)
Add source within ixO for iws: w=w+qdt*S[wCT].
Definition mod_source.t:132
Subroutines for TVD-MUSCL schemes.
Definition mod_tvd.t:2
subroutine, public tvdlimit2(method, qdt, ixil, ixicl, ixol, idims, wl, wr, wnew, x, fc, dxs)
Definition mod_tvd.t:39
Module with all the methods that users can customize in AMRVAC.
integer nwflux
Number of flux variables.