85 c ------------------------------------------------------------------
87 subroutine bc2amr(val,aux,nrow,ncol,meqn,naux,
88 1 hx, hy, level, time,
89 2 xlo_patch, xhi_patch, ylo_patch,yhi_patch)
94 use amr_module, only: mthbc,xlower,ylower,xupper,yupper
95 use amr_module, only: xperdom,yperdom,spheredom
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
109 if (xperdom .and. (yperdom .or. spheredom)) go to 499
112 c-------------------------------------------------------
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.
122 c
# number of grid cells from this patch lying outside physical domain:
123 nxl = (xlower+hxmarg-xlo_patch)/hx
125 go to(100,110,120,130) mthbc(1)+1
128 c
# user-specified boundary conditions go here in place of error output
130 &
'*** ERROR *** mthbc(1)=0 and no BCs specified in bc2amr'
135 c
# zero-order extrapolation:
139 val(m,i,j) = val(m,nxl+1,j)
144 c
# periodic: handled elsewhere in amr
148 c
# solid wall (assumes 2'nd component is velocity or momentum in x):
152 val(m,i,j) = val(m,2*nxl+1-i,j)
154 c
# negate the normal velocity:
157 val(2,i,j) = -val(2,i,j)
163 c-------------------------------------------------------
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.
173 c
# number of grid cells lying outside physical domain:
174 nxr = (xhi_patch - xupper + hxmarg)/hx
175 ibeg = max0(nrow-nxr+1, 1)
177 go to(200,210,220,230) mthbc(2)+1
180 c
# user-specified boundary conditions go here in place of error output
182 &
'*** ERROR *** mthbc(2)=0 and no BCs specified in bc2amr'
187 c
# zero-order extrapolation:
191 val(m,i,j) = val(m,ibeg-1,j)
196 c
# periodic: handled elsewhere in amr
200 c
# solid wall (assumes 2'nd component is velocity or momentum in x):
204 val(m,i,j) = val(m,2*ibeg-1-i,j)
206 c
# negate the normal velocity:
209 val(2,i,j) = -val(2,i,j)
215 c-------------------------------------------------------
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.
225 c
# number of grid cells lying outside physical domain:
226 nyb = (ylower+hymarg-ylo_patch)/hy
228 go to(300,310,320,330) mthbc(3)+1
231 c
# user-specified boundary conditions go here in place of error output
233 &
'*** ERROR *** mthbc(3)=0 and no BCs specified in bc2amr'
238 c
# zero-order extrapolation:
242 val(m,i,j) = val(m,i,nyb+1)
247 c
# periodic: handled elsewhere in amr
251 c
# solid wall (assumes 3'rd component is velocity or momentum in y):
255 val(m,i,j) = val(m,i,2*nyb+1-j)
257 c
# negate the normal velocity:
260 val(3,i,j) = -val(3,i,j)
266 c-------------------------------------------------------
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.
276 c
# number of grid cells lying outside physical domain:
277 nyt = (yhi_patch - yupper + hymarg)/hy
278 jbeg = max0(ncol-nyt+1, 1)
280 go to(400,410,420,430) mthbc(4)+1
283 c
# user-specified boundary conditions go here in place of error output
285 &
'*** ERROR *** mthbc(4)=0 and no BCs specified in bc2amr'
290 c
# zero-order extrapolation:
294 val(m,i,j) = val(m,i,jbeg-1)
299 c
# periodic: handled elsewhere in amr
303 c
# solid wall (assumes 3'rd component is velocity or momentum in y):
307 val(m,i,j) = val(m,i,2*jbeg-1-j)
309 c
# negate the normal velocity:
312 val(3,i,j) = -val(3,i,j)
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...