2D AMRCLAW
Functions/Subroutines
intcopy.f File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine intcopy (val, mitot, mjtot, nvar, ilo, ihi, jlo, jhi, level, iputst, jputst)
 For a rectangle that is on level level, described by ilo, ihi, jlo, jhi and made up by mitot $ \times$ mjtot cells, copy solution from all (OLD) grids on level level to val, which stores solution on this rectangle. More...
 

Function/Subroutine Documentation

◆ intcopy()

subroutine intcopy ( dimension(nvar,mitot,mjtot)  val,
  mitot,
  mjtot,
  nvar,
  ilo,
  ihi,
  jlo,
  jhi,
  level,
  iputst,
  jputst 
)

For a rectangle that is on level level, described by ilo, ihi, jlo, jhi and made up by mitot $ \times$ mjtot cells, copy solution from all (OLD) grids on level level to val, which stores solution on this rectangle.

Some portions of the rectangle may not be filled since they do not overlap with any level level grids.

Definition at line 14 of file intcopy.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, and amr_module::store1.

Referenced by filval(), and preintcopy().

14 
15  use amr_module
16  implicit double precision (a-h, o-z)
17 
18  dimension val(nvar,mitot,mjtot)
19 
20 
21 c OLD INDEXING
22 c iadd(i,j,ivar) = loc + i - 1 + mi*((ivar-1)*mj+j-1)
23 c NEW INDEXING ORDER SWITCHED
24  iadd(ivar,i,j) = loc + ivar-1 + nvar*((j-1)*mi+i-1)
25 
26 c ::::::::::::::::::::::::::: INTCOPY :::::::::::::::::::::::::::::::
27 c
28 c find intersecting grids at the same level. copy data from
29 c old grid to val.
30 c old grid has "nghost" ghost cells - passed in nodal common block.
31 c new grid has no ghost cells - indices describe entire patch.
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 c
36 c :::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
37 
38 
39  mptr = lstart(level)
40 
41  10 if (mptr .eq. 0) go to 99
42  iglo = node(ndilo,mptr)
43  ighi = node(ndihi,mptr)
44  jglo = node(ndjlo,mptr)
45  jghi = node(ndjhi,mptr)
46 
47 c # does it intersect?
48  ixlo = max(iglo,ilo)
49  ixhi = min(ighi,ihi)
50  jxlo = max(jglo,jlo)
51  jxhi = min(jghi,jhi)
52 
53  if (ixlo .le. ixhi .and. jxlo .le. jxhi) then
54  loc = node(store1,mptr)
55  nx = ighi - iglo + 1
56  ny = jghi - jglo + 1
57  mi = nx + 2*nghost
58  mj = ny + 2*nghost
59  do 20 j = jxlo, jxhi
60  do 20 ivar = 1, nvar
61  do 30 i = ixlo, ixhi
62  val(ivar,iputst+i-ilo,jputst+j-jlo) =
63  1 alloc(iadd(ivar,i-iglo+nghost+1,j-jglo+nghost+1))
64  30 continue
65  20 continue
66  endif
67  mptr = node(levelptr, mptr)
68  go to 10
69 
70  99 return
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
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
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
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: