8 logical :: fix_conserve_at_step = .true.
22 integer,
intent(in) :: iit
24 integer :: iigrid, igrid, idimsplit
33 call advect(idimsplit,idimsplit)
38 do idimsplit=
ndim,1,-1
39 call advect(idimsplit,idimsplit)
63 integer,
intent(in) :: idim^LIM
64 integer :: iigrid, igrid
71 do iigrid=1,igridstail; igrid=igrids(iigrid);
72 ps1(igrid)%w=ps(igrid)%w
92 do iigrid=1,igridstail; igrid=igrids(iigrid);
93 ps1(igrid)%w=ps(igrid)%w
100 call mpistop(
"unkown onestep time_integrator in advect")
108 fix_conserve_at_step = .false.
118 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
119 ps(igrid)%w = ps(igrid)%w+
rk_b1*(ps1(igrid)%w-ps(igrid)%w)/
rk_a21
129 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
130 ps(igrid)%w = half*ps(igrid)%w+half*ps1(igrid)%w
131 if(
stagger_grid) ps(igrid)%ws = half*ps(igrid)%ws+half*ps1(igrid)%ws
138 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
139 ps2(igrid)%w = ps(igrid)%w
146 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
147 ps(igrid)%w = ps(igrid)%w+2.0d0*(ps2(igrid)%w-ps1(igrid)%w)
148 if(
stagger_grid) ps(igrid)%ws = ps(igrid)%ws+2.0d0*(ps2(igrid)%ws-ps1(igrid)%ws)
156 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
157 ps2(igrid)%w = half*(ps(igrid)%w+ps1(igrid)%w)
158 if(
stagger_grid) ps2(igrid)%ws = half*(ps(igrid)%ws+ps1(igrid)%ws)
163 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
164 ps1(igrid)%w = ps1(igrid)%w+half*
dt*ps(igrid)%w
165 if(
stagger_grid) ps1(igrid)%ws = ps1(igrid)%ws+half*
dt*ps(igrid)%ws
169 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
170 ps(igrid)%w = ps2(igrid)%w+half*
dt*ps(igrid)%w
171 if(
stagger_grid) ps(igrid)%ws = ps2(igrid)%ws+half*
dt*ps(igrid)%ws
177 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
178 ps(igrid)%w = ps(igrid)%w+ps2(igrid)%w-ps1(igrid)%w
179 if(
stagger_grid) ps(igrid)%ws = ps(igrid)%ws+ps2(igrid)%ws-ps1(igrid)%ws
193 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
194 ps2(igrid)%w = ps(igrid)%w
207 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
213 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
214 ps(igrid)%w = half*(ps(igrid)%w + ps1(igrid)%w + ps2(igrid)%w)
215 if(
stagger_grid) ps(igrid)%ws = half*(ps(igrid)%ws + ps1(igrid)%ws + ps2(igrid)%ws)
219 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
220 ps1(igrid)%w = ps1(igrid)%w + (1.0d0 - 2.0d0*
imex222_lambda)*ps2(igrid)%w
228 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
229 ps2(igrid)%w = 2.0d0*ps2(igrid)%w - ps1(igrid)%w -
imex222_lambda*ps2(igrid)%w
238 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
239 ps(igrid)%w = ps(igrid)%w + (ps2(igrid)%w - ps1(igrid)%w) / (2.0d0 *
imex222_lambda)
247 call mpistop(
"unkown twostep time_integrator in advect")
256 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
263 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
274 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
275 ps3(igrid)%w=(ps1(igrid)%w-ps(igrid)%w)/
rk3_a21
280 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
281 ps2(igrid)%w=ps(igrid)%w+
rk3_a31*ps3(igrid)%w
287 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
288 ps(igrid)%w=ps(igrid)%w+
rk3_b1*ps3(igrid)%w &
291 ps(igrid)%ws=ps(igrid)%ws+
rk3_b1*ps3(igrid)%ws &
302 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
303 ps4(igrid)%w=(ps1(igrid)%w-ps(igrid)%w)/
ars_gamma
309 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
310 ps1(igrid)%w=(ps2(igrid)%w-ps1(igrid)%w)/
ars_gamma
315 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
316 ps3(igrid)%w=ps(igrid)%w+(
ars_gamma-1.0d0)*ps4(igrid)%w+(1.0d0-2.0d0*
ars_gamma)*ps1(igrid)%w
318 ps3(igrid)%ws=ps(igrid)%ws+(
ars_gamma-1.0d0)*ps4(igrid)%ws+(1.0d0-2.0d0*
ars_gamma)*ps1(igrid)%ws
324 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
325 ps2(igrid)%w=ps1(igrid)%w+(ps3(igrid)%w-(ps(igrid)%w+ &
328 ps2(igrid)%ws=ps1(igrid)%ws+(ps3(igrid)%ws-(ps(igrid)%ws+ &
335 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
336 ps(igrid)%w=ps(igrid)%w+half*ps2(igrid)%w &
337 +half*(ps4(igrid)%w-ps3(igrid)%w)/
ars_gamma
339 ps(igrid)%ws=ps(igrid)%ws+half*ps2(igrid)%ws &
340 +half*(ps4(igrid)%ws-ps3(igrid)%ws)/
ars_gamma
350 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
351 ps4(igrid)%w=(ps1(igrid)%w-ps(igrid)%w)/
imex_a21
352 ps3(igrid)%w=ps(igrid)%w
354 ps4(igrid)%ws=(ps1(igrid)%ws-ps(igrid)%ws)/
imex_a21
355 ps3(igrid)%ws=ps(igrid)%ws
361 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
362 ps1(igrid)%w=ps1(igrid)%w+
imex_ha21*
dt*ps3(igrid)%w
369 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
370 ps(igrid)%w=ps(igrid)%w+
imex_a31*ps4(igrid)%w &
373 ps(igrid)%ws=ps(igrid)%ws+
imex_a31*ps4(igrid)%ws &
379 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
380 ps3(igrid)%w=ps1(igrid)%w-
imex_a21*ps4(igrid)%w &
383 ps3(igrid)%ws=ps1(igrid)%ws-
imex_a21*ps4(igrid)%ws &
390 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
394 ps2(igrid)%ws=(ps(igrid)%ws-ps3(igrid)%ws-
imex_a31*ps4(igrid)%ws)/
imex_a32 &
400 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
401 ps1(igrid)%w=ps3(igrid)%w+
imex_b1*ps4(igrid)%w+
imex_b2*ps2(igrid)%w
403 ps1(igrid)%ws=ps3(igrid)%ws+
imex_b1*ps4(igrid)%ws+
imex_b2*ps2(igrid)%ws
409 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
410 ps(igrid)%w=ps1(igrid)%w+ps2(igrid)%w-ps(igrid)%w
412 ps(igrid)%ws=ps1(igrid)%ws+ps2(igrid)%ws-ps(igrid)%ws
427 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
428 ps3(igrid)%w = ps(igrid)%w +
imex_a32/
imex_a22 * (ps2(igrid)%w - ps1(igrid)%w)
429 ps(igrid)%w = ps(igrid)%w +
imex_b2 /
imex_a22 * (ps2(igrid)%w - ps1(igrid)%w)
430 ps1(igrid)%w = ps3(igrid)%w
438 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
439 ps(igrid)%w = ps(igrid)%w +
imex_b2 /
imex_ha32 * (ps3(igrid)%w - ps1(igrid)%w)
444 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
445 ps(igrid)%w = ps(igrid)%w +
imex_b3 /
imex_a33 * (ps1(igrid)%w - ps3(igrid)%w)
452 call mpistop(
"unkown threestep time_integrator in advect")
465 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
472 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
479 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
489 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
490 ps2(igrid)%w=ps(igrid)%w
491 ps3(igrid)%w=ps(igrid)%w
493 ps2(igrid)%ws=ps(igrid)%ws
494 ps3(igrid)%ws=ps(igrid)%ws
502 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
503 ps(igrid)%w=(1.0d0/3.0d0)*(-ps(igrid)%w+ps1(igrid)%w+2.0d0*ps2(igrid)%w+ps3(igrid)%w)
505 *(-ps(igrid)%ws+ps1(igrid)%ws+2.0d0*ps2(igrid)%ws+ps3(igrid)%ws)
511 call mpistop(
"unkown fourstep time_integrator in advect")
521 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
528 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
535 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
541 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
548 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
549 ps(igrid)%w=ps3(igrid)%w+
rk_alfa55*ps2(igrid)%w &
552 ps(igrid)%ws=ps3(igrid)%ws+
rk_alfa55*ps2(igrid)%ws &
561 call mpistop(
"unkown fivestep time_integrator in advect")
565 call mpistop(
"unkown time_stepper in advect")
576 type(state),
target :: psa(max_blocks)
577 type(state),
target :: psb(max_blocks)
578 double precision,
intent(in) :: qdt
579 double precision,
intent(in) :: qtC
580 double precision,
intent(in) :: dtfactor
582 integer :: iigrid, igrid
586 do iigrid=1,igridstail; igrid=igrids(iigrid);
587 psa(igrid)%w = psb(igrid)%w
604 type(state),
target :: psa(max_blocks)
605 double precision,
intent(in) :: qtC
614 subroutine advect1(method,dtfactor,idim^LIM,qtC,psa,qt,psb)
620 integer,
intent(in) :: idim^LIM
621 type(state),
target :: psa(max_blocks)
622 type(state),
target :: psb(max_blocks)
623 double precision,
intent(in) :: dtfactor
624 double precision,
intent(in) :: qtC
625 double precision,
intent(in) :: qt
626 integer,
intent(in) :: method(nlevelshi)
629 double precision :: fC(ixG^T,1:nwflux,1:ndim)
631 double precision :: fE(ixG^T,sdim:3)
632 double precision :: qdt
633 integer :: iigrid, igrid
644 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
649 qtc,psa(igrid),qt,psb(igrid),fc,fe,
rnode(rpdx1_:
rnodehi,igrid),ps(igrid)%x)
676 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
689 subroutine advect1_grid(method,qdt,dtfactor,ixI^L,idim^LIM,qtC,sCT,qt,s,fC,fE,dxs,x)
700 integer,
intent(in) :: method
701 integer,
intent(in) :: ixI^L, idim^LIM
702 double precision,
intent(in) :: qdt, dtfactor, qtC, qt, dxs(ndim), x(ixI^S,1:ndim)
703 type(state),
target :: sCT, s
704 double precision :: fC(ixI^S,1:nwflux,1:ndim), wprim(ixI^S,1:nw)
705 double precision :: fE(ixI^S,sdim:3)
709 ixo^l=ixi^l^lsubnghostcells;
712 call finite_volume(method,qdt,dtfactor,ixi^l,ixo^l,idim^lim,qtc,sct,qt,s,fc,fe,dxs,x)
714 call centdiff(method,qdt,dtfactor,ixi^l,ixo^l,idim^lim,qtc,sct,qt,s,fc,fe,dxs,x)
716 call hancock(qdt,dtfactor,ixi^l,ixo^l,idim^lim,qtc,sct,qt,s,dxs,x)
718 call fd(qdt,dtfactor,ixi^l,ixo^l,idim^lim,qtc,sct,qt,s,fc,fe,dxs,x)
720 call centdiff(
fs_cd,qdt,dtfactor,ixi^l,ixo^l,idim^lim,qtc,sct,qt,s,fc,fe,dxs,x)
721 call tvdlimit(method,qdt,ixi^l,ixo^l,idim^lim,sct,qt+qdt,s,fc,dxs,x)
725 call addsource2(qdt*dble(idimmax-idimmin+1)/dble(ndim),&
726 dtfactor*dble(idimmax-idimmin+1)/dble(ndim),&
727 ixi^l,ixo^l,1,nw,qtc,sct%w,wprim,qt,s%w,x,.false.)
731 call mpistop(
"unknown flux scheme in advect1_grid")
745 integer,
intent(in) :: iit
746 double precision,
intent(in):: qt
748 integer:: iigrid, igrid
756 do iigrid=1,igridstail; igrid=igrids(iigrid);
761 qt,ps(igrid)%w,ps(igrid)%x)
779 integer,
intent(in) :: iit
780 double precision,
intent(in):: qt
782 integer:: iigrid, igrid
790 do iigrid=1,igridstail; igrid=igrids(iigrid);
796 qt,ps(igrid)%w,ps(igrid)%x)
Module containing all the time stepping schemes.
subroutine, public process_advanced(iit, qt)
process_advanced is user entry in time loop, just after advance allows to modify solution,...
subroutine, public advance(iit)
Advance all the grids over one time step, including all sources.
subroutine advect1(method, dtfactor, idimLIM, qtC, psa, qt, psb)
Integrate all grids by one partial step.
subroutine, public process(iit, qt)
process is a user entry in time loop, before output and advance allows to modify solution,...
subroutine evaluate_implicit(qtC, psa)
Evaluate Implicit part in place, i.e. psa==>F_im(psa)
subroutine global_implicit_update(dtfactor, qdt, qtC, psa, psb)
Implicit global update step within IMEX schemes, advance psa=psb+dtfactor*qdt*F_im(psa)
subroutine advect1_grid(method, qdt, dtfactor, ixIL, idimLIM, qtC, sCT, qt, s, fC, fE, dxs, x)
Advance a single grid over one partial time step.
subroutine advect(idimLIM)
Advance all grids over one time step, but without taking dimensional splitting or split source terms ...
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Module with finite difference methods for fluxes.
subroutine, public centdiff(method, qdt, dtfactor, ixIL, ixOL, idimsLIM, qtC, sCT, qt, s, fC, fE, dxs, x)
subroutine, public fd(qdt, dtfactor, ixIL, ixOL, idimsLIM, qtC, sCT, qt, snew, fC, fE, dxs, x)
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 hancock(qdt, dtfactor, ixIL, ixOL, idimsLIM, qtC, sCT, qt, snew, dxs, x)
The non-conservative Hancock predictor for TVDLF.
Module for flux conservation near refinement boundaries.
subroutine, public recvflux(idimLIM)
subroutine, public fix_edges(psuse, idimLIM)
subroutine, public sendflux(idimLIM)
subroutine, public fix_conserve(psb, idimLIM, nw0, nwfluxin)
subroutine, public store_edge(igrid, ixIL, fE, idimLIM)
subroutine, public init_comm_fix_conserve(idimLIM, nwfluxin)
subroutine, public store_flux(igrid, fC, idimLIM, nwfluxin)
update ghost cells of all blocks including physical boundaries
subroutine getbc(time, qdt, psb, nwstart, nwbc, req_diag)
do update ghost cells of all blocks including physical boundaries
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, parameter ssprk4
integer, dimension(:), allocatable typepred1
The spatial discretization for the predictor step when using a two step PC method.
integer, parameter fs_tvdlf
double precision rk_alfa41
double precision imex_a32
integer, parameter fs_hlld
integer, parameter imex_euler
integer, parameter fs_tvd
double precision rk_beta55
double precision global_time
The global simulation time.
double precision rk_alfa22
integer istep
Index of the sub-step in a multi-step time integrator.
integer, parameter threestep
integer, parameter imex_232
double precision rk_alfa55
double precision rk_beta22
integer, parameter fivestep
double precision imex_ha32
double precision ars_gamma
IMEX_ARS3 parameter ars_gamma.
integer, parameter fs_hllcd
integer, parameter ndim
Number of spatial dimensions for grid variables.
integer, parameter ssprk5
logical stagger_grid
True for using stagger grid.
integer, parameter imex_ars3
logical use_particles
Use particles module or not.
integer, parameter fs_tvdmu
double precision rk_alfa54
integer, parameter imex_trapezoidal
integer, parameter plevel_
double precision, dimension(:), allocatable, parameter d
double precision rk_alfa21
double precision dt
global time step
double precision imex222_lambda
IMEX-222(lambda) one-parameter family of schemes.
integer, dimension(:), allocatable flux_method
Which flux scheme of spatial discretization to use (per grid level)
integer, parameter fs_nul
double precision rk_beta54
integer, parameter rk2_alf
integer, parameter fs_hll
flux schemes
double precision rk_beta11
logical time_advance
do time evolving
double precision, dimension(:,:), allocatable rnode
Corner coordinates.
double precision rk_alfa31
integer, parameter imex_sp
double precision rk_alfa44
double precision imex_ha22
double precision imex_a31
integer, parameter twostep
integer, parameter onestep
double precision rk_a21
RK2(alfa) method parameters from Butcher tableau.
integer, parameter predictor_corrector
integer, parameter forward_euler
logical fix_conserve_global
Whether to apply flux conservation at refinement boundaries.
character(len=std_len) typedimsplit
double precision rk_alfa53
double precision imex_ha21
integer, parameter ssprk3
double precision, dimension(^nd) dxlevel
store unstretched cell size of current level
integer t_stepper
time stepper type
integer, parameter fs_cd4
integer, parameter rk3_bt
integer, parameter fourstep
double precision imex_a21
integer, parameter rnodehi
grid location info (corner coordinates and grid spacing)
integer, parameter fs_hllc
double precision imex_a22
IMEX_CB3a extra parameters.
integer, parameter imex_cb3a
integer, parameter imex_222
double precision imex_a33
integer t_integrator
time integrator method
integer, parameter fs_hancock
integer, dimension(:,:), allocatable node
integer, parameter ssprk2
double precision rk_alfa33
integer, parameter fs_source
integer, parameter imex_midpoint
double precision rk_beta33
double precision rk_beta44
Module containing all the particle routines.
This module defines the procedures of a physics module. It contains function pointers for the various...
procedure(sub_convert), pointer phys_to_primitive
procedure(sub_evaluate_implicit), pointer phys_evaluate_implicit
logical phys_req_diagonal
Whether the physics routines require diagonal ghost cells, for example for computing a curl.
procedure(sub_implicit_update), pointer phys_implicit_update
procedure(sub_face_to_center), pointer phys_face_to_center
procedure(sub_special_advance), pointer phys_special_advance
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].
subroutine, public add_split_source(prior)
Subroutines for TVD-MUSCL schemes.
subroutine, public tvdlimit(method, qdt, ixIL, ixOL, idimLIM, s, qt, snew, fC, dxs, x)
Module with all the methods that users can customize in AMRVAC.
procedure(process_grid), pointer usr_process_grid
procedure(process_adv_grid), pointer usr_process_adv_grid
procedure(process_global), pointer usr_process_global
procedure(process_adv_global), pointer usr_process_adv_global