2D AMRCLAW
Functions/Subroutines
prefilp_orig.f90 File Reference

Go to the source code of this file.

Functions/Subroutines

recursive subroutine prefilrecur (level, nvar, valbig, aux, naux, time, mitot, mjtot, nrowst, ncolst, ilo, ihi, jlo, jhi, fullGrid)
 
integer pure function iadd (n, i, j)
 
integer pure function iaddscratch (n, i, j)
 

Function/Subroutine Documentation

◆ iadd()

integer pure function prefilrecur::iadd ( integer, intent(in)  n,
integer, intent(in)  i,
integer, intent(in)  j 
)

Definition at line 188 of file prefilp_orig.f90.

References iadd().

188  implicit none
189  integer, intent(in) :: n, i, j
190  iadd = locflip + n-1 + nvar*((j-1)*nr+i-1)
integer pure function iadd(ivar, i, j)
Definition: intfil.f90:294
Here is the call graph for this function:

◆ iaddscratch()

integer pure function prefilrecur::iaddscratch ( integer, intent(in)  n,
integer, intent(in)  i,
integer, intent(in)  j 
)

Definition at line 194 of file prefilp_orig.f90.

References iaddscratch().

194  implicit none
195  integer, intent(in) :: n, i, j
196  iaddscratch = n + nvar*((j-1)*nr+i-1) ! no subtract 1
integer pure function iaddscratch(n, i, j)
Definition: prefilp.f90:207
Here is the call graph for this function:

◆ prefilrecur()

recursive subroutine prefilrecur ( integer, intent(in)  level,
integer, intent(in)  nvar,
real(kind=8), dimension(nvar,mitot,mjtot), intent(inout)  valbig,
real(kind=8), dimension(naux,mitot,mjtot), intent(inout)  aux,
integer, intent(in)  naux,
real(kind=8), intent(in)  time,
integer, intent(in)  mitot,
integer, intent(in)  mjtot,
integer, intent(in)  nrowst,
integer, intent(in)  ncolst,
integer, intent(in)  ilo,
integer, intent(in)  ihi,
integer, intent(in)  jlo,
integer, intent(in)  jhi,
logical  fullGrid 
)

Definition at line 20 of file prefilp_orig.f90.

References filrecur(), amr_module::hxposs, amr_module::hyposs, iaddscratch(), amr_module::iregsz, amr_module::jregsz, amr_module::needs_to_be_set, amr_module::nghost, setaux(), amr_module::spheredom, amr_module::xlower, amr_module::xperdom, amr_module::ylower, and amr_module::yperdom.

20 
21 
22 
25 
26  !for setaux timing
27  use amr_module, only: timesetaux, timesetauxcpu
28 
29  implicit none
30 
31  ! Input
32  integer, intent(in) :: level, nvar, naux, mitot, mjtot, nrowst, ncolst
33  integer, intent(in) :: ilo,ihi,jlo,jhi
34  real(kind=8), intent(in) :: time
35  logical :: fullgrid ! true first time called, false for recursive coarse sub-patches
36 
37  ! Output
38  real(kind=8), intent(in out) :: valbig(nvar,mitot,mjtot)
39  real(kind=8), intent(in out) :: aux(naux,mitot,mjtot)
40 
41  ! Local storage
42  integer :: i, j, ii, jj, ivar, nr, nc, ng, i1, i2, j1, j2, iputst, jputst
43  integer :: jbump, iwrap1, iwrap2, jwrap1, tmp, locflip, rect(4)
44  real(kind=8) :: xlwrap, ybwrap
45 
46  integer :: ist(3), iend(3), jst(3), jend(3), ishift(3), jshift(3)
47  real(kind=8) :: scratch(max(mitot,mjtot)*nghost*nvar)
48  real(kind=8) :: scratchaux(max(mitot,mjtot)*nghost*naux)
49 
50  !for timings
51  integer :: clock_start, clock_finish, clock_rate
52  real(kind=8) :: cpu_start, cpu_finish
53 
54 ! # will divide patch into 9 possibilities (some empty):
55 ! x sticks out left, x interior, x sticks out right
56 ! same for y. for example, the max. would be
57 ! i from (ilo,-1), (0,iregsz(level)-1), (iregsz(level),ihi)
58 
59  if (xperdom) then
60  ist(1) = ilo
61  ist(2) = 0
62  ist(3) = iregsz(level)
63  iend(1) = -1
64  iend(2) = iregsz(level)-1
65  iend(3) = ihi
66  ishift(1) = iregsz(level)
67  ishift(2) = 0
68  ishift(3) = -iregsz(level)
69  else ! if not periodic, set vals to only have nonnull intersection for interior regoin
70  ist(1) = iregsz(level)
71  ist(2) = ilo
72  ist(3) = iregsz(level)
73  iend(1) = -1
74  iend(2) = ihi
75  iend(3) = -1
76  ishift(1) = 0
77  ishift(2) = 0
78  ishift(3) = 0
79  endif
80 
81  if (yperdom .or. spheredom) then
82  jst(1) = jlo
83  jst(2) = 0
84  jst(3) = jregsz(level)
85  jend(1) = -1
86  jend(2) = jregsz(level)-1
87  jend(3) = jhi
88  jshift(1) = jregsz(level)
89  jshift(2) = 0
90  jshift(3) = -jregsz(level)
91  else
92  jst(1) = jregsz(level)
93  jst(2) = jlo
94  jst(3) = jregsz(level)
95  jend(1) = -1
96  jend(2) = jhi
97  jend(3) = -1
98  jshift(1) = 0
99  jshift(2) = 0
100  jshift(3) = 0
101  endif
102 
103 ! ## loop over the 9 regions (in 2D) of the patch - the interior is i=j=2 plus
104 ! ## the ghost cell regions. If any parts stick out of domain and are periodic
105 ! ## map indices periodically, but stick the values in the correct place
106 ! ## in the orig grid (indicated by iputst,jputst.
107 ! ## if a region sticks out of domain but is not periodic, not handled in (pre)-icall
108 ! ## but in setaux/bcamr (not called from here).
109  do 20 i = 1, 3
110  i1 = max(ilo, ist(i))
111  i2 = min(ihi, iend(i))
112  if (i1 .gt. i2) go to 20
113  do 10 j = 1, 3
114  j1 = max(jlo, jst(j))
115  j2 = min(jhi, jend(j))
116 
117  ! part of patch in this region
118  if (j1 <= j2) then
119 
120  ! there is something to fill. j=2 case is interior, no special
121  ! mapping needed even if spherical bc
122  if (.not. spheredom .or. j == 2 ) then
123  iputst = (i1 - ilo) + nrowst
124  jputst = (j1 - jlo) + ncolst
125 
126  call filrecur(level,nvar,valbig,aux,naux,time,mitot,mjtot, &
127  iputst,jputst,i1+ishift(i),i2+ishift(i),j1+jshift(j),j2+jshift(j),.false.)
128  else
129  nr = i2 - i1 + 1
130  nc = j2 - j1 + 1
131  ng = 0 ! no ghost cells in this little patch. fill everything.
132 
133  jbump = 0
134  if (j1 < 0) jbump = abs(j1) ! bump up so new bottom is 0
135  if (j2 >= jregsz(level)) jbump = -(j2+1-jregsz(level)) ! bump down
136 
137  ! next 2 lines would take care of periodicity in x
138  iwrap1 = i1 + ishift(i)
139  iwrap2 = i2 + ishift(i)
140  ! next 2 lines take care of mapped sphere bcs
141  iwrap1 = iregsz(level) - iwrap1 -1
142  iwrap2 = iregsz(level) - iwrap2 -1
143  ! swap so that smaller one is left index, etc since mapping reflects
144  tmp = iwrap1
145  iwrap1 = iwrap2
146  iwrap2 = tmp
147 
148  jwrap1 = j1 + jbump
149  xlwrap = xlower + iwrap1*hxposs(level)
150  ybwrap = ylower + jwrap1*hyposs(level)
151 
152  if (naux>0) then
153  scratchaux = needs_to_be_set !flag all cells with signal since dimensioned strangely
154 
155  call system_clock(clock_start,clock_rate)
156  call cpu_time(cpu_start)
157  call setaux(ng,nr,nc,xlwrap,ybwrap,hxposs(level),hyposs(level),naux,scratchaux)
158  call system_clock(clock_finish,clock_rate)
159  call cpu_time(cpu_finish)
160  timesetaux = timesetaux + clock_finish - clock_start
161  timesetauxcpu = timesetauxcpu + cpu_finish - cpu_start
162  endif
163 
164  rect = [iwrap1,iwrap2,j1+jbump,j2+jbump]
165  call filrecur(level,nvar,scratch,scratchaux,naux,time,nr, &
166  nc,1,1,iwrap1,iwrap2,j1+jbump,j2+jbump,.false.)
167 
168  ! copy back using weird mapping for spherical folding
169  do ii = i1, i2
170  do jj = j1, j2
171  ! write(dbugunit,'(" filling loc ",2i5," with ",2i5)') nrowst+ii-fill_indices(1),ncolst+jj-fill_indices(3),nr-(ii-i1),nc-jj+j1
172 
173  do ivar = 1, nvar
174  valbig(ivar,nrowst+(ii-ilo),ncolst+(jj-jlo)) = &
175  scratch(iaddscratch(ivar,nr-(ii-i1),nc-(jj-j1)))
176  end do
177  ! write(dbugunit,'(" new val is ",4e15.7)')(valbig(ivar,nrowst+(ii-fill_indices(1)),ncolst+(jj-fill_indices(3))),ivar=1,nvar)
178  end do
179  end do
180  endif
181  endif
182  10 continue
183  20 continue
184 
185 contains
186 
187  integer pure function iadd(n,i,j)
188  implicit none
189  integer, intent(in) :: n, i, j
190  iadd = locflip + n-1 + nvar*((j-1)*nr+i-1)
191  end function iadd
192 
193  integer pure function iaddscratch(n,i,j)
194  implicit none
195  integer, intent(in) :: n, i, j
196  iaddscratch = n + nvar*((j-1)*nr+i-1) ! no subtract 1
197  end function iaddscratch
198 
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
integer pure function iaddscratch(n, i, j)
Definition: prefilp.f90:207
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
real(kind=8) ylower
Definition: amr_module.f90:231
real(kind=8), parameter needs_to_be_set
Definition: amr_module.f90:168
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
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
subroutine setaux(mbc, mx, my, xlower, ylower, dx, dy, maux, aux)
Definition: setaux.f90:2
Here is the call graph for this function: