MPI-AMRVAC 3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
Loading...
Searching...
No Matches
mod_boundary_conditions.t
Go to the documentation of this file.
2
3 implicit none
4 private
5
6 public :: bc_phys
7 public :: getintbc
8
9contains
10
11 !> fill ghost cells at a physical boundary
12 subroutine bc_phys(iside,idims,time,qdt,s,ixG^L,ixB^L,nwhead,nwtail)
14 use mod_bc_data, only: bc_data_set
16 use mod_physics
18
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
22
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
26
27 associate(x=>s%x,w=>s%w,ws=>s%ws)
28 select case (idims)
29 {case (^d)
30 if (iside==2) then
31 ! maximal boundary
32 ib=2*^d
33 ixomin^dd=ixbmax^d+1-nghostcells^d%ixOmin^dd=ixbmin^dd;
34 ixomax^dd=ixbmax^dd;
35 ! cont/symm/asymm types
36 do iw=nwhead,nwtail
37 select case (typeboundary(iw,ib))
38 case (bc_special)
39 ! skip it here, do AFTER all normal type boundaries are set
40 case (bc_cont)
41 do ix^d=ixomin^d,ixomax^d
42 w(ix^d^d%ixO^s,iw) = w(ixomin^d-1^d%ixO^s,iw)
43 end do
44 case (bc_symm)
45 w(ixo^s,iw) = w(ixomin^d-1:ixomin^d-nghostcells:-1^d%ixO^s,iw)
46 case (bc_asymm)
47 w(ixo^s,iw) =-w(ixomin^d-1:ixomin^d-nghostcells:-1^d%ixO^s,iw)
48 case (bc_periodic)
49 ! skip it here, periodic bc info should come from neighbors
50 case(bc_noinflow)
51 if (iw==1+^d)then
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)
54 end do
55 else
56 do ix^d=ixomin^d,ixomax^d
57 w(ix^d^d%ixO^s,iw) = w(ixomin^d-1^d%ixO^s,iw)
58 end do
59 end if
60 case (bc_aperiodic)
61 !this just multiplies the variables with (-), they have been set from neighbors just like periodic.
62 w(ixo^s,iw) = - w(ixo^s,iw)
63 case (bc_data)
64 ! skip it here, do AFTER all normal type boundaries are set
65 case (bc_icarus)
66 ! skip it here, do AFTER all normal type boundaries are set
67 case (bc_character)
68 ! skip it here, do AFTER all normal type boundaries are set
69 case default
70 write (unitterm,*) "Undefined boundarytype found in bc_phys", &
71 "for variable iw=",iw," and side iB=",ib
72 end select
73 end do
74 if(stagger_grid) then
75 do idir=1,nws
76 ! At this stage, extrapolation is applied only to the tangential components
77 if(idir==^d) cycle
78 ixosmax^dd=ixomax^dd;
79 ixosmin^dd=ixomin^dd-kr(^dd,idir);
80 select case(typeboundary(iw_mag(idir),ib))
81 case (bc_special)
82 ! skip it here, do AFTER all normal type boundaries are set
83 case (bc_cont)
84 do ix^d=ixosmin^d,ixosmax^d
85 ws(ix^d^d%ixOs^s,idir) = ws(ixosmin^d-1^d%ixOs^s,idir)
86 end do
87 case (bc_symm)
88 ws(ixos^s,idir) = ws(ixosmin^d-1:ixosmin^d-nghostcells:-1^d%ixOs^s,idir)
89 case (bc_asymm)
90 ws(ixos^s,idir) =-ws(ixosmin^d-1:ixosmin^d-nghostcells:-1^d%ixOs^s,idir)
91 case (bc_periodic)
92 case default
93 write (unitterm,*) "Undefined boundarytype found in bc_phys", &
94 "for variable iw=",iw," and side iB=",ib
95 end select
96 end do
97 ! Now that the tangential components are set,
98 ! fill the normal components using a prescription for the divergence.
99 ! This prescription is given by the typeB for the normal component.
100 do idir=1,nws
101 ! Consider only normal direction
102 if (idir/=^d) cycle
103 ixos^l=ixo^l;
104 select case(typeboundary(iw_mag(idir),ib))
105 case(bc_cont)
106 hxo^l=ixo^l-nghostcells*kr(^dd,^d);
107 ! Calculate divergence and partial divergence
108 call get_divb(s%w,ixg^l,hxo^l,q)
109 ws(ixos^s,idir)=zero
110 do ix^d=0,nghostcells-1
111 call get_divb(s%w,ixg^l,ixo^l,qp)
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)
116 end do
117 ! Fill cell averages
118 call phys_face_to_center(ixo^l,s)
119 case(bc_symm)
120 ! (a)symmetric normal B ensures symmetric divB
121 ws(ixos^s,idir)= ws(ixosmin^d-2:ixosmin^d-nghostcells-1:-1^d%ixOs^s,idir)
122 ! Fill cell averages
123 call phys_face_to_center(ixo^l,s)
124 case(bc_asymm)
125 ! (a)symmetric normal B ensures symmetric divB
126 ws(ixos^s,idir)=-ws(ixosmin^d-2:ixosmin^d-nghostcells-1:-1^d%ixOs^s,idir)
127 ! Fill cell averages
128 call phys_face_to_center(ixo^l,s)
129 case(bc_periodic)
130 ! Fill cell averages
131 call phys_face_to_center(ixo^l,s)
132 end select
133 end do
134 end if
135 else
136 ! minimal boundary
137 ib=2*^d-1
138 ixomin^dd=ixbmin^dd;
139 ixomax^dd=ixbmin^d-1+nghostcells^d%ixOmax^dd=ixbmax^dd;
140 ! cont/symm/asymm types
141 do iw=nwhead,nwtail
142 select case (typeboundary(iw,ib))
143 case (bc_special)
144 ! skip it here, do AFTER all normal type boundaries are set
145 case (bc_cont)
146 do ix^d=ixomin^d,ixomax^d
147 w(ix^d^d%ixO^s,iw) = w(ixomax^d+1^d%ixO^s,iw)
148 end do
149 case (bc_symm)
150 w(ixo^s,iw) = w(ixomax^d+nghostcells:ixomax^d+1:-1^d%ixO^s,iw)
151 case (bc_asymm)
152 w(ixo^s,iw) =-w(ixomax^d+nghostcells:ixomax^d+1:-1^d%ixO^s,iw)
153 case (bc_periodic)
154 ! skip it here, periodic bc info should come from neighbors
155 case(bc_noinflow)
156 if (iw==1+^d)then
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)
159 end do
160 else
161 do ix^d=ixomin^d,ixomax^d
162 w(ix^d^d%ixO^s,iw) = w(ixomax^d+1^d%ixO^s,iw)
163 end do
164 end if
165 case (bc_aperiodic)
166 !this just multiplies the variables with (-), they have been set from neighbors just like periodic.
167 w(ixo^s,iw) = - w(ixo^s,iw)
168 case (bc_data)
169 ! skip it here, do AFTER all normal type boundaries are set
170 case (bc_icarus)
171 ! skip it here, do AFTER all normal type boundaries are set
172 case (bc_character)
173 ! skip it here, do AFTER all normal type boundaries are set
174 case default
175 write (unitterm,*) "Undefined boundarytype found in bc_phys", &
176 "for variable iw=",iw," and side iB=",ib
177 end select
178 end do
179 if(stagger_grid) then
180 do idir=1,nws
181 ! At this stage, extrapolation is applied only to the tangential components
182 if(idir==^d) cycle
183 ixosmax^dd=ixomax^dd;
184 ixosmin^dd=ixomin^dd-kr(^dd,idir);
185 select case(typeboundary(iw_mag(idir),ib))
186 case (bc_special)
187 ! skip it here, periodic bc info should come from neighbors
188 case (bc_cont)
189 do ix^d=ixosmin^d,ixosmax^d
190 ws(ix^d^d%ixOs^s,idir) = ws(ixosmax^d+1^d%ixOs^s,idir)
191 end do
192 case (bc_symm)
193 ws(ixos^s,idir) = ws(ixosmax^d+nghostcells:ixosmax^d+1:-1^d%ixOs^s,idir)
194 case (bc_asymm)
195 ws(ixos^s,idir) =-ws(ixosmax^d+nghostcells:ixosmax^d+1:-1^d%ixOs^s,idir)
196 case (bc_periodic)
197 case default
198 write (unitterm,*) "Undefined boundarytype in bc_phys", &
199 "for variable iw=",iw," and side iB=",ib
200 end select
201 end do
202 ! Now that the tangential components are set,
203 ! fill the normal components using a prescription for the divergence.
204 ! This prescription is given by the typeB for the normal component.
205 do idir=1,nws
206 ! Consider only normal direction
207 if (idir/=^d) cycle
208 ixos^l=ixo^l-kr(^dd,^d);
209 select case(typeboundary(iw_mag(idir),ib))
210 case(bc_cont)
211 jxo^l=ixo^l+nghostcells*kr(^dd,^d);
212 ! Calculate divergence and partial divergence
213 call get_divb(s%w,ixg^l,jxo^l,q)
214 ws(ixos^s,idir)=zero
215 do ix^d=0,nghostcells-1
216 call get_divb(s%w,ixg^l,ixo^l,qp)
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)
221 end do
222 ! Fill cell averages
223 call phys_face_to_center(ixo^l,s)
224 case(bc_symm)
225 ! (a)symmetric normal B ensures symmetric divB
226 ws(ixos^s,idir)= ws(ixosmax^d+nghostcells+1:ixosmax^d+2:-1^d%ixOs^s,idir)
227 ! Fill cell averages
228 call phys_face_to_center(ixo^l,s)
229 case(bc_asymm)
230 ! (a)symmetric normal B ensures symmetric divB
231 ws(ixos^s,idir)=-ws(ixosmax^d+nghostcells+1:ixosmax^d+2:-1^d%ixOs^s,idir)
232 ! Fill cell averages
233 call phys_face_to_center(ixo^l,s)
234 case(bc_periodic)
235 ! Fill cell averages
236 call phys_face_to_center(ixo^l,s)
237 end select
238 end do
239 end if
240 end if \}
241 end select
242
243 ! do user defined special boundary conditions
244 if (any(typeboundary(1:nwflux+nwaux,ib)==bc_special)) then
245 if (.not. associated(usr_special_bc)) &
246 call mpistop("usr_special_bc not defined")
247 call usr_special_bc(time,ixg^l,ixo^l,ib,w,x)
248 end if
249
250 ! fill boundary conditions from external data vtk files
251 if (any(typeboundary(1:nwflux+nwaux,ib)==bc_data)) then
252 call bc_data_set(time,ixg^l,ixo^l,ib,w,x)
253 end if
254
255 ! fill boundary conditions from external data vtk files and do user defined special boundary conditions
256 if (any(typeboundary(1:nwflux+nwaux,ib)==bc_icarus)) then
257 call bc_data_set(time,ixg^l,ixo^l,ib,w,x)
258 if (.not. associated(usr_special_bc)) &
259 call mpistop("usr_special_bc not defined")
260 call usr_special_bc(time,ixg^l,ixo^l,ib,w,x)
261 end if
262
263 end associate
264 end subroutine bc_phys
265
266 !> fill inner boundary values
267 subroutine getintbc(time,ixG^L)
270
271 double precision, intent(in) :: time
272 integer, intent(in) :: ixg^l
273
274 integer :: iigrid, igrid, ixo^l
275
276 ixo^l=ixg^l^lsubnghostcells;
277
278 !$OMP PARALLEL DO SCHEDULE(dynamic) PRIVATE(igrid)
279 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
280 block=>ps(igrid)
281 ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
282
283 if (associated(usr_internal_bc)) then
284 call usr_internal_bc(node(plevel_,igrid),time,ixg^l,ixo^l,ps(igrid)%w,ps(igrid)%x)
285 end if
286 end do
287 !$OMP END PARALLEL DO
288
289 end subroutine getintbc
290
Module to set boundary conditions from user data.
Definition mod_bc_data.t:2
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, nth_in)
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, 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.
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 nghostcells
Number of ghost cells surrounding a grid.
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...
Definition mod_physics.t:4
procedure(sub_face_to_center), pointer phys_face_to_center
Definition mod_physics.t:81
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