12 subroutine fd(qdt,dtfactor,ixI^L,ixO^L,idims^LIM,qtC,sCT,qt,snew,fC,fE,dxs,x)
19 double precision,
intent(in) :: qdt, dtfactor, qtc, qt, dxs(
ndim)
20 integer,
intent(in) :: ixi^
l, ixo^
l, idims^lim
21 double precision,
dimension(ixI^S,1:ndim),
intent(in) :: x
23 type(state) :: sct, snew
24 double precision,
dimension(ixI^S,1:nwflux,1:ndim),
intent(out) :: fc
25 double precision,
dimension(ixI^S,sdim:3) :: fe
27 double precision,
dimension(ixI^S,1:nwflux) :: fct
28 double precision,
dimension(ixI^S,1:nw) :: fm, fp, fmr, fpl, wprim
30 double precision,
dimension(ixI^S,1:nw) :: wlc, wrc
32 double precision,
dimension(ixI^S,1:nw) :: wlp, wrp
33 double precision,
dimension(ixI^S) :: cmaxc, cminc
34 double precision,
dimension(ixI^S) :: hspeed
35 double precision,
dimension(ixO^S) :: inv_volume
36 double precision,
dimension(1:ndim) :: dxinv
37 logical :: transport, active
38 integer :: idims, iw, ixc^
l, ix^
l, hxo^
l, kxc^
l, kxr^
l
39 type(ct_velocity) :: vcts
41 associate(wct=>sct%w,wnew=>snew%w)
56 hxo^
l=ixo^
l-
kr(idims,^
d);
62 kxcmin^
d=iximin^
d; kxcmax^
d=iximax^
d-
kr(idims,^
d);
63 kxr^
l=kxc^
l+
kr(idims,^
d);
68 wrp(kxc^s,1:nw)=wprim(kxr^s,1:nw)
69 wlp(kxc^s,1:nw)=wprim(kxc^s,1:nw)
72 ixcmax^
d=ixomax^
d; ixcmin^
d=hxomin^
d;
84 call reconstructl(ixi^
l,ixc^
l,idims,fp,fpl)
85 call reconstructr(ixi^
l,ixc^
l,idims,fm,fmr)
87 fc(ixc^s,iwstart:nwflux,idims) = fpl(ixc^s,iwstart:nwflux) + fmr(ixc^s,iwstart:nwflux)
91 call reconstruct_lr(ixi^
l,ixc^
l,ixc^
l,idims,wprim,wlc,wrc,wlp,wrp,x,dxs(idims))
95 call phys_get_cbounds(wlc,wrc,wlp,wrp,x,ixi^
l,ixc^
l,idims,hspeed,cmaxc,cminc)
107 hxo^
l=ixo^
l-
kr(idims,^
d);
110 fc(ixi^s,iw,idims) = -
block%dt(ixi^s)*dtfactor/dxs(idims) * fc(ixi^s,iw,idims)
111 wnew(ixo^s,iw)=wnew(ixo^s,iw)+(fc(ixo^s,iw,idims)-fc(hxo^s,iw,idims))
115 fc(ixi^s,iw,idims) = dxinv(idims) * fc(ixi^s,iw,idims)
116 wnew(ixo^s,iw)=wnew(ixo^s,iw)+(fc(ixo^s,iw,idims)-fc(hxo^s,iw,idims))
121 inv_volume=1.d0/
block%dvolume(ixo^s)
123 hxo^
l=ixo^
l-
kr(idims,^
d);
126 fc(ixi^s,iw,idims)=-
block%dt(ixi^s)*dtfactor*fc(ixi^s,iw,idims)*
block%surfaceC(ixi^s,idims)
127 wnew(ixo^s,iw)=wnew(ixo^s,iw)+ &
128 (fc(ixo^s,iw,idims)-fc(hxo^s,iw,idims))*inv_volume(ixo^s)
132 fc(ixi^s,iw,idims)=-qdt*fc(ixi^s,iw,idims)*
block%surfaceC(ixi^s,idims)
133 wnew(ixo^s,iw)=wnew(ixo^s,iw)+ &
134 (fc(ixo^s,iw,idims)-fc(hxo^s,iw,idims))*inv_volume(ixo^s)
155 dtfactor*dble(idimsmax-idimsmin+1)/dble(
ndim), &
156 ixi^
l,ixo^
l,1,nw,qtc,wct,wprim,qt,wnew,x,.false.,active)
162 subroutine reconstructl(ixI^L,iL^L,idims,w,wLC)
168 integer,
intent(in) :: ixi^
l, il^
l, idims
169 double precision,
intent(in) :: w(ixi^s,1:nw)
171 double precision,
intent(out) :: wlc(ixi^s,1:nw)
173 double precision :: ldw(ixi^s), dwc(ixi^s)
174 integer :: jxr^
l, ixc^
l, jxc^
l, kxc^
l, iw
180 call weno5limiterl(ixi^
l,il^
l,idims,w,wlc,1)
182 call weno5nmlimiterl(ixi^
l,il^
l,idims,w,wlc,1)
184 call weno5limiterl(ixi^
l,il^
l,idims,w,wlc,2)
186 call weno5nmlimiterl(ixi^
l,il^
l,idims,w,wlc,2)
188 call weno5limiterl(ixi^
l,il^
l,idims,w,wlc,3)
190 call weno5nmlimiterl(ixi^
l,il^
l,idims,w,wlc,3)
193 kxcmin^
d=iximin^
d; kxcmax^
d=iximax^
d-
kr(idims,^
d);
195 wlc(kxc^s,iwstart:nwflux) = w(kxc^s,iwstart:nwflux)
197 jxr^
l=il^
l+
kr(idims,^
d);
199 ixcmax^
d=jxrmax^
d; ixcmin^
d=ilmin^
d-
kr(idims,^
d);
200 jxc^
l=ixc^
l+
kr(idims,^
d);
203 dwc(ixc^s)=w(jxc^s,iw)-w(ixc^s,iw)
207 wlc(il^s,iw)=wlc(il^s,iw)+half*ldw(il^s)
211 end subroutine reconstructl
213 subroutine reconstructr(ixI^L,iL^L,idims,w,wRC)
219 integer,
intent(in) :: ixi^
l, il^
l, idims
220 double precision,
intent(in) :: w(ixi^s,1:nw)
222 double precision,
intent(out) :: wrc(ixi^s,1:nw)
224 double precision :: rdw(ixi^s), dwc(ixi^s)
225 integer :: jxr^
l, ixc^
l, jxc^
l, kxc^
l, kxr^
l, iw
231 call weno5limiterr(ixi^
l,il^
l,idims,w,wrc,1)
233 call weno5nmlimiterr(ixi^
l,il^
l,idims,w,wrc,1)
235 call weno5limiterr(ixi^
l,il^
l,idims,w,wrc,2)
237 call weno5nmlimiterr(ixi^
l,il^
l,idims,w,wrc,2)
239 call weno5limiterr(ixi^
l,il^
l,idims,w,wrc,3)
241 call weno5nmlimiterr(ixi^
l,il^
l,idims,w,wrc,3)
244 kxcmin^
d=iximin^
d; kxcmax^
d=iximax^
d-
kr(idims,^
d);
245 kxr^
l=kxc^
l+
kr(idims,^
d);
247 wrc(kxc^s,iwstart:nwflux)=w(kxr^s,iwstart:nwflux)
249 jxr^
l=il^
l+
kr(idims,^
d);
250 ixcmax^
d=jxrmax^
d; ixcmin^
d=ilmin^
d-
kr(idims,^
d);
251 jxc^
l=ixc^
l+
kr(idims,^
d);
254 dwc(ixc^s)=w(jxc^s,iw)-w(ixc^s,iw)
258 wrc(il^s,iw)=wrc(il^s,iw)-half*rdw(jxr^s)
262 end subroutine reconstructr
264 subroutine centdiff(method,qdt,dtfactor,ixI^L,ixO^L,idims^LIM,qtC,sCT,qt,s,fC,fE,dxs,x)
279 integer,
intent(in) :: method
280 integer,
intent(in) :: ixi^
l, ixo^
l, idims^lim
281 double precision,
intent(in) :: qdt, dtfactor, qtc, qt, dxs(
ndim)
282 type(state) :: sct, s
283 double precision,
intent(in) :: x(ixi^s,1:
ndim)
285 double precision :: fe(ixi^s,
sdim:3)
287 double precision :: v(ixi^s,
ndim), f(ixi^s,
nwflux)
289 double precision,
dimension(ixI^S,1:nw) :: wprim, wlc, wrc
291 double precision,
dimension(ixI^S,1:nw) :: wlp, wrp
292 double precision,
dimension(ixI^S) :: vlc, cmaxlc, cmaxrc
293 double precision,
dimension(ixI^S,1:number_species) :: cmaxc
294 double precision,
dimension(ixI^S,1:number_species) :: cminc
295 double precision,
dimension(ixI^S) :: hspeed
296 double precision,
dimension(ixO^S) :: inv_volume
298 double precision :: dxinv(1:
ndim)
299 integer :: idims, iw, ix^
l, hxo^
l, ixc^
l, jxc^
l, hxc^
l, kxc^
l, kkxc^
l, kkxr^
l
300 type(ct_velocity) :: vcts
301 logical :: transport, new_cmax, patchw(ixi^s), active
303 associate(wct=>sct%w,w=>s%w)
307 ix^
l=ix^
l^ladd2*
kr(idims,^
d);
310 if (ixi^
l^ltix^
l|.or.|.or.)
then
311 call mpistop(
"Error in centdiff: Non-conforming input limits")
322 ix^
l=ixo^
l^ladd2*
kr(idims,^
d);
323 hxo^
l=ixo^
l-
kr(idims,^
d);
333 ixcmax^
d=ixomax^
d; ixcmin^
d=hxomin^
d;
335 hxc^
l=ixc^
l-
kr(idims,^
d);
336 jxc^
l=ixc^
l+
kr(idims,^
d);
337 kxc^
l=ixc^
l+2*
kr(idims,^
d);
339 kkxcmin^
d=iximin^
d; kkxcmax^
d=iximax^
d-
kr(idims,^
d);
340 kkxr^
l=kkxc^
l+
kr(idims,^
d);
345 call reconstruct_lr(ixi^
l,ixc^
l,ixc^
l,idims,wprim,wlc,wrc,wlp,wrp,x,dxs(idims))
351 call phys_get_cbounds(wlc,wrc,wlp,wrp,x,ixi^
l,ixc^
l,idims,hspeed,cmaxc,cminc)
359 vlc(ixc^s)=max(cmaxrc(ixc^s),cmaxlc(ixc^s))
364 if(method==
fs_cd)
then
373 fc(ixc^s,iw,idims)=fc(ixc^s,iw,idims)-
tvdlfeps*half*vlc(ixc^s) &
374 *(wrc(ixc^s,iw)-wlc(ixc^s,iw))
386 hxo^
l=ixo^
l-
kr(idims,^
d);
389 fc(ixi^s,iw,idims)=-
block%dt(ixi^s)*dtfactor/dxs(idims)*fc(ixi^s,iw,idims)
391 w(ixo^s,iw)=w(ixo^s,iw)+(fc(ixo^s,iw,idims)-fc(hxo^s,iw,idims))
395 fc(ixi^s,iw,idims)=dxinv(idims)*fc(ixi^s,iw,idims)
397 w(ixo^s,iw)=w(ixo^s,iw)+(fc(ixo^s,iw,idims)-fc(hxo^s,iw,idims))
402 inv_volume=1.d0/
block%dvolume
404 hxo^
l=ixo^
l-
kr(idims,^
d);
407 fc(ixi^s,iw,idims)=-
block%dt(ixi^s)*dtfactor*
block%surfaceC(ixi^s,idims)*fc(ixi^s,iw,idims)
408 w(ixo^s,iw)=w(ixo^s,iw)+ &
409 (fc(ixo^s,iw,idims)-fc(hxo^s,iw,idims))*inv_volume(ixo^s)
413 fc(ixi^s,iw,idims)=-qdt*
block%surfaceC(ixi^s,idims)*fc(ixi^s,iw,idims)
414 w(ixo^s,iw)=w(ixo^s,iw)+ &
415 (fc(ixo^s,iw,idims)-fc(hxo^s,iw,idims))*inv_volume(ixo^s)
431 dtfactor*dble(idimsmax-idimsmin+1)/dble(
ndim), &
432 ixi^
l,ixo^
l,1,
nw,qtc,wct,wprim,qt,w,x,.false.,active)
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Module with finite difference methods for fluxes.
subroutine, public fd(qdt, dtfactor, ixil, ixol, idimslim, qtc, sct, qt, snew, fc, fe, dxs, x)
subroutine, public centdiff(method, qdt, dtfactor, ixil, ixol, idimslim, qtc, sct, qt, s, fc, fe, dxs, x)
Module with finite volume methods for fluxes.
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.
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.
logical h_correction
If true, do H-correction to fix the carbuncle problem at grid-aligned shocks.
double precision tvdlfeps
integer, dimension(3, 3) kr
Kronecker delta tensor.
integer, parameter ndim
Number of spatial dimensions for grid variables.
logical stagger_grid
True for using stagger grid.
double precision cmax_global
global fastest wave speed needed in fd scheme and glm method
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.
integer, dimension(:), allocatable type_limiter
Type of slope limiter used for reconstructing variables on cell edges.
integer nghostcells
Number of ghost cells surrounding a grid.
integer, parameter sdim
starting dimension for electric field
logical fix_small_values
fix small values with average or replace methods
logical crash
Save a snapshot before crash a run met unphysical values.
logical slab_uniform
uniform Cartesian geometry or not (stretched Cartesian)
Module with slope/flux limiters.
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
integer, parameter limiter_weno5
integer, parameter limiter_wenoz5
integer, parameter limiter_wenoz5nm
integer, parameter limiter_wenozp5nm
integer, parameter limiter_mp5
integer, parameter limiter_weno5nm
Module containing the MP5 (fifth order) flux scheme.
subroutine, public mp5limiterl(ixil, ill, idims, w, wlc)
subroutine, public mp5limiterr(ixil, ill, idims, w, wrc)
This module defines the procedures of a physics module. It contains function pointers for the various...
procedure(sub_get_ct_velocity), pointer phys_get_ct_velocity
procedure(sub_convert), pointer phys_to_primitive
procedure(sub_small_values), pointer phys_handle_small_values
procedure(sub_get_flux), pointer phys_get_flux
procedure(sub_get_cbounds), pointer phys_get_cbounds
procedure(sub_add_source_geom), pointer phys_add_source_geom
procedure(sub_update_faces), pointer phys_update_faces
procedure(sub_face_to_center), pointer phys_face_to_center
procedure(sub_get_h_speed), pointer phys_get_h_speed
procedure(sub_get_cmax), pointer phys_get_cmax
Module for handling split source terms (split from the fluxes)
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].
Module with all the methods that users can customize in AMRVAC.
integer nw
Total number of variables.
integer nwflux
Number of flux variables.