2D AMRCLAW
Functions/Subroutines
preintcopy.f File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine preintcopy (val, mitot, mjtot, nvar, ilo, ihi, jlo, jhi, level, fliparray)
 

Function/Subroutine Documentation

◆ preintcopy()

subroutine preintcopy ( dimension(nvar,mitot,mjtot)  val,
  mitot,
  mjtot,
  nvar,
  ilo,
  ihi,
  jlo,
  jhi,
  level,
dimension((mitot+mjtot)*nghost*nvar)  fliparray 
)

Definition at line 6 of file preintcopy.f.

References amr_module::hxposs, amr_module::hyposs, iadd(), intcopy(), amr_module::iregsz, amr_module::jregsz, amr_module::nghost, amr_module::spheredom, amr_module::xlower, amr_module::xperdom, amr_module::ylower, and amr_module::yperdom.

Referenced by filval().

6 c
7  use amr_module
8  implicit double precision (a-h,o-z)
9 
10  dimension fliparray((mitot+mjtot)*nghost*nvar)
11  dimension val(nvar,mitot,mjtot)
12  dimension ist(3), iend(3), jst(3), jend(3), ishift(3), jshift(3)
13 
14 c NEW INDEXING ORDER SWITCHED
15  iadd(ivar,i,j) = locflip + ivar-1 + nvar*((j-1)*nr+i-1)
16 c
17 c :::::::::::::: PREINTCOPY :::::::::::::::::::::::::::::::::::::::::::
18 c For periodic boundary conditions more work needed to initialize a
19 c new grid that sticks out. This routine was
20 c called because the patch sticks out of the domain,
21 c and has periodic bc.s preprocess the patch before calling
22 c intcopy to shift the patch periodically back into the domain.
23 c
24 c Inputs to this routine:
25 c ilo,ihi,jlo,jhi = the location in index space of this patch.
26 c
27 c Outputs from this routine:
28 c The values of the grid are inserted
29 c directly into the enlarged val array for this piece.
30 c
31 c :::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
32 c
33  locflip = 1
34 
35 c
36 c # will divide patch into 9 possibilities (some empty):
37 c x sticks out left, x interior, x sticks out right
38 c same for y. for example, the max. would be
39 c i from (ilo,-1), (0,iregsz(level)-1), (iregsz(level),ihi)
40 
41  if (xperdom) then
42  ist(1) = ilo
43  ist(2) = 0
44  ist(3) = iregsz(level)
45  iend(1) = -1
46  iend(2) = iregsz(level)-1
47  iend(3) = ihi
48  ishift(1) = iregsz(level)
49  ishift(2) = 0
50  ishift(3) = -iregsz(level)
51  else ! if not periodic, set vals to only have nonnull intersection for interior region
52  ist(1) = iregsz(level)
53  ist(2) = ilo
54  ist(3) = iregsz(level)
55  iend(1) = -1
56  iend(2) = ihi
57  iend(3) = -1
58  ishift(1) = 0
59  ishift(2) = 0
60  ishift(3) = 0
61  endif
62 
63 
64  if (yperdom .or. spheredom) then
65  jst(1) = jlo
66  jst(2) = 0
67  jst(3) = jregsz(level)
68  jend(1) = -1
69  jend(2) = jregsz(level)-1
70  jend(3) = jhi
71  jshift(1) = jregsz(level)
72  jshift(2) = 0
73  jshift(3) = -jregsz(level)
74  else
75  jst(1) = jregsz(level)
76  jst(2) = jlo
77  jst(3) = jregsz(level)
78  jend(1) = -1
79  jend(2) = jhi
80  jend(3) = -1
81  jshift(1) = 0
82  jshift(2) = 0
83  jshift(3) = 0
84  endif
85 
86 
87  do 20 i = 1, 3
88  i1 = max(ilo, ist(i))
89  i2 = min(ihi, iend(i))
90  if (i1 .gt. i2) go to 20
91  do 10 j = 1, 3
92  j1 = max(jlo, jst(j))
93  j2 = min(jhi, jend(j))
94 
95 
96  if (j1 .le. j2) then ! part of patch in this region
97 c
98 c check if special mapping needed for spherical bc.
99 c(j=2 is interior,nothing special needed)
100  if (.not. spheredom .or. j .eq. 2) then
101  iputst = (i1 - ilo) + 1
102  jputst = (j1 - jlo) + 1
103  call intcopy(val,mitot,mjtot,nvar,
104  2 i1+ishift(i),i2+ishift(i),
105  3 j1+jshift(j),j2+jshift(j),level,
106  4 iputst,jputst)
107  else
108  nr = i2 - i1 + 1
109  nc = j2 - j1 + 1
110  ng = 0 ! no ghost cells in this little patch. fill everything.
111 
112  jbump = 0
113  if (j1 < 0) jbump = abs(j1) ! bump up so new bottom is 0
114  if (j2 >= jregsz(level)) jbump = -(j2+1-jregsz(level)) ! bump down
115 
116 c next 2 lines would take care of periodicity in x
117  iwrap1 = i1 + ishift(i)
118  iwrap2 = i2 + ishift(i)
119 c next 2 lines take care of mapped sphere bcs
120  iwrap1 = iregsz(level) - iwrap1 -1
121  iwrap2 = iregsz(level) - iwrap2 -1
122 c swap so that smaller one is left index, etc since mapping reflects
123  tmp = iwrap1
124  iwrap1 = iwrap2
125  iwrap2 = tmp
126 
127  jwrap1 = j1 + jbump
128  xlwrap = xlower + iwrap1*hxposs(level)
129  ybwrap = ylower + jwrap1*hyposs(level)
130  jwrap2 = j2 + jbump
131 
132 c write(dbugunit,101) i1,i2,j1,j2
133 c write(dbugunit6,102) iwrap1,iwrap2,j1+jbump,j2+jbump
134  101 format(" actual patch from i:",2i5," j :",2i5)
135  102 format(" intcopy called w i:",2i5," j :",2i5)
136  call intcopy(fliparray,nr,nc,nvar,
137  1 iwrap1,iwrap2,jwrap1,jwrap2,level,1,1)
138 
139 c copy back using weird mapping for spherical folding
140  nrowst = 1 ! start filling up val at (1,1) - no additional offset
141  ncolst = 1
142  do 15 ii = i1, i2
143  do 15 jj = j1, j2
144 c write(dbugunit6,100)nrowst+ii-ilo,ncolst+jj-jlo,nr-(ii-i1),
145 c 1 nc-jj+j1
146  100 format(" filling loc ",2i5," with ",2i5)
147 
148  do 15 ivar = 1, nvar
149  iindex = nr-(ii-i1)
150  jindex = nc-(jj-j1)
151  index = iadd(ivar,nr-(ii-i1),nc-(jj-j1))
152  val(ivar,nrowst+(ii-ilo),ncolst+(jj-jlo)) =
153  1 fliparray(iadd(ivar,nr-(ii-i1),nc-(jj-j1)))
154  15 continue
155 
156 
157  endif
158 
159  endif
160 
161  10 continue
162  20 continue
163 
164 
165  return
integer, dimension(maxlv) iregsz
Definition: amr_module.f90:198
real(kind=8), dimension(maxlv) hyposs
Definition: amr_module.f90:193
real(kind=8), dimension(maxlv) hxposs
Definition: amr_module.f90:193
integer, dimension(maxlv) jregsz
Definition: amr_module.f90:198
real(kind=8) xlower
Definition: amr_module.f90:231
integer pure function iadd(ivar, i, j)
Definition: intfil.f90:294
logical yperdom
Definition: amr_module.f90:230
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 mjtot c...
Definition: intcopy.f:14
real(kind=8) ylower
Definition: amr_module.f90:231
logical spheredom
Definition: amr_module.f90:230
logical xperdom
Definition: amr_module.f90:230
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
Here is the call graph for this function:
Here is the caller graph for this function: