2D AMRCLAW
regrid.f
Go to the documentation of this file.
1 c
2 c -----------------------------------------------------------
3 c
4  subroutine regrid (nvar,lbase,cut,naux,start_time)
5 c
6  use amr_module
7  implicit double precision (a-h,o-z)
8  integer newnumgrids(maxlv)
9  integer clock_start2, clock_finish, clock_rate
10 c
11 c :::::::::::::::::::::::::::: REGRID :::::::::::::::::::::::::::::::
12 
20 c
21 c input parameters:
22 c lbase = highest level that stays fixed during regridding
23 c cutoff = criteria for measuring goodness of rect. fit.
24 
25 c local variables:
26 c lcheck = the level being examined.
27 c lfnew = finest grid to be. will replace lfine.
28 
29 c global
30 c mstart = start of very coarsest grids.
31 c :::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
32 c
34  lcheck = min0(lfine,mxnest-1)
35  lfnew = lbase
36  do 10 i = 1, mxnest
37  newnumgrids(i) = 0
38  10 newstl(i) = 0
39  time = rnode(timemult, lstart(lbase))
40 c
41  20 if (lcheck .lt. lbase) go to 50
42  call grdfit(lbase,lcheck,nvar,naux,cut,time,start_time)
43  if (newstl(lcheck+1) .eq. 0) go to 40
44  lfnew = max0(lcheck + 1,lfnew)
45  40 continue
46  lcheck = lcheck - 1
47 c
48  go to 20
49  50 continue
50 c
51 c end of level loop
52 c
53 c remaining tasks left in regridding:
54 c 1. count number of new grids at each level
55  maxnumnewgrids = 0 ! max over all levels. needed for dimensioning
56  do lev = lbase+1,lfnew
57  ngridcount = 0
58  mptr = newstl(lev)
59  52 if (mptr .eq. 0) go to 55
60  ngridcount = ngridcount + 1
61  mptr = node(levelptr,mptr)
62  go to 52
63 
64  55 newnumgrids(lev) = ngridcount
65  maxnumnewgrids = max(maxnumnewgrids,ngridcount)
66  end do
67 c
68 c 2. interpolate storage for the new grids. the starting pointers
69 c for each level are in newstl. also reclaim some space before new
70 c allocations.
71  call system_clock(clock_start2,clock_rate)
72  call gfixup(lbase,lfnew,nvar,naux,newnumgrids,maxnumnewgrids)
73  call system_clock(clock_finish,clock_rate)
74  timegrdfit2 = timegrdfit2 + clock_finish - clock_start2
75 c
76 c 3. merge data structures (newstl and lstart )
77 c finish storage allocation, reclaim space, etc. set up boundary
78 c flux conservation arrays
79 c
80  do 60 level = lbase, lfine-1
81  call prepf(level+1,nvar,naux)
82  call prepc(level,nvar)
83  60 continue
84 c
85 c reset numgrids per level, needed for omp parallelization.
86 c note that grids may have disappeared, so next loop resets to 0
87 c if there are no grids from lfine+1 to mxnest
88 c
89  do 72 levnew = lbase+1, mxnest
90  mptr = lstart(levnew)
91  ngridcount = 0
92  ncells = 0
93  do while (mptr .gt. 0)
94  ngridcount = ngridcount + 1
95  ncells = ncells + (node(ndihi,mptr)-node(ndilo,mptr)+1)
96  . * (node(ndjhi,mptr)-node(ndjlo,mptr)+1)
97  mptr = node(levelptr, mptr)
98  end do
99  numgrids(levnew) = ngridcount
100  numcells(levnew) = ncells
101  avenumgrids(levnew) = avenumgrids(levnew) + ngridcount
102  iregridcount(levnew) = iregridcount(levnew) + 1
103 c sort grids to first ones are the most work. this helps load
104 c balancing, but doesn't help locality
105  if (ngridcount .gt. 1) call arrangegrids(levnew,ngridcount)
106 
107  if (verbosity_regrid .ge. levnew) then
108  write(*,100) ngridcount,ncells,levnew
109  write(outunit,100) ngridcount,ncells,levnew
110  100 format("there are ",i6," grids with ",i10,
111  & " cells at level ", i3)
112  endif
113 72 continue
114 c
115 c set up array of grids instead of recomputing at each step
116  call makegridlist(lbase)
117  do levnew = lbase+1, lfine
118  call makebndrylist(levnew) ! does one level at a time
119  end do
120 
121  return
122  end
123 c
124 c -------------------------------------------------------------------
125 c
131 
132  subroutine arrangegrids(level, numg)
133 c
134  use amr_module
135  implicit double precision (a-h,o-z)
136  integer listgrids(numg), cost(numg), index(numg), prevptr
137 c
138 c slow sort for now, putting most expensive grids first on lstart list
139 c measure cost by number of cells
140 c
141  mptr = lstart(level)
142  do i = 1, numg
143  listgrids(i) = mptr
144  cost(i) = (node(ndihi,mptr)-node(ndilo,mptr)+3) *
145  1 (node(ndjhi,mptr)-node(ndjlo,mptr)+2)
146  index(i) = i
147  mptr = node(levelptr, mptr)
148  end do
149 c
150 c write(*,*)" before sorting"
151 c write(*,*) index
152 c
153  call qsorti(index, numg, cost)
154 
155 c write(*,*)"after sorting"
156 c write(*,*) index
157 
158 c qsort returns in ascending order, repack in descending order
159 c grids can stay in place, just their levelptrs need to change
160  lstart(level) = listgrids(index(numg)) ! last grid is most expensive
161  prevptr = listgrids(index(numg))
162  do i = 1, numg-1
163  node(levelptr, prevptr) = listgrids(index(numg-i))
164  prevptr = listgrids(index(numg-i))
165  end do
166  node(levelptr,prevptr) = null !signal the last grid
167 
168  return
169  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, dimension(7) method
Definition: amr_module.f90:253
integer, parameter timemult
current simulation time on this grid
Definition: amr_module.f90:151
integer timegrdfit2
Definition: amr_module.f90:240
subroutine regrid(nvar, lbase, cut, naux, start_time)
Definition: regrid.f:5
integer, dimension(maxlv) numcells
Definition: amr_module.f90:198
subroutine prepc(level, nvar)
This routine is called because regridding just changed the fine grids.
Definition: prepc.f:12
integer, parameter maxlv
Definition: amr_module.f90:174
subroutine qsorti(ORD, N, A)
Definition: quick_sort1.f:23
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, dimension(maxlv) newstl
Definition: amr_module.f90:198
integer, parameter ndihi
global i index of right border of this grid
Definition: amr_module.f90:111
integer, dimension(nsize, maxgr) node
Definition: amr_module.f90:198
integer, dimension(maxlv) numgrids
Definition: amr_module.f90:198
real(kind=8), dimension(rsize, maxgr) rnode
Definition: amr_module.f90:193
integer, dimension(maxlv) iregridcount
Definition: amr_module.f90:238
real(kind=8), dimension(maxlv) avenumgrids
Definition: amr_module.f90:237
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
subroutine grdfit(lbase, lcheck, nvar, naux, cut, time, start_time)
This subroutine is called by setgrd and regrid to actually fit the new grids on each level...
Definition: grdfit2.f:22
subroutine gfixup(lbase, lfnew, nvar, naux, newnumgrids, maxnumnewgrids)
Interpolate initial values for the newly created grids, whose levels start from level lbase+1...
Definition: gfixup.f:8
integer, parameter outunit
Definition: amr_module.f90:290
integer, parameter ndjhi
global j index of upper border of this grid
Definition: amr_module.f90:117
subroutine arrangegrids(level, numg)
Sort all grids at level level.
Definition: regrid.f:133
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 mxnest
Definition: amr_module.f90:198
integer verbosity_regrid
Definition: amr_module.f90:259
The module contains the definition of a "node descriptor" as well as other global variables used duri...
Definition: amr_module.f90:21
integer lfine
Definition: amr_module.f90:198
subroutine prepf(level, nvar, naux)
For new fine grids (grids on level level), allocate space for saving fluxes at boundary (even if it's...
Definition: prepf.f:19