2D AMRCLAW
Functions/Subroutines
bound.f90 File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine bound (time, nvar, ng, valbig, mitot, mjtot, mptr, aux, naux)
 This routine sets the boundary values for a given grid at level level. More...
 

Function/Subroutine Documentation

◆ bound()

subroutine bound ( real(kind=8), intent(in)  time,
integer, intent(in)  nvar,
integer, intent(in)  ng,
real(kind=8), dimension(nvar,mitot,mjtot), intent(inout)  valbig,
integer, intent(in)  mitot,
integer, intent(in)  mjtot,
integer, intent(in)  mptr,
real(kind=8), dimension(naux,mitot,mjtot), intent(inout)  aux,
integer, intent(in)  naux 
)

This routine sets the boundary values for a given grid at level level.

It sets the values for a strip ng zones wide all the way around the border, in 4 rectangular strips.

This routine calls the routine filpatch() which for any block of mesh points on a given level, intersects that block with all grids on that level, copies the values into the appropriate intersecting regions, and interpolates the remaining cells from coarser grids as required. After calling filpatch for four regions around the border, call bc2amr() to set any regions that are outside physical boundary.

Output:

The values around the border of the grid are inserted directly into the enlarged valbig array.

Parameters
timethe grid mptr needs ghost cell values at time t
nvarnumber of equations for the system
ngwidth of ghost regions on each side
valbigdata array for solution $q $ (cover the whole grid msrc)
mitotnumber of cells in i direction on grid mptr
mjtotnumber of cells in j direction on grid mptr
mptrindex of the grid to be processed
auxdata array for auxiliary variables
nauxnumber of auxiliary variables

Definition at line 52 of file bound.f90.

References bc2amr(), amr_module::cornxhi, amr_module::cornxlo, amr_module::cornyhi, amr_module::cornylo, filrecur(), amr_module::hxposs, amr_module::hyposs, amr_module::ndihi, amr_module::ndilo, amr_module::ndjhi, amr_module::ndjlo, amr_module::nestlevel, amr_module::node, prefilrecur(), amr_module::rnode, amr_module::spheredom, amr_module::xlower, amr_module::xperdom, amr_module::xupper, amr_module::ylower, amr_module::yperdom, and amr_module::yupper.

Referenced by advanc(), flagger(), prepbigstep(), and spest2().

52 
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 
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:88
integer, parameter cornxlo
x-coordinate of the left border of this grid
Definition: amr_module.f90:143
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.
Definition: prefilp.f90:20
real(kind=8), dimension(maxlv) hyposs
Definition: amr_module.f90:193
real(kind=8), dimension(maxlv) hxposs
Definition: amr_module.f90:193
integer, parameter ndihi
global i index of right border of this grid
Definition: amr_module.f90:111
integer, dimension(nsize, maxgr) node
Definition: amr_module.f90:198
integer, parameter nestlevel
AMR level of the grid.
Definition: amr_module.f90:44
real(kind=8) xupper
Definition: amr_module.f90:231
real(kind=8), dimension(rsize, maxgr) rnode
Definition: amr_module.f90:193
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
integer, parameter ndilo
global i index of left border of this grid
Definition: amr_module.f90:108
integer, parameter ndjlo
global j index of lower border of this grid
Definition: amr_module.f90:114
logical spheredom
Definition: amr_module.f90:230
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:73
integer, parameter ndjhi
global j index of upper border of this grid
Definition: amr_module.f90:117
logical xperdom
Definition: amr_module.f90:230
integer, parameter cornylo
y-coordinate of the lower border of this grid
Definition: amr_module.f90:145
integer, parameter cornxhi
x-coordinate of the right border of this grid
Definition: amr_module.f90:147
The module contains the definition of a "node descriptor" as well as other global variables used duri...
Definition: amr_module.f90:21
integer, parameter cornyhi
y-coordinate of the upper border of this grid
Definition: amr_module.f90:149
Here is the call graph for this function:
Here is the caller graph for this function: