2D AMRCLAW
bound.f90
Go to the documentation of this file.
1 !
2 ! :::::::::::::: BOUND :::::::::::::::::::::::::::::::::::::::::::
3 ! This routine sets the boundary values for a given grid
4 ! at level level.
5 ! We are setting the values for a strip ng zones wide all
6 ! the way around the border, in 4 rectangular strips.
7 !
8 ! Outputs from this routine:
9 ! The values around the border of the grid are inserted
10 ! directly into the enlarged valbig array.
11 !
12 ! This routine calls the routine filpatch
13 ! which for any block of mesh points on a given level,
14 ! intersects that block with all grids on that level and with
15 ! the physical boundaries, copies the values into the
16 ! appropriate intersecting regions, and interpolates the remaining
17 ! cells from coarser grids as required.
18 ! :::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
19 !
20 ! Below are comments for Doxygen
51 subroutine bound(time,nvar,ng,valbig,mitot,mjtot,mptr,aux,naux)
52 
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
57 
58  implicit none
59 
60  ! Input
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)
65 
66  ! Locals
67  integer :: ilo, ihi, jlo, jhi, level ! rect(4)
68  real(kind=8) :: xleft, xright, ybot, ytop, hx, hy, xl, xr, yb, yt
69  real(kind=8) :: xlowithghost, xhiwithghost, ylowithghost, yhiwithghost
70  logical :: patchonly
71 
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)
81  hx = hxposs(level)
82  hy = hyposs(level)
83 
84  xlowithghost = xleft - ng*hx
85  xhiwithghost = xright + ng*hx
86  ylowithghost = ybot - ng*hy
87  yhiwithghost = ytop + ng*hy
88  ! used in filaptch for bc2amr: for patches it is called. for full grids called from bound below
89  patchonly = .false.
90 
91 
92 
93 
94  ! left boundary
95  xl = xleft - ng*hx
96  xr = xleft
97  yb = ybot
98  yt = ytop
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)
102  else
103  call filrecur(level,nvar,valbig,aux,naux,time,mitot,mjtot,1,ng+1,ilo-ng,ilo-1,jlo,jhi,patchonly,mptr)
104  endif
105 
106  ! right boundary
107  xl = xright
108  xr = xright + ng*hx
109  yb = ybot
110  yt = ytop
111 
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)
115  else
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)
118  endif
119 
120  ! bottom boundary
121  xl = xleft - ng*hx
122  xr = xright + ng*hx
123  yb = ybot - ng*hy
124  yt = ybot
125 
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)
131  else
132  call filrecur(level,nvar,valbig,aux,naux,time,mitot,mjtot,1,1,ilo-ng,ihi+ng,jlo-ng,jlo-1,patchonly,mptr)
133  endif
134 
135 
136  ! top boundary
137  xl = xleft - ng*hx
138  xr = xright + ng*hx
139  yb = ytop
140  yt = ytop + ng*hy
141 
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)
146  else
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)
149  endif
150 
151 
152  ! set all exterior (physical) boundary conditions for this grid at once
153  ! used to be done from filpatch, but now only for recursive calls with new patch
154  ! where the info matches. more efficient to do whole grid at once, and avoid copying
155  call bc2amr(valbig,aux,mitot,mjtot,nvar,naux,hx,hy,level,time, &
156  xlowithghost,xhiwithghost,ylowithghost,yhiwithghost)
157 
158 end subroutine bound
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...
Definition: bc2amr.f:87
recursive subroutine prefilrecur(level, nvar, valbig, auxbig, naux, time, mitot, mjtot, nrowst, ncolst, ilo, ihi, jlo, jhi, iglo, ighi, jglo, jghi, patchOnly)
Definition: prefilp.f90:18
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:
Definition: filpatch.f90:71
subroutine bound(time, nvar, ng, valbig, mitot, mjtot, mptr, aux, naux)
This routine sets the boundary values for a given grid at level level.
Definition: bound.f90:51