12 subroutine bc_phys(iside,idims,time,qdt,s,ixG^L,ixB^L,nwhead,nwtail)
19 integer,
intent(in) :: iside, idims, ixg^
l,ixb^
l,nwhead,nwtail
20 double precision,
intent(in) :: time,qdt
21 type(state),
intent(inout) :: s
23 double precision :: q(ixg^s),qp(ixg^s)
24 integer :: idir, ixos^
l,hxo^
l,jxo^
l
25 integer :: iw, ib, ix^
d, ixo^
l
27 associate(x=>s%x,w=>s%w,ws=>s%ws)
41 do ix^
d=ixomin^
d,ixomax^
d
42 w(ix^
d^
d%ixO^s,iw) = w(ixomin^
d-1^
d%ixO^s,iw)
52 do ix^
d=ixomin^
d,ixomax^
d
53 w(ix^
d^
d%ixO^s,iw) = max(w(ixomin^
d-1^
d%ixO^s,iw),zero)
56 do ix^
d=ixomin^
d,ixomax^
d
57 w(ix^
d^
d%ixO^s,iw) = w(ixomin^
d-1^
d%ixO^s,iw)
62 w(ixo^s,iw) = - w(ixo^s,iw)
70 write (
unitterm,*)
"Undefined boundarytype found in bc_phys", &
71 "for variable iw=",iw,
" and side iB=",ib
79 ixosmin^dd=ixomin^dd-
kr(^dd,idir);
84 do ix^
d=ixosmin^
d,ixosmax^
d
85 ws(ix^
d^
d%ixOs^s,idir) = ws(ixosmin^
d-1^
d%ixOs^s,idir)
88 ws(ixos^s,idir) = ws(ixosmin^
d-1:ixosmin^
d-
nghostcells:-1^
d%ixOs^s,idir)
90 ws(ixos^s,idir) =-ws(ixosmin^
d-1:ixosmin^
d-
nghostcells:-1^
d%ixOs^s,idir)
93 write (
unitterm,*)
"Undefined boundarytype found in bc_phys", &
94 "for variable iw=",iw,
" and side iB=",ib
112 ws(ixosmin^
d+ix^
d^
d%ixOs^s,idir)=&
113 (q(hxomax^
d^
d%hxO^s)*s%dvolume(hxomax^
d^
d%hxO^s)&
114 -qp(ixomin^
d+ix^
d^
d%ixO^s)*s%dvolume(ixomin^
d+ix^
d^
d%ixO^s))&
115 /s%surfaceC(ixosmin^
d+ix^
d^
d%ixOs^s,^
d)
121 ws(ixos^s,idir)= ws(ixosmin^
d-2:ixosmin^
d-
nghostcells-1:-1^
d%ixOs^s,idir)
126 ws(ixos^s,idir)=-ws(ixosmin^
d-2:ixosmin^
d-
nghostcells-1:-1^
d%ixOs^s,idir)
146 do ix^
d=ixomin^
d,ixomax^
d
147 w(ix^
d^
d%ixO^s,iw) = w(ixomax^
d+1^
d%ixO^s,iw)
157 do ix^
d=ixomin^
d,ixomax^
d
158 w(ix^
d^
d%ixO^s,iw) = min(w(ixomax^
d+1^
d%ixO^s,iw),zero)
161 do ix^
d=ixomin^
d,ixomax^
d
162 w(ix^
d^
d%ixO^s,iw) = w(ixomax^
d+1^
d%ixO^s,iw)
167 w(ixo^s,iw) = - w(ixo^s,iw)
175 write (
unitterm,*)
"Undefined boundarytype found in bc_phys", &
176 "for variable iw=",iw,
" and side iB=",ib
183 ixosmax^dd=ixomax^dd;
184 ixosmin^dd=ixomin^dd-
kr(^dd,idir);
189 do ix^
d=ixosmin^
d,ixosmax^
d
190 ws(ix^
d^
d%ixOs^s,idir) = ws(ixosmax^
d+1^
d%ixOs^s,idir)
193 ws(ixos^s,idir) = ws(ixosmax^
d+
nghostcells:ixosmax^
d+1:-1^
d%ixOs^s,idir)
195 ws(ixos^s,idir) =-ws(ixosmax^
d+
nghostcells:ixosmax^
d+1:-1^
d%ixOs^s,idir)
198 write (
unitterm,*)
"Undefined boundarytype in bc_phys", &
199 "for variable iw=",iw,
" and side iB=",ib
208 ixos^
l=ixo^
l-
kr(^dd,^
d);
217 ws(ixosmax^
d-ix^
d^
d%ixOs^s,idir)=&
218 -(q(jxomin^
d^
d%jxO^s)*s%dvolume(jxomin^
d^
d%jxO^s)&
219 -qp(ixomax^
d-ix^
d^
d%ixO^s)*s%dvolume(ixomax^
d-ix^
d^
d%ixO^s))&
220 /s%surfaceC(ixosmax^
d-ix^
d^
d%ixOs^s,^
d)
226 ws(ixos^s,idir)= ws(ixosmax^
d+
nghostcells+1:ixosmax^
d+2:-1^
d%ixOs^s,idir)
231 ws(ixos^s,idir)=-ws(ixosmax^
d+
nghostcells+1:ixosmax^
d+2:-1^
d%ixOs^s,idir)
246 call mpistop(
"usr_special_bc not defined")
259 call mpistop(
"usr_special_bc not defined")
271 double precision,
intent(in) :: time
272 integer,
intent(in) :: ixg^
l
274 integer :: iigrid, igrid, ixo^
l
276 ixo^
l=ixg^
l^lsubnghostcells;
279 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, nwhead, nwtail)
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_
double precision, dimension(:), allocatable, parameter d
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
double precision, dimension(^nd) dxlevel
store unstretched cell size of current level
integer, parameter bc_aperiodic
integer, dimension(:,:), allocatable node
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