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