2D AMRCLAW
Functions/Subroutines
bc2amr.f File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine bc2amr (val, aux, nrow, ncol, meqn, naux,
 Take a grid patch with mesh widths hx,hy, of dimensions nrow by ncol, and set the values of any piece of of the patch which extends outside the physical domain using the boundary conditions. More...
 

Function/Subroutine Documentation

◆ bc2amr()

subroutine bc2amr ( real*8, dimension(meqn,nrow,ncol)  val,
real*8, dimension(naux,nrow,ncol)  aux,
integer  nrow,
integer  ncol,
integer  meqn,
integer  naux 
)

Take a grid patch with mesh widths hx,hy, of dimensions nrow by ncol, and set the values of any piece of of the patch which extends outside the physical domain using the boundary conditions.

Standard boundary condition choices for amr2ez in clawpack

At each boundary k = 1 (left), 2 (right), 3 (bottom), 4 (top):

mthbc(k) =

  • 0 for user-supplied BC's (must be inserted!)
  • 1 for zero-order extrapolation
  • 2 for periodic boundary conditions
  • 3 for solid walls, assuming this can be implemented by reflecting the data about the boundary and then negating the 2'nd (for k=1,2) or 3'rd (for k=3,4) component of q.
  • 5 for sphere bcs (left half maps to right half of same side, and vice versa), as if domain folded in half

The corners of the grid patch are at (xlo_patch,ylo_patch) – lower left corner (xhi_patch,yhi_patch) – upper right corner

The physical domain itself is a rectangle bounded by (xlower,ylower) – lower left corner (xupper,yupper) – upper right corner Any cells that lie outside the physical domain are ghost cells whose values should be set in this routine. This is tested for by comparing xlo_patch with xlower to see if values need to be set at the left and similarly at the other boundaries. Patches are guaranteed to have at least 1 row of cells filled with interior values so it is possible to extrapolate. Fix trimbd() if you want more than 1 row pre-set.

Make sure the order the boundaries are specified is correct so that diagonal corner cells are also properly taken care of.

Periodic boundaries are set before calling this routine, so if the domain is periodic in one direction only you can safely extrapolate in the other direction.

Don't overwrite ghost cells in periodic directions!

Parameters
valdata array for solution $q $ (cover the whole grid msrc)
auxdata array for auxiliary variables
nrownumber of cells in i direction on this grid
ncolnumber of cells in j direction on this grid
meqnnumber of equations for the system
nauxnumber of auxiliary variables
hxspacing (mesh size) in i direction
hyspacing (mesh size) in j direction
levelAMR level of this grid
timesetting ghost cell values at time time
xlo_patchleft bound of the input grid
xhi_patchright bound of the input grid
ylo_patchlower bound of the input grid
yhi_patchupper bound of the input grid

Definition at line 88 of file bc2amr.f.

References amr_module::mthbc, amr_module::spheredom, amr_module::xlower, amr_module::xperdom, amr_module::xupper, amr_module::ylower, amr_module::yperdom, and amr_module::yupper.

Referenced by bound(), filrecur(), filval(), and saveqc().

88  1 hx, hy, level, time,
89  2 xlo_patch, xhi_patch, ylo_patch,yhi_patch)
90 
91 c
92 c
93 
96 
97  implicit none
98 
99  real*8 val(meqn,nrow,ncol), aux(naux,nrow,ncol)
100  integer nrow,ncol,meqn,naux,level
101  real*8 hx,hy,time, hxmarg, hymarg
102  real*8 xlo_patch,xhi_patch,ylo_patch,yhi_patch
103  integer nxl,nxr,ibeg,nyb,nyt,jbeg,i,j,m
104 
105 
106  hxmarg = hx*.01
107  hymarg = hy*.01
108 
109  if (xperdom .and. (yperdom .or. spheredom)) go to 499
110 c
111 c
112 c-------------------------------------------------------
113 c # left boundary:
114 c-------------------------------------------------------
115  if (xlo_patch .ge. xlower-hxmarg) then
116 c # not a physical boundary -- no cells at this edge lies
117 c # outside the physical bndry.
118 c # values are set elsewhere in amr code.
119  go to 199
120  endif
121 c
122 c # number of grid cells from this patch lying outside physical domain:
123  nxl = (xlower+hxmarg-xlo_patch)/hx
124 c
125  go to (100,110,120,130) mthbc(1)+1
126 c
127  100 continue
128 c # user-specified boundary conditions go here in place of error output
129  write(6,*)
130  & '*** ERROR *** mthbc(1)=0 and no BCs specified in bc2amr'
131  stop
132  go to 199
133 c
134  110 continue
135 c # zero-order extrapolation:
136  do 115 j = 1,ncol
137  do 115 i=1,nxl
138  do 115 m=1,meqn
139  val(m,i,j) = val(m,nxl+1,j)
140  115 continue
141  go to 199
142 
143  120 continue
144 c # periodic: handled elsewhere in amr
145  go to 199
146 
147  130 continue
148 c # solid wall (assumes 2'nd component is velocity or momentum in x):
149  do 135 j = 1,ncol
150  do 135 i=1,nxl
151  do 135 m=1,meqn
152  val(m,i,j) = val(m,2*nxl+1-i,j)
153  135 continue
154 c # negate the normal velocity:
155  do 136 j = 1,ncol
156  do 136 i=1,nxl
157  val(2,i,j) = -val(2,i,j)
158  136 continue
159  go to 199
160 
161  199 continue
162 c
163 c-------------------------------------------------------
164 c # right boundary:
165 c-------------------------------------------------------
166  if (xhi_patch .le. xupper+hxmarg) then
167 c # not a physical boundary -- no cells at this edge lies
168 c # outside the physical bndry.
169 c # values are set elsewhere in amr code.
170  go to 299
171  endif
172 c
173 c # number of grid cells lying outside physical domain:
174  nxr = (xhi_patch - xupper + hxmarg)/hx
175  ibeg = max0(nrow-nxr+1, 1)
176 c
177  go to (200,210,220,230) mthbc(2)+1
178 c
179  200 continue
180 c # user-specified boundary conditions go here in place of error output
181  write(6,*)
182  & '*** ERROR *** mthbc(2)=0 and no BCs specified in bc2amr'
183  stop
184  go to 299
185 
186  210 continue
187 c # zero-order extrapolation:
188  do 215 j = 1,ncol
189  do 215 i=ibeg,nrow
190  do 215 m=1,meqn
191  val(m,i,j) = val(m,ibeg-1,j)
192  215 continue
193  go to 299
194 
195  220 continue
196 c # periodic: handled elsewhere in amr
197  go to 299
198 
199  230 continue
200 c # solid wall (assumes 2'nd component is velocity or momentum in x):
201  do 235 j = 1,ncol
202  do 235 i=ibeg,nrow
203  do 235 m=1,meqn
204  val(m,i,j) = val(m,2*ibeg-1-i,j)
205  235 continue
206 c # negate the normal velocity:
207  do 236 j = 1,ncol
208  do 236 i=ibeg,nrow
209  val(2,i,j) = -val(2,i,j)
210  236 continue
211  go to 299
212 
213  299 continue
214 c
215 c-------------------------------------------------------
216 c # bottom boundary:
217 c-------------------------------------------------------
218  if (ylo_patch .ge. ylower-hymarg) then
219 c # not a physical boundary -- no cells at this edge lies
220 c # outside the physical bndry.
221 c # values are set elsewhere in amr code.
222  go to 399
223  endif
224 c
225 c # number of grid cells lying outside physical domain:
226  nyb = (ylower+hymarg-ylo_patch)/hy
227 c
228  go to (300,310,320,330) mthbc(3)+1
229 c
230  300 continue
231 c # user-specified boundary conditions go here in place of error output
232  write(6,*)
233  & '*** ERROR *** mthbc(3)=0 and no BCs specified in bc2amr'
234  stop
235  go to 399
236 c
237  310 continue
238 c # zero-order extrapolation:
239  do 315 j=1,nyb
240  do 315 i=1,nrow
241  do 315 m=1,meqn
242  val(m,i,j) = val(m,i,nyb+1)
243  315 continue
244  go to 399
245 
246  320 continue
247 c # periodic: handled elsewhere in amr
248  go to 399
249 
250  330 continue
251 c # solid wall (assumes 3'rd component is velocity or momentum in y):
252  do 335 j=1,nyb
253  do 335 i=1,nrow
254  do 335 m=1,meqn
255  val(m,i,j) = val(m,i,2*nyb+1-j)
256  335 continue
257 c # negate the normal velocity:
258  do 336 j=1,nyb
259  do 336 i=1,nrow
260  val(3,i,j) = -val(3,i,j)
261  336 continue
262  go to 399
263 
264  399 continue
265 c
266 c-------------------------------------------------------
267 c # top boundary:
268 c-------------------------------------------------------
269  if (yhi_patch .le. yupper+hymarg) then
270 c # not a physical boundary -- no cells at this edge lies
271 c # outside the physical bndry.
272 c # values are set elsewhere in amr code.
273  go to 499
274  endif
275 c
276 c # number of grid cells lying outside physical domain:
277  nyt = (yhi_patch - yupper + hymarg)/hy
278  jbeg = max0(ncol-nyt+1, 1)
279 c
280  go to (400,410,420,430) mthbc(4)+1
281 c
282  400 continue
283 c # user-specified boundary conditions go here in place of error output
284  write(6,*)
285  & '*** ERROR *** mthbc(4)=0 and no BCs specified in bc2amr'
286  stop
287  go to 499
288 
289  410 continue
290 c # zero-order extrapolation:
291  do 415 j=jbeg,ncol
292  do 415 i=1,nrow
293  do 415 m=1,meqn
294  val(m,i,j) = val(m,i,jbeg-1)
295  415 continue
296  go to 499
297 
298  420 continue
299 c # periodic: handled elsewhere in amr
300  go to 499
301 
302  430 continue
303 c # solid wall (assumes 3'rd component is velocity or momentum in y):
304  do 435 j=jbeg,ncol
305  do 435 i=1,nrow
306  do 435 m=1,meqn
307  val(m,i,j) = val(m,i,2*jbeg-1-j)
308  435 continue
309 c # negate the normal velocity:
310  do 436 j=jbeg,ncol
311  do 436 i=1,nrow
312  val(3,i,j) = -val(3,i,j)
313  436 continue
314  go to 499
315 
316  499 continue
317 
318  return
real(kind=8) xupper
Definition: amr_module.f90:231
real(kind=8) xlower
Definition: amr_module.f90:231
logical yperdom
Definition: amr_module.f90:230
real(kind=8) yupper
Definition: amr_module.f90:231
real(kind=8) ylower
Definition: amr_module.f90:231
logical spheredom
Definition: amr_module.f90:230
logical xperdom
Definition: amr_module.f90:230
integer, dimension(4) mthbc
Definition: amr_module.f90:232
The module contains the definition of a "node descriptor" as well as other global variables used duri...
Definition: amr_module.f90:21
Here is the caller graph for this function: