2D AMRCLAW
prepc.f
Go to the documentation of this file.
1 c
9 c ----------------------------------------------------------
10 c
11  subroutine prepc(level,nvar)
12 c
13  use amr_module
14  implicit double precision (a-h,o-z)
15 
16 c
17 c :::::::::::::::::::: PREPC ::::::::::::::::::::::::::::::::::::::
18 c
19 c this routine called because regridding just changed the fine grids.
20 c modify coarse grid boundary lists to store fluxes in appropriate
21 c fine grids lists.
22 c assume new fine grids have node(cfluxptr) initialized to null
23 c
24 c first compute max. possible number of list cells. allocate
25 c initially so that one pass through is enough.
26 c
27 c ::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
28 c
29  maxsp = 0
30  mkid = lstart(level+1)
31  10 if (mkid .eq. 0) go to 20
32  ikeep = (node(ndihi,mkid)-node(ndilo,mkid)+1)/intratx(level)
33  jkeep = (node(ndjhi,mkid)-node(ndjlo,mkid)+1)/intraty(level)
34  maxsp = maxsp + 2*(ikeep+jkeep)
35  mkid = node(levelptr,mkid)
36  go to 10
37  20 listsp(level) = maxsp
38  if (maxsp .eq. 0) go to 99
39 c
40  hxpar = hxposs(level)
41  hypar = hyposs(level)
42  hxkid = hxposs(level+1)
43  hykid = hyposs(level+1)
44  imax = iregsz(level) - 1
45  jmax = jregsz(level) - 1
46 
47  mpar = lstart(level)
48  30 if (mpar .eq. 0) go to 99
49 c
50  ispot = 0
51  ilo = node(ndilo,mpar)
52  jlo = node(ndjlo,mpar)
53  ihi = node(ndihi,mpar)
54  jhi = node(ndjhi,mpar)
55  locbc = igetsp(5*maxsp)
56 c # initialize list to 0 (0 terminator indicates end of bc list)
57  do 35 i = 1,5*maxsp
58  35 alloc(locbc+i-1) = 0.d0
59  node(cfluxptr,mpar) = locbc
60 c
61  mkid = lstart(level+1)
62  40 if (mkid .eq. 0) go to 60
63 
64  iclo = node(ndilo,mkid)/intratx(level)
65  jclo = node(ndjlo,mkid)/intraty(level)
66  ichi = node(ndihi,mkid)/intratx(level)
67  jchi = node(ndjhi,mkid)/intraty(level)
68 
69  iplo = max(ilo,iclo)
70  jplo = max(jlo,jclo)
71  iphi = min(ihi,ichi)
72  jphi = min(jhi,jchi)
73 
74 c regular intersections (will check in setuse that no duplication)
75 c this first call is only interior interfaces.
76 
77  if (iplo .le. iphi+1 .and. jplo .le. jphi+1) then
78  kflag = 1 ! interior stuff, no mappings
79  call setuse(alloc(locbc),maxsp,ispot,mkid,
80  2 ilo,ihi,jlo,jhi,iclo,ichi,jclo,jchi,kflag)
81  endif
82 
83 c for fine grids touching periodic boundary on right
84  if (xperdom .and. ilo .eq. 0 .and. ichi .eq. imax) then
85  kflag = 1 ! periodic in x
86  call setuse(alloc(locbc),maxsp,ispot,mkid,
87  2 ilo,ihi,jlo,jhi,iclo-iregsz(level),ichi-iregsz(level),
88  3 jclo,jchi,kflag)
89  endif
90 
91 c for fine grids touching periodic boundary on left
92  if (xperdom .and. iclo .eq. 0 .and. ihi .eq. imax) then
93  kflag = 1
94  call setuse(alloc(locbc),maxsp,ispot,mkid,
95  2 ilo,ihi,jlo,jhi,iclo+iregsz(level),ichi+iregsz(level),
96  3 jclo,jchi,kflag)
97  endif
98 
99 c for fine grids touching periodic boundary on top
100  if (yperdom .and. jlo .eq. 0 .and. jchi .eq. jmax) then
101  kflag = 1
102  call setuse(alloc(locbc),maxsp,ispot,mkid,
103  2 ilo,ihi,jlo,jhi,iclo,ichi,
104  3 jclo-jregsz(level),jchi-jregsz(level),kflag)
105  endif
106 
107 c for fine grids touching periodic boundary on bottom
108  if (yperdom .and. jclo .eq. 0 .and. jhi .eq. jmax) then
109  kflag = 1
110  call setuse(alloc(locbc),maxsp,ispot,mkid,
111  2 ilo,ihi,jlo,jhi,iclo,ichi,
112  3 jclo+jregsz(level),jchi+jregsz(level),kflag)
113  endif
114 
115 c for fine grids touching boundary on top in spherically mapped case
116 c and coarse grid touches top too. see if (mapped) x extent overlap.
117  if (spheredom .and. jhi .eq. jmax .and. jchi .eq. jmax) then
118  kflag = 2
119 c write(dbugunit,*)" for coarse grid ",mpar
120  iwrap2 = iregsz(level) - iclo - 1 !higher mapped index
121  iwrap1 = iregsz(level) - ichi - 1 !lower mapped index
122  if (max(ilo,iwrap1) .le. min(ihi,iwrap2)) then
123  call setuse(alloc(locbc),maxsp,ispot,mkid,
124  1 ilo,ihi,jlo,jhi,iclo,ichi,
125  2 jclo,jchi,kflag)
126  endif
127  endif
128 
129 c fine grids touching boundary on bottom for spherically mapped case
130 c coarse grid touches bottom too. see if (mapped) x extents overlap
131  if (spheredom .and. jclo .eq. 0 .and. jlo .eq. 0) then
132  kflag = 3
133  iwrap2 = iregsz(level) - iclo - 1 !higher mapped index
134  iwrap1 = iregsz(level) - ichi - 1 !lower mapped index
135  if (max(ilo,iwrap1) .le. min(ihi,iwrap2)) then
136  call setuse(alloc(locbc),maxsp,ispot,mkid,
137  1 ilo,ihi,jlo,jhi,iclo,ichi,
138  2 jclo,jchi,kflag)
139  endif
140  endif
141 
142  50 mkid = node(levelptr,mkid)
143  go to 40
144 c
145 c done with subgrid cycle. if no cells would need fixing, all done
146 c else cycle through again to set up list with info. for bc processing
147 c
148  60 continue
149 c
150 c for now, leave unused space allocated to the grid. alternative is to
151 c return (maxsp-ispot) amt starting at loc node(cfluxptr,mpar)+ispot.
152 c
153  mpar = node(levelptr,mpar)
154  go to 30
155 c
156  99 return
157  end
function igetsp(nwords)
Allocate contiguous space of length nword in main storage array alloc.
Definition: igetsp.f:9
subroutine prepc(level, nvar)
This routine is called because regridding just changed the fine grids.
Definition: prepc.f:12
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, parameter ndihi
global i index of right border of this grid
Definition: amr_module.f90:111
integer, dimension(maxlv) jregsz
Definition: amr_module.f90:198
integer, dimension(nsize, maxgr) node
Definition: amr_module.f90:198
logical yperdom
Definition: amr_module.f90:230
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
integer, dimension(maxlv) lstart
Definition: amr_module.f90:198
integer, dimension(maxlv) intraty
Definition: amr_module.f90:198
integer, parameter ndjhi
global j index of upper border of this grid
Definition: amr_module.f90:117
integer, dimension(maxlv) listsp
Definition: amr_module.f90:198
logical xperdom
Definition: amr_module.f90:230
integer, dimension(maxlv) intratx
Definition: amr_module.f90:198
integer, parameter levelptr
node number (index) of next grid on the same level
Definition: amr_module.f90:35
subroutine setuse(listbc, maxsp, ispot, mkid, ilo, ihi, jlo, jhi, iclo, ichi, jclo, jchi, kflag)
Add intersection information between grid mptr and a finer grid mkid to the boundary list...
Definition: setuse.f:64
The module contains the definition of a "node descriptor" as well as other global variables used duri...
Definition: amr_module.f90:21
integer, parameter cfluxptr
Pointer to an 5 by maxsp array, which has boundary information for this grid.
Definition: amr_module.f90:80
real(kind=8), dimension(:), allocatable alloc
Definition: amr_module.f90:218