12 subroutine bc_phys(iside,idims,time,qdt,s,ixG^L,ixB^L)
19 integer,
intent(in) :: iside, idims, ixg^
l,ixb^
l
20 double precision,
intent(in) :: time,qdt
21 type(state),
intent(inout) :: s
22 double precision :: wtmp(ixg^s,1:nwflux)
25 integer :: ixos^
l,hxo^
l,jxo^
l
26 double precision :: q(ixg^s),qp(ixg^s)
27 integer :: iw, ib, ix^
d, ixo^
l,
ixm^
l, nghostcellsi,iib^
d
28 logical :: isphysbound
30 associate(x=>s%x,w=>s%w,ws=>s%ws)
44 do ix^
d=ixomin^
d,ixomax^
d
45 w(ix^
d^
d%ixO^s,iw) = w(ixomin^
d-1^
d%ixO^s,iw)
55 do ix^
d=ixomin^
d,ixomax^
d
56 w(ix^
d^
d%ixO^s,iw) = max(w(ixomin^
d-1^
d%ixO^s,iw),zero)
59 do ix^
d=ixomin^
d,ixomax^
d
60 w(ix^
d^
d%ixO^s,iw) = w(ixomin^
d-1^
d%ixO^s,iw)
65 w(ixo^s,iw) = - w(ixo^s,iw)
73 write (
unitterm,*)
"Undefined boundarytype found in bc_phys", &
74 "for variable iw=",iw,
" and side iB=",ib
82 ixosmin^dd=ixomin^dd-
kr(^dd,idir);
87 do ix^
d=ixosmin^
d,ixosmax^
d
88 ws(ix^
d^
d%ixOs^s,idir) = ws(ixosmin^
d-1^
d%ixOs^s,idir)
91 ws(ixos^s,idir) = ws(ixosmin^
d-1:ixosmin^
d-
nghostcells:-1^
d%ixOs^s,idir)
93 ws(ixos^s,idir) =-ws(ixosmin^
d-1:ixosmin^
d-
nghostcells:-1^
d%ixOs^s,idir)
96 write (
unitterm,*)
"Undefined boundarytype found in bc_phys", &
97 "for variable iw=",iw,
" and side iB=",ib
115 ws(ixosmin^
d+ix^
d^
d%ixOs^s,idir)=&
116 (q(hxomax^
d^
d%hxO^s)*s%dvolume(hxomax^
d^
d%hxO^s)&
117 -qp(ixomin^
d+ix^
d^
d%ixO^s)*s%dvolume(ixomin^
d+ix^
d^
d%ixO^s))&
118 /s%surfaceC(ixosmin^
d+ix^
d^
d%ixOs^s,^
d)
124 ws(ixos^s,idir)= ws(ixosmin^
d-2:ixosmin^
d-
nghostcells-1:-1^
d%ixOs^s,idir)
129 ws(ixos^s,idir)=-ws(ixosmin^
d-2:ixosmin^
d-
nghostcells-1:-1^
d%ixOs^s,idir)
149 do ix^
d=ixomin^
d,ixomax^
d
150 w(ix^
d^
d%ixO^s,iw) = w(ixomax^
d+1^
d%ixO^s,iw)
160 do ix^
d=ixomin^
d,ixomax^
d
161 w(ix^
d^
d%ixO^s,iw) = min(w(ixomax^
d+1^
d%ixO^s,iw),zero)
164 do ix^
d=ixomin^
d,ixomax^
d
165 w(ix^
d^
d%ixO^s,iw) = w(ixomax^
d+1^
d%ixO^s,iw)
170 w(ixo^s,iw) = - w(ixo^s,iw)
178 write (
unitterm,*)
"Undefined boundarytype found in bc_phys", &
179 "for variable iw=",iw,
" and side iB=",ib
186 ixosmax^dd=ixomax^dd;
187 ixosmin^dd=ixomin^dd-
kr(^dd,idir);
192 do ix^
d=ixosmin^
d,ixosmax^
d
193 ws(ix^
d^
d%ixOs^s,idir) = ws(ixosmax^
d+1^
d%ixOs^s,idir)
196 ws(ixos^s,idir) = ws(ixosmax^
d+
nghostcells:ixosmax^
d+1:-1^
d%ixOs^s,idir)
198 ws(ixos^s,idir) =-ws(ixosmax^
d+
nghostcells:ixosmax^
d+1:-1^
d%ixOs^s,idir)
201 write (
unitterm,*)
"Undefined boundarytype in bc_phys", &
202 "for variable iw=",iw,
" and side iB=",ib
211 ixos^
l=ixo^
l-
kr(^dd,^
d);
220 ws(ixosmax^
d-ix^
d^
d%ixOs^s,idir)=&
221 -(q(jxomin^
d^
d%jxO^s)*s%dvolume(jxomin^
d^
d%jxO^s)&
222 -qp(ixomax^
d-ix^
d^
d%ixO^s)*s%dvolume(ixomax^
d-ix^
d^
d%ixO^s))&
223 /s%surfaceC(ixosmax^
d-ix^
d^
d%ixOs^s,^
d)
229 ws(ixos^s,idir)= ws(ixosmax^
d+
nghostcells+1:ixosmax^
d+2:-1^
d%ixOs^s,idir)
234 ws(ixos^s,idir)=-ws(ixosmax^
d+
nghostcells+1:ixosmax^
d+2:-1^
d%ixOs^s,idir)
249 call mpistop(
"usr_special_bc not defined")
262 call mpistop(
"usr_special_bc not defined")
274 double precision,
intent(in) :: time
275 integer,
intent(in) :: ixg^
l
277 integer :: iigrid, igrid, ixo^
l
279 ixo^
l=ixg^
l^lsubnghostcells;
282 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
Module to set boundary conditions from user data.
subroutine, public bc_data_set(qt, ixIL, ixOL, iB, w, x)
Set boundary conditions according to user data.
subroutine, public getintbc(time, ixGL)
fill inner boundary values
subroutine, public bc_phys(iside, idims, time, qdt, s, ixGL, ixBL)
fill ghost cells at a physical boundary
subroutine, public get_divb(w, ixIL, ixOL, divb, fourthorder)
Calculate div B within ixO.
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 bc_noinflow
integer, parameter bc_asymm
integer, dimension(3, 3) kr
Kronecker delta tensor.
integer, parameter bc_character
integer, dimension(:, :), allocatable typeboundary
Array indicating the type of boundary condition per variable and per physical boundary.
logical stagger_grid
True for using stagger grid.
integer, parameter bc_data
integer, parameter plevel_
integer, dimension(:), allocatable, parameter d
integer ixm
the mesh range of a physical block without ghost cells
integer, parameter bc_periodic
integer, parameter bc_special
boundary condition types
integer, parameter unitterm
Unit for standard output.
double precision, dimension(:,:), allocatable rnode
Corner coordinates.
integer, parameter bc_cont
integer, parameter bc_icarus
integer nghostcells
Number of ghost cells surrounding a grid.
integer, parameter bc_symm
integer, parameter bc_aperiodic
integer, dimension(:,:), allocatable node
double precision, dimension(ndim) dxlevel
This module defines the procedures of a physics module. It contains function pointers for the various...
procedure(sub_face_to_center), pointer phys_face_to_center
Module with all the methods that users can customize in AMRVAC.
procedure(special_bc), pointer usr_special_bc
procedure(internal_bc), pointer usr_internal_bc