2D AMRCLAW
Functions/Subroutines
grdfit2.f File Reference

Go to the source code of this file.

Functions/Subroutines

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. More...
 

Function/Subroutine Documentation

◆ grdfit()

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.

lcheck is the level being error estimated so that lcheck+1 will be the level of the new grids.

It generates all grids on level lcheck+1 and ensure proper nesting.

Parameters
lbasebase AMR level for current refinement, which stays fixed. Note that lbase is always less or equal to lcheck
lcheckthe level being error estimated so that lcheck+1 will be the level of the new grids.
nvarnumber of equations for the system
nauxnumber of auxiliary variables
cutefficiency threshold for the clustering algorithm
timecurrent simulation time
start_timestart time of current simulation

Definition at line 22 of file grdfit2.f.

References amr_module::alloc, birect(), amr_module::cornxhi, amr_module::cornxlo, amr_module::cornyhi, amr_module::cornylo, flglvl2(), amr_module::gprint, amr_module::hxposs, amr_module::hyposs, amr_module::iinfinity, amr_module::intratx, amr_module::intraty, amr_module::iregend, amr_module::iregst, amr_module::iregsz, amr_module::jregend, amr_module::jregst, amr_module::jregsz, amr_module::levelptr, amr_module::maxcl, moment(), amr_module::ndihi, amr_module::ndilo, amr_module::ndjhi, amr_module::ndjlo, amr_module::nestlevel, amr_module::newstl, amr_module::node, nodget(), amr_module::nsize, amr_module::null, amr_module::outunit, reclam(), amr_module::rnode, amr_module::rsize, smartbis(), amr_module::timeflglvl, amr_module::timegrdfitall, amr_module::timemult, amr_module::xlower, and amr_module::ylower.

Referenced by moment(), regrid(), and setgrd().

22 c
23  use amr_module
24  implicit double precision (a-h,o-z)
25  integer clock_start, clock_finish, clock_rate
26  integer clock_start1
27 c
28  dimension corner(nsize,maxcl)
29  integer numptc(maxcl), prvptr
30  logical fit, cout
31  logical fit2, nestck2
32  data cout/.false./
33 c
34 c
35  call system_clock(clock_start,clock_rate)
36 
37 
38 c ### initialize region start and end indices for new level grids
39  iregst(lcheck+1) = iinfinity
40  jregst(lcheck+1) = iinfinity
41  iregend(lcheck+1) = -1
42  jregend(lcheck+1) = -1
43 
44 c ## flag all grids at given level based on error ests.
45 c ## npts is number of points actually colated - some
46 c ## flagged points turned off due to proper nesting requirement.
47 c ## (storage based on nptmax calculation however).
48 
49  call system_clock(clock_start1,clock_rate)
50  call flglvl2(nvar,naux,lcheck,nptmax,index,lbase,npts,start_time)
51  call system_clock(clock_finish,clock_rate)
52  timeflglvl = timeflglvl + clock_finish - clock_start1
53 
54  if (npts .eq. 0) go to 99
55 c
56  levnew = lcheck + 1
57  hxfine = hxposs(levnew)
58  hyfine = hyposs(levnew)
59 c
60 c ## call smart_bisect grid gen. to make the clusters
61 c till each cluster ok. needs scratch space.
62 c
63  idim = iregsz(lcheck)
64  jdim = jregsz(lcheck)
65 c lociscr = igetsp(idim+jdim)
66 c locjscr = lociscr + idim
67  call smartbis(alloc(index),npts,cut,numptc,nclust,lbase,
68  2 corner,idim,jdim)
69 c 2 corner,alloc(lociscr),alloc(locjscr),idim,jdim)
70 c call reclam(lociscr,idim+jdim)
71 
72  if (gprint) then
73  write(outunit,103) nclust
74  write(outunit,104) (icl, numptc(icl),icl=1,nclust)
75  103 format(' ',i4,' clusters after bisect')
76  104 format(' cluster ',i5,' has points: ',i8)
77  endif
78 c
79 c ## for each cluster, fit the actual grid, set up some data structures
80 c
81  50 ibase = 0
82  icl = 1
83  prvptr = null
84 c
85  70 mnew = nodget()
86 c if (lcheck .eq. 2 .and. (mnew .ne. 6 .and. mnew .ne. 7)) go to 69
87 c if (lcheck .eq. 1 .and. (mnew .ne. 3 .and. mnew .ne. 2 )) go to 69
88  75 call moment(node(1,mnew),alloc(index+2*ibase),numptc(icl),usage)
89 
90  if (gprint) write(outunit,100) icl,mnew,usage,numptc(icl)
91 100 format(' cluster ',i5,' new rect.',i5,
92  1 ' usage ',e12.5,' with ',i5,' pts.')
93 
94  node(ndilo,mnew) = node(ndilo,mnew)*intratx(lcheck)
95  node(ndjlo,mnew) = node(ndjlo,mnew)*intraty(lcheck)
96  node(ndihi,mnew) = (node(ndihi,mnew)+1)*intratx(lcheck) - 1
97  node(ndjhi,mnew) = (node(ndjhi,mnew)+1)*intraty(lcheck) - 1
98  rnode(cornxlo,mnew) = node(ndilo,mnew)*hxfine + xlower
99  rnode(cornylo,mnew) = node(ndjlo,mnew)*hyfine + ylower
100  rnode(cornxhi,mnew) = (node(ndihi,mnew)+1)*hxfine + xlower
101  rnode(cornyhi,mnew) = (node(ndjhi,mnew)+1)*hyfine + ylower
102  node(nestlevel,mnew) = levnew
103  rnode(timemult,mnew) = time
104 
105  if (gprint) write(outunit,101) (node(i,mnew),i=1,nsize),
106  & (rnode(i,mnew),i=1,rsize)
107  101 format(4i5,4i15,/,4i15,5i15,/,2i15,/,5e15.7)
108 c
109 c ## if new grid doesn't fit in base grid, nestck bisect it
110 c ## and returns 2 clusters where there used to be 1.
111 c
112 c 2/28/02 : Added naux to argument list; needed by call to outtre in nestck
113 
114  fit2 = nestck2(mnew,lbase,alloc(index+2*ibase),numptc(icl),numptc,
115  1 icl,nclust,nvar, naux)
116  if (.not. fit2) go to 75
117 c
118 c ## grid accepted. put in list.
119  if (newstl(levnew) .eq. null) then
120  newstl(levnew) = mnew
121  else
122  node(levelptr,prvptr) = mnew
123  endif
124  prvptr = mnew
125 c # keep track of min and max location of grids at this level
126  iregst(levnew) = min(iregst(levnew), node(ndilo,mnew))
127  jregst(levnew) = min(jregst(levnew), node(ndjlo,mnew))
128  iregend(levnew) = max(iregend(levnew),node(ndihi,mnew))
129  jregend(levnew) = max(jregend(levnew),node(ndjhi,mnew))
130 
131 c ## on to next cluster
132  69 ibase = ibase + numptc(icl)
133  icl = icl + 1
134  if (icl .le. nclust) go to 70
135 c
136 c ## clean up. for all grids check final size.
137  call birect(newstl(levnew))
138 
139  99 continue
140 c ## may have npts 0 but array was allocated due to initially flagged points
141 c ## that were not allowed for proper nesting or other reasons. in this case
142 c ## the array was still allocated, so need to test further to see if colating
143 c ## array space needs to be reclaimed
144  if (nptmax .gt. 0) call reclam(index, 2*nptmax)
145 c
146  call system_clock(clock_finish,clock_rate)
147  timegrdfitall = timegrdfitall + clock_finish - clock_start
148 
149  return
integer, parameter timemult
current simulation time on this grid
Definition: amr_module.f90:151
integer, parameter rsize
Definition: amr_module.f90:30
integer, parameter cornxlo
x-coordinate of the left border of this grid
Definition: amr_module.f90:143
integer timegrdfitall
Definition: amr_module.f90:240
subroutine flglvl2(nvar, naux, lcheck, nxypts, index, lbase, npts, start_time)
Controls the error estimation/flagging bad pts.
Definition: flglvl2.f:26
subroutine reclam(index, nwords)
Definition: reclam.f:5
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) 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(maxlv) jregsz
Definition: amr_module.f90:198
integer, dimension(nsize, maxgr) node
Definition: amr_module.f90:198
integer, dimension(maxlv) iregst
Definition: amr_module.f90:198
integer, parameter nestlevel
AMR level of the grid.
Definition: amr_module.f90:44
logical gprint
Definition: amr_module.f90:297
real(kind=8), dimension(rsize, maxgr) rnode
Definition: amr_module.f90:193
real(kind=8) xlower
Definition: amr_module.f90:231
integer, parameter maxcl
maximum number of clusters (grids) on each grid level
Definition: amr_module.f90:177
integer, parameter nsize
Definition: amr_module.f90:31
real(kind=8) ylower
Definition: amr_module.f90:231
integer, dimension(maxlv) jregst
Definition: amr_module.f90:198
subroutine moment(intrect, badpts, npt, usage)
Compute enclosing rectangle around flagged points.
Definition: moment.f:17
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, parameter iinfinity
Definition: amr_module.f90:170
logical function nestck2(mnew, lbase, badpts, npts, numptc, icl, nclust, nvar, naux)
Check that the potential grid mnew is completely contained in the (coarser) finest grid which stays f...
Definition: nestck2.f:39
integer, dimension(maxlv) intraty
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, parameter outunit
Definition: amr_module.f90:290
subroutine smartbis(badpts, npts, cutoff, numptc, nclust, lbase, intcorn, idim, jdim)
Smart bisect rectangles until cutoff reached for each.
Definition: smartbis.f:27
integer, dimension(maxlv) jregend
Definition: amr_module.f90:198
integer, parameter ndjhi
global j index of upper border of this grid
Definition: amr_module.f90:117
subroutine birect(mptr1)
Check each grid, starting with mptr1 (either newstl or lstart) to see that it has no more than max1d ...
Definition: birect.f:11
integer, parameter cornylo
y-coordinate of the lower border of this grid
Definition: amr_module.f90:145
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
integer, parameter null
Definition: amr_module.f90:155
integer, parameter cornxhi
x-coordinate of the right border of this grid
Definition: amr_module.f90:147
integer timeflglvl
Definition: amr_module.f90:240
The module contains the definition of a "node descriptor" as well as other global variables used duri...
Definition: amr_module.f90:21
integer, dimension(maxlv) iregend
Definition: amr_module.f90:198
integer, parameter cornyhi
y-coordinate of the upper border of this grid
Definition: amr_module.f90:149
real(kind=8), dimension(:), allocatable alloc
Definition: amr_module.f90:218
Here is the call graph for this function:
Here is the caller graph for this function: