2D AMRCLAW
Functions/Subroutines
icall.f File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine icall (val, aux, nrow, ncol, nvar, naux, ilo, ihi, jlo, jhi, level, iputst, jputst)
 For a rectangle defined on level level and bound by ilo, ihi, jlo, jhi, find intersecting grids at the same level and copy data from intersecting grids to solution and auxiliary variables on the rectangle, stored in val and aux arrays. More...
 

Function/Subroutine Documentation

◆ icall()

subroutine icall ( dimension(nvar,nrow,ncol)  val,
dimension(naux,nrow,ncol)  aux,
  nrow,
  ncol,
  nvar,
  naux,
  ilo,
  ihi,
  jlo,
  jhi,
  level,
  iputst,
  jputst 
)

For a rectangle defined on level level and bound by ilo, ihi, jlo, jhi, find intersecting grids at the same level and copy data from intersecting grids to solution and auxiliary variables on the rectangle, stored in val and aux arrays.

Definition at line 10 of file icall.f.

References amr_module::alloc, iadd(), amr_module::levelptr, amr_module::lstart, amr_module::ndihi, amr_module::ndilo, amr_module::ndjhi, amr_module::ndjlo, amr_module::nghost, amr_module::node, amr_module::store1, and amr_module::storeaux.

Referenced by filval(), preicall(), and saveqc().

10 
11  use amr_module
12  implicit double precision (a-h, o-z)
13 
14  dimension val(nvar,nrow,ncol)
15  dimension aux(naux,nrow,ncol)
16 
17  logical sticksout
18 
19 
20 c new index ordering
21  iadd(ivar,i,j) = loc + ivar-1 + nvar*((j-1)*mitot+i-1)
22  iaddaux(ivar,i,j) = locaux + ivar-1 + naux*((j-1)*mitot+i-1)
23 
24 c ::::::::::::::::::::::::::: icall :::::::::::::::::::::::::::::::
25 c
26 c find intersecting grids at the same level. copy data from
27 c intersecting grids to both val and aux arrays.
28 c
29 c use larger definition of grids here - boundary data already in.
30 c aux arrays also enlarged size.
31 c
32 c iputst, jputst: where to copy values into. may not be in
33 c location corresponding to ilo,ihi,etc. if
34 c the patch has been periodically wrapped.
35 
36 c
37 c :::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
38 
39 
40  mptr = lstart(level)
41 
42  10 if (mptr .eq. 0) go to 99
43  iglo = node(ndilo,mptr)
44  ighi = node(ndihi,mptr)
45  jglo = node(ndjlo,mptr)
46  jghi = node(ndjhi,mptr)
47 
48 c # does it intersect?
49 c$$$ ixlo = max(iglo-nghost,ilo)
50 c$$$ ixhi = min(ighi+nghost,ihi)
51 c$$$ jxlo = max(jglo-nghost,jlo)
52 c$$$ jxhi = min(jghi+nghost,jhi)
53 c how did ghost cells get in the allowable region? they are not filled
54 c(since we may be interpolating from newly filled grids, not just grids
55 c that have been primed with bcs to be advanced.
56  ixlo = max(iglo,ilo)
57  ixhi = min(ighi,ihi)
58  jxlo = max(jglo,jlo)
59  jxhi = min(jghi,jhi)
60 
61 
62  if (ixlo .le. ixhi .and. jxlo .le. jxhi) then
63  loc = node(store1,mptr)
64  locaux = node(storeaux,mptr)
65  nx = ighi - iglo + 1
66  ny = jghi - jglo + 1
67  mitot = nx + 2*nghost
68  mjtot = ny + 2*nghost
69  do 30 j = jxlo, jxhi
70  do 30 i = ixlo, ixhi
71  do 20 ivar = 1, nvar
72  ialloc = iadd(ivar,i-iglo+nghost+1,j-jglo+nghost+1)
73  val(ivar,i-ilo+iputst,j-jlo+jputst) = alloc(ialloc)
74  20 continue
75  do 25 iaux = 1, naux
76  ialloc = iaddaux(iaux,i-iglo+nghost+1,j-jglo+nghost+1)
77  aux(iaux,i-ilo+iputst,j-jlo+jputst) = alloc(ialloc)
78  25 continue
79  30 continue
80  endif
81  mptr = node(levelptr, mptr)
82  go to 10
83 
84  99 continue
85 
86 c if cells stick out of domain but not periodic then set elsewhere
87 c either setaux and bc2amr. (called from routine that called this, e.g.
88 c saveqc or filval)
89 
90  return
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
subroutine saveqc(level, nvar, naux)
Prepare new fine grids to save fluxes after each integration step for future conservative fix-up on c...
Definition: saveqc.f:14
subroutine filval(val, mitot, mjtot, dx, dy, level, time, mic, mjc, xleft, xright, ybot, ytop, nvar, mptr, ilo, ihi, jlo, jhi, aux, naux)
Fill grid mptr on level level by copying values from OLD level level grids if available, otherwise by interpolating values from coarser grids.
Definition: filval.f90:6
logical pure function sticksout(iplo, iphi, jplo, jphi)
Definition: filpatch.f90:349
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 pure function iadd(ivar, i, j)
Definition: intfil.f90:294
subroutine icall(val, aux, nrow, ncol, nvar, naux, ilo, ihi, jlo, jhi, level, iputst, jputst)
For a rectangle defined on level level and bound by ilo, ihi, jlo, jhi, find intersecting grids at th...
Definition: icall.f:10
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
integer, parameter store1
pointer to the address of memory storing the first copy of solution data on this grid, usually for storing new solution
Definition: amr_module.f90:101
integer, dimension(maxlv) lstart
Definition: amr_module.f90:198
integer, parameter ndjhi
global j index of upper border of this grid
Definition: amr_module.f90:117
integer, parameter levelptr
node number (index) of next grid on the same level
Definition: amr_module.f90:35
subroutine domain(nvar, vtime, nx, ny, naux, start_time)
Definition: domain.f:5
integer nghost
Definition: amr_module.f90:232
The module contains the definition of a "node descriptor" as well as other global variables used duri...
Definition: amr_module.f90:21
subroutine setaux(mbc, mx, my, xlower, ylower, dx, dy, maux, aux)
Definition: setaux.f90:2
integer, parameter storeaux
pointer to the address of memory storing auxiliary data on this grid
Definition: amr_module.f90:120
real(kind=8), dimension(:), allocatable alloc
Definition: amr_module.f90:218
Here is the call graph for this function:
Here is the caller graph for this function: