14 integer :: igrid,iigrid
19 do iigrid=1,igridstail; igrid=igrids(iigrid);
37 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
54 integer :: igrid,idim^LIM
57 double precision :: dx^D
58 double precision :: fC(ixG^T,1:nwflux,1:ndim)
59 double precision :: fE(ixG^T,sdim:3)
76 double precision :: fC(ixI^S,1:nwflux,1:ndim)
77 double precision :: fE(ixI^S,sdim:3)
78 double precision :: dx^D
80 double precision :: A(s%ixGs^S,1:3)
81 integer :: ixIs^L,ixO^L,idir
83 associate(ws=>s%ws,x=>s%x)
89 ixo^l=ixi^l^lsubnghostcells;
96 fe(ixi^s,idir) =-a(ixi^s,idir)
107 integer,
intent(in) :: ixIs^L, ixI^L, ixO^L
108 double precision,
intent(inout) :: ws(ixIs^S,1:nws),A(ixIs^S,1:3)
109 double precision,
intent(in) :: x(ixI^S,1:ndim)
111 double precision :: xC(ixIs^S,1:ndim),xCC(ixIs^S,1:ndim)
112 double precision :: circ(ixIs^S,1:ndim)
113 integer :: ixC^L, hxC^L, idim1, idim2, idir
119 xcc(ixi^s,1:ndim)=x(ixi^s,1:ndim)
121 xcc(ixismin^
d^
d%ixI^s,1:ndim)=x(iximin^
d^
d%ixI^s,1:ndim)
122 xcc(ixismin^
d^
d%ixIs^s,^
d)=x({iximin^dd,},^
d)-
block%dx({iximin^dd,},^
d)
125 xcc(iximin1:iximax1,ixismin2,ixismin3,1)=x(iximin1:iximax1,iximin2,iximin3,1)
126 xcc(ixismin1,iximin2:iximax2,ixismin3,2)=x(iximin1,iximin2:iximax2,iximin3,2)
127 xcc(ixismin1,ixismin2,iximin3:iximax3,3)=x(iximin1,iximin2,iximin3:iximax3,3)
132 ixcmin^
d=ixomin^
d-1+
kr(idir,^
d);
135 if (idim1/=idir)
then
136 xc(ixc^s,idim1)=xcc(ixc^s,idim1)+half*
block%dx(ixc^s,idim1)
138 xc(ixc^s,idim1)=xcc(ixc^s,idim1)
143 a(ixc^s,idir)=a(ixc^s,idir)*
block%dsC(ixc^s,idir)
151 ixcmin^
d=ixomin^
d-
kr(idim1,^
d);
154 if(
lvc(idim1,idim2,idir)==0) cycle
156 hxc^l=ixc^l-
kr(idim2,^
d);
158 circ(ixc^s,idim1)=circ(ixc^s,idim1)&
159 +
lvc(idim1,idim2,idir)* &
165 where(
block%surfaceC(ixc^s,idim1)==0)
166 circ(ixc^s,idim1)=zero
168 circ(ixc^s,idim1)=circ(ixc^s,idim1)/
block%surfaceC(ixc^s,idim1)
170 ws(ixc^s,idim1) = circ(ixc^s,idim1)
181 integer,
intent(in) :: ixI^L, ixC^L, idir
182 double precision,
intent(in) :: q(ixI^S)
183 double precision,
intent(out) :: qL(ixI^S), qR(ixI^S)
185 double precision :: qC(ixI^S)
186 double precision,
dimension(ixI^S) :: dqC,ldq,rdq
187 integer :: ixO^L,jxC^L,gxC^L,hxC^L
189 jxc^l=ixc^l+
kr(idir,^
d);
190 gxcmin^
d=ixcmin^
d-
kr(idir,^
d);gxcmax^
d=jxcmax^
d;
191 hxc^l=gxc^l+
kr(idir,^
d);
199 ixomin^
d=ixcmin^
d+
kr(idir,^
d);
201 call ppmlimitervar(ixi^l,ixo^l,idir,q,q,ql,qr)
203 call mp5limitervar(ixi^l,ixc^l,idir,q,ql,qr)
205 call weno3limiter(ixi^l,ixc^l,idir,
dxlevel(idir),q,ql,qr,1)
207 call weno3limiter(ixi^l,ixc^l,idir,
dxlevel(idir),q,ql,qr,2)
209 call weno5limiter(ixi^l,ixc^l,idir,
dxlevel(idir),q,ql,qr,1)
211 call weno5nmlimiter(ixi^l,ixc^l,idir,
dxlevel(idir),q,ql,qr,1)
213 call weno5limiter(ixi^l,ixc^l,idir,
dxlevel(idir),q,ql,qr,2)
215 call weno5nmlimiter(ixi^l,ixc^l,idir,
dxlevel(idir),q,ql,qr,2)
217 call weno5limiter(ixi^l,ixc^l,idir,
dxlevel(idir),q,ql,qr,3)
219 call weno5nmlimiter(ixi^l,ixc^l,idir,
dxlevel(idir),q,ql,qr,3)
221 call weno5cu6limiter(ixi^l,ixc^l,idir,q,ql,qr)
223 call teno5adlimiter(ixi^l,ixc^l,idir,
dxlevel(idir),q,ql,qr)
225 call weno7limiter(ixi^l,ixc^l,idir,q,ql,qr,1)
227 call weno7limiter(ixi^l,ixc^l,idir,q,ql,qr,2)
229 call venklimiter(ixi^l,ixc^l,idir,
dxlevel(idir),q,ql,qr)
231 dqc(gxc^s)= qr(gxc^s)-ql(gxc^s)
233 ql(ixc^s) = ql(ixc^s) + half*ldq(ixc^s)
234 qr(ixc^s) = qr(ixc^s) - half*rdq(jxc^s)
subroutine recalculateb
re-calculate the magnetic field from the vector potential in a completely divergency free way
subroutine fake_update(ixIL, s, fC, fE, dxD)
fake update magnetic field from vector potential
subroutine b_from_vector_potentiala(ixIsL, ixIL, ixOL, ws, x, A)
calculate magnetic field from vector potential A at cell edges
subroutine reconstruct(ixIL, ixCL, idir, q, qL, qR)
Reconstruct scalar q within ixO^L to 1/2 dx in direction idir Return both left and right reconstructe...
subroutine fake_advance(igrid, idimLIM, s)
fake advance a step to calculate magnetic field
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, dimension(3, 3, 3) lvc
Levi-Civita tensor.
double precision global_time
The global simulation time.
integer, dimension(3, 3) kr
Kronecker delta tensor.
integer, parameter ndim
Number of spatial dimensions for grid variables.
double precision, dimension(:), allocatable, parameter d
double precision, dimension(:,:), allocatable rnode
Corner coordinates.
integer, dimension(:), allocatable type_limiter
Type of slope limiter used for reconstructing variables on cell edges.
integer, parameter sdim
starting dimension for electric field
double precision, dimension(^nd) dxlevel
store unstretched cell size of current level
Module with slope/flux limiters.
integer, parameter limiter_weno5cu6
integer, parameter limiter_mpweno7
integer, parameter limiter_teno5ad
integer, parameter limiter_weno3
integer, parameter limiter_ppm
subroutine dwlimiter2(dwC, ixIL, ixCL, idims, typelim, ldw, rdw, a2max)
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_weno7
integer, parameter limiter_wenozp5nm
integer, parameter limiter_mp5
integer, parameter limiter_venk
integer, parameter limiter_weno5nm
integer, parameter limiter_wenoyc3
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_convert), pointer phys_to_conserved
procedure(sub_face_to_center), pointer phys_face_to_center
Module with all the methods that users can customize in AMRVAC.
procedure(init_vector_potential), pointer usr_init_vector_potential