2D AMRCLAW
nodget.f
Go to the documentation of this file.
1 c
2 c ------------------------------------------------------------
3 c
10  integer function nodget()
11 c
12  use amr_module
13  implicit double precision (a-h,o-z)
14 
15 c
16 c ::::::::::::::::: NODGET ::::::::::::::::::::::::::::::::::::;
17 c nodget = get first free node of the linked list kept in node
18 c array. adjust pointers accordingly.
19 c :::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::;
20 c
21  if (ndfree .ne. null) go to 10
22  write(outunit,100) maxgr
23  write(*,100) maxgr
24 100 format(' out of nodal space - allowed ',i8,' grids')
25  do level = 1, lfine
26  write(*,101) level,numgrids(level)
27  101 format(" level ",i4," has ",i7,' grids')
28  end do
29  write(*,*)" Could need twice as many grids as on any given"
30  write(*,*)" level if regridding/birecting"
31  stop
32 c
33 c update pointers
34 c
35  10 nodget = ndfree
37 c
38 c initialize new block
39 c
40  do 20 i = 1, nsize
41  node(i,nodget) = 0
42  20 continue
43 c
44  do 30 i = 1, rsize
45  rnode(i,nodget) = 0.0d0
46  30 continue
47 c
48  return
49  end
50 c
51 c ------------------------------------------------------------
52 c
53  integer function nodget_bnd()
54 c
55  use amr_module
56  implicit double precision (a-h,o-z)
57 
58 c
59 c ::::::::::::::::: NODGET_BND ::::::::::::::::::::::::::::::::::::;
60 c nodget_bnd = same as above but for bndry list
61 c :::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::;
62 c
63  if (ndfree_bnd .ne. null) go to 10
64  write(outunit,100) bndlistsize
65  write(*,100) bndlistsize
66 100 format(' out of bndry space - allowed ',i5,' bndry grids')
67  ! calc average number of bndry nbors per grid
68  nbortotal = 0
69  numgridstotal = 0
70  do lev = 1, lfine
71  numgridstotal = numgridstotal + numgrids(lev)
72  do mptr = 1, numgrids(lev)
73  nbortotal = nbortotal + node(bndlistnum,mptr)
74  end do
75  end do
76  avgnbors = float(nbortotal)/numgridstotal
77  write(*,101) numgridstotal,nbortotal,avgnbors
78  101 format(" There are ",i8," total grids", i10," bndry nbors",
79  . " average num/grid ",f10.3)
80 
81  stop
82 c
83 c ## adjust pointers
84 c
87 c
88 c ## initialize to 0
89 c
90  bndlist(nodget_bnd,1) = 0
91  bndlist(nodget_bnd,2) = 0
92 c
93  return
94  end
95 c
96 c -----------------------------------------------------------------
97 c
108  subroutine makegridlist(lbase)
109 c
110  use amr_module
111  implicit none
112 
113  integer lbase, levSt, lev, mptr, n
114 
115 c :::::::::::::::::::::::::::: make_gridList :::::::::::::::::::::::::
116 c make array of grid numbers (after sorting them so in decreasing
117 c order of workload, done in arrangeGrid and put back into linked
118 c list. Done every time there is regridding, initial gridding,
119 c or restarting. Most often finest level is regridded, so
120 c put it last in array. lbase is the level that didnt change, so
121 c only redo from lbase+1 to lfine.
122 c :::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
123 
124  !write(*,*)"mgl: lbase,lfine",lbase,lfine
125  do lev = lbase+1, lfine
126  levst = liststart(lev)
127  !write(*,*)"mgl: level ",lev," starts at ",levSt
128  mptr = lstart(lev)
129 c traverse linked list into array. list already sorted by arrangegrids
130  do n = 1, numgrids(lev)
131  listofgrids(levst+n-1) = mptr
132  mptr = node(levelptr,mptr)
133  end do
134 c
135 c next level starts one after where this one ends.
136 c Using a sentinel in dimension of
137 c listStart so no need to test if level = mxnest
138 c
139  liststart(lev+1) = levst + numgrids(lev)
140  end do
141 
142  return
143  end
144 c
145 c -----------------------------------------------------------------
146 c
147  subroutine initbndrylist()
149  use amr_module
150  implicit none
151 
152  integer i
153 c
154 c need to manage the boundary List too
155 c
156  do i = 1, bndlistsize
157  bndlist(i,nextfree) = i+1
158  end do
159 
161  ndfree_bnd = 1
162 
163  end
164 c
165 c -----------------------------------------------------------------
166 c
182  subroutine makebndrylist(level)
183 c
184  use amr_module
185  implicit none
186 
187  integer level, n, levSt, k, nborCount
188  integer nodget_bnd, nextSpot, prevNbor, msrc, mptr
189  integer imin, imax, jmin, jmax
190  integer imlo, imhi, jmlo, jmhi
191  integer ixlo, ixhi, jxlo, jxhi
192 
193 c :::::::::::::::::::::::::::: makeBndryList :::::::::::::::::::::::::
194 c preprocess each grid to have linked list of other grids at
195 c same level that supply ghost cells.
196 c :::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
197 
198 c traverse linked list into array. list already sorted by arrangegrids
199  levst = liststart(level)
200  do n = 1, numgrids(level)
201  mptr = listofgrids(levst+n-1)
202  imin = node(ndilo,mptr) - nghost ! ghost cells included since
203  imax = node(ndihi,mptr) + nghost ! this is what you want to fill
204  jmin = node(ndjlo,mptr) - nghost ! may also use for filval stuff,
205  jmax = node(ndjhi,mptr) + nghost ! change nghost to mbuff, etc
206  nborcount = 0
207 
208  do k = 1, numgrids(level) ! loop over all other grids once to find touching ones
209  if (k .eq. n) cycle ! dont count yourself as source grid
210  msrc = listofgrids(levst+k-1)
211 
212  ! Check if grid mptr and patch intersect
213  imlo = node(ndilo, msrc)
214  jmlo = node(ndjlo, msrc)
215  imhi = node(ndihi, msrc)
216  jmhi = node(ndjhi, msrc)
217 
218  ixlo = max(imlo,imin)
219  ixhi = min(imhi,imax)
220  jxlo = max(jmlo,jmin)
221  jxhi = min(jmhi,jmax)
222 
223  if (ixlo .le. ixhi .and. jxlo .le. jxhi) then ! put on bnd list for mptr
224  nborcount = nborcount + 1
225  nextspot = nodget_bnd()
226  bndlist(nextspot,gridnbor) = msrc
227  ! get spot in bnd list. insert next grid at front to avoid traversing
228  bndlist(nextspot,nextfree) = node(bndlistst,mptr)
229  node(bndlistst,mptr) = nextspot
230  endif
231 
232  end do
233 
234 ! save final count
235  node(bndlistnum,mptr) = nborcount
236  end do
237 
238  return
239  end
240 c
241 c -----------------------------------------------------------------
242 c
246  subroutine freebndrylist(mold)
247 c
248  use amr_module
249  implicit none
250 
251  integer nborCount, mold,nextSpot, i, nextnext
252 
253 c :::::::::::::::::::::::::::: freeBndryList :::::::::::::::::::::::::
254 c free the linked list of intersecting "boundary" grids for grid 'mold'
255 c that is no longer active
256 c :::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
257 
258  nborcount = node(bndlistnum,mold) ! count for this grid
259  nextspot = node(bndlistst,mold) ! first index of this grids nbors
260  do i = 1, nborcount
261  nextnext = bndlist(nextspot,nextfree)
262  call putnod_bnd(nextspot)
263  nextspot = nextnext
264  end do
265 
266  return
267  end
subroutine makebndrylist(level)
Preprocess each grid on level level to have a linked list of other grids at the same level that suppl...
Definition: nodget.f:183
integer, parameter rsize
Definition: amr_module.f90:30
integer ndfree
Definition: amr_module.f90:198
integer, parameter bndlistnum
number of grids (on the same level) that border this grid
Definition: amr_module.f90:139
integer, parameter maxgr
Definition: amr_module.f90:173
subroutine makegridlist(lbase)
Make (modify) array of grid numbers, listOfGrids, (after sorting them in the linked list so they are ...
Definition: nodget.f:109
integer, parameter ndihi
global i index of right border of this grid
Definition: amr_module.f90:111
subroutine freebndrylist(mold)
Free the linked list of intersecting "boundary" grids for grid 'mold' that is no longer active...
Definition: nodget.f:247
integer, dimension(nsize, maxgr) node
Definition: amr_module.f90:198
integer, dimension(0:maxlv+1) liststart
Definition: amr_module.f90:189
integer, dimension(maxlv) numgrids
Definition: amr_module.f90:198
real(kind=8), dimension(rsize, maxgr) rnode
Definition: amr_module.f90:193
integer ndfree_bnd
Definition: amr_module.f90:198
integer, dimension(maxgr) listofgrids
Definition: amr_module.f90:189
integer, parameter gridnbor
Definition: amr_module.f90:158
integer, parameter bndlistsize
Definition: amr_module.f90:190
integer, parameter bndlistst
pointer (actually it's an index in the bndList array) to the first node of a linked list...
Definition: amr_module.f90:136
integer, parameter nsize
Definition: amr_module.f90:31
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, dimension(maxlv) lstart
Definition: amr_module.f90:198
integer function nodget()
Get first free node of the linked list kept in node array.
Definition: nodget.f:11
integer, dimension(bndlistsize, 2) bndlist
Definition: amr_module.f90:191
subroutine putnod_bnd(mcell)
Definition: putnod.f:26
integer, parameter outunit
Definition: amr_module.f90:290
integer, parameter nextfree
Definition: amr_module.f90:154
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, parameter null
Definition: amr_module.f90:155
integer function nodget_bnd()
Definition: nodget.f:54
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 initbndrylist()
Definition: nodget.f:148
integer lfine
Definition: amr_module.f90:198