51 subroutine bound(time,nvar,ng,valbig,mitot,mjtot,mptr,aux,naux)
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)
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)
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...
integer, parameter cornxlo
x-coordinate of the left border of this grid
recursive subroutine prefilrecur(level, nvar, valbig, auxbig, naux, time, mitot, mjtot, nrowst, ncolst, ilo, ihi, jlo, jhi, iglo, ighi, jglo, jghi, patchOnly)
For periodic boundary conditions more work needed to fill the piece of the boundary.
real(kind=8), dimension(maxlv) hyposs
real(kind=8), dimension(maxlv) hxposs
integer, parameter ndihi
global i index of right border of this grid
integer, dimension(nsize, maxgr) node
integer, parameter nestlevel
AMR level of the grid.
real(kind=8), dimension(rsize, maxgr) rnode
integer, parameter ndilo
global i index of left border of this grid
integer, parameter ndjlo
global j index of lower border of this grid
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:
integer, parameter ndjhi
global j index of upper border of this grid
subroutine bound(time, nvar, ng, valbig, mitot, mjtot, mptr, aux, naux)
This routine sets the boundary values for a given grid at level level.
integer, parameter cornylo
y-coordinate of the lower border of this grid
integer, parameter cornxhi
x-coordinate of the right border of this grid
The module contains the definition of a "node descriptor" as well as other global variables used duri...
integer, parameter cornyhi
y-coordinate of the upper border of this grid