51 subroutine bound(time,nvar,ng,valbig,mitot,mjtot,mptr,aux,naux)
53 use amr_module, only: rnode, node, hxposs, hyposs, cornxlo, cornxhi
54 use amr_module, only: cornylo, cornyhi, ndilo, ndihi, ndjlo, ndjhi
55 use amr_module, only: nestlevel, xlower, xupper, ylower, yupper
56 use amr_module, only: xperdom, yperdom, spheredom
61 integer,
intent(in) :: nvar, ng, mitot, mjtot, mptr, naux
62 real(kind=8),
intent(in) :: time
63 real(kind=8),
intent(in out) :: valbig(nvar,mitot,mjtot)
64 real(kind=8),
intent(in out) :: aux(naux,mitot,mjtot)
67 integer :: ilo, ihi, jlo, jhi, level
68 real(kind=8) :: xleft, xright, ybot, ytop, hx, hy, xl, xr, yb, yt
69 real(kind=8) :: xlowithghost, xhiwithghost, ylowithghost, yhiwithghost
72 xleft = rnode(cornxlo, mptr)
73 xright = rnode(cornxhi, mptr)
74 ybot = rnode(cornylo, mptr)
75 ytop = rnode(cornyhi, mptr)
76 ilo = node(ndilo, mptr)
77 ihi = node(ndihi, mptr)
78 jlo = node(ndjlo, mptr)
79 jhi = node(ndjhi, mptr)
80 level = node(nestlevel, mptr)
84 xlowithghost = xleft - ng*hx
85 xhiwithghost = xright + ng*hx
86 ylowithghost = ybot - ng*hy
87 yhiwithghost = ytop + ng*hy
99 if ((xl < xlower) .and. xperdom)
then
100 call
prefilrecur(level,nvar,valbig,aux,naux,time,mitot,mjtot,1,ng+1, &
101 ilo-ng,ilo-1,jlo,jhi,ilo-ng,ihi+ng,jlo-ng,jhi+ng,patchonly)
103 call
filrecur(level,nvar,valbig,aux,naux,time,mitot,mjtot,1,ng+1,ilo-ng,ilo-1,jlo,jhi,patchonly,mptr)
112 if ((xr .gt. xupper) .and. xperdom)
then
113 call
prefilrecur(level,nvar,valbig,aux,naux,time,mitot,mjtot, &
114 mitot-ng+1,ng+1,ihi+1,ihi+ng,jlo,jhi,ilo-ng,ihi+ng,jlo-ng,jhi+ng,patchonly)
116 call
filrecur(level,nvar,valbig,aux,naux,time,mitot,mjtot, &
117 mitot-ng+1,ng+1,ihi+1,ihi+ng,jlo,jhi,patchonly,mptr)
126 if ( ((yb < ylower) .and. (yperdom .or. spheredom)) .or. &
127 ( ((xl < xlower) .or. (xr > xupper)) .and. xperdom) )
then
128 call
prefilrecur(level,nvar,valbig,aux,naux,time,mitot,mjtot, &
129 1,1,ilo-ng,ihi+ng,jlo-ng,jlo-1, &
130 ilo-ng,ihi+ng,jlo-ng,jhi+ng,patchonly)
132 call
filrecur(level,nvar,valbig,aux,naux,time,mitot,mjtot,1,1,ilo-ng,ihi+ng,jlo-ng,jlo-1,patchonly,mptr)
142 if ( ((yt .gt. yupper) .and. (yperdom .or. spheredom)) .or. &
143 (((xl .lt. xlower) .or. (xr .gt. xupper)) .and. xperdom) )
then
144 call
prefilrecur(level,nvar,valbig,aux,naux,time,mitot,mjtot, &
145 1,mjtot-ng+1,ilo-ng,ihi+ng,jhi+1,jhi+ng,ilo-ng,ihi+ng,jlo-ng,jhi+ng,patchonly)
147 call
filrecur(level,nvar,valbig,aux,naux,time,mitot,mjtot, &
148 1,mjtot-ng+1,ilo-ng,ihi+ng,jhi+1,jhi+ng,patchonly,mptr)
155 call
bc2amr(valbig,aux,mitot,mjtot,nvar,naux,hx,hy,level,time, &
156 xlowithghost,xhiwithghost,ylowithghost,yhiwithghost)
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...
recursive subroutine prefilrecur(level, nvar, valbig, auxbig, naux, time, mitot, mjtot, nrowst, ncolst, ilo, ihi, jlo, jhi, iglo, ighi, jglo, jghi, patchOnly)
recursive subroutine filrecur(level, nvar, valbig, aux, naux, t, mitot, mjtot, nrowst, ncolst, ilo, ihi, jlo, jhi, patchOnly, msrc)
Fill a region (patch) described by:
subroutine bound(time, nvar, ng, valbig, mitot, mjtot, mptr, aux, naux)
This routine sets the boundary values for a given grid at level level.