2D AMRCLAW
Functions/Subroutines
gfixup.f File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine gfixup (lbase, lfnew, nvar, naux, newnumgrids, maxnumnewgrids)
 Interpolate initial values for the newly created grids, whose levels start from level lbase+1. More...
 
subroutine prepnewgrids (listnewgrids, num, level)
 

Function/Subroutine Documentation

◆ gfixup()

subroutine gfixup (   lbase,
  lfnew,
  nvar,
  naux,
integer, dimension(maxlv)  newnumgrids,
  maxnumnewgrids 
)

Interpolate initial values for the newly created grids, whose levels start from level lbase+1.

Definition at line 8 of file gfixup.f.

References amr_module::alloc, amr_module::cornxhi, amr_module::cornxlo, amr_module::cornyhi, amr_module::cornylo, filval(), freebndrylist(), amr_module::hxposs, amr_module::hyposs, igetsp(), amr_module::intratx, amr_module::intraty, amr_module::levelptr, amr_module::lfine, amr_module::lstart, amr_module::mxnest, amr_module::ndihi, amr_module::ndilo, amr_module::ndjhi, amr_module::ndjlo, amr_module::newstl, amr_module::nghost, amr_module::node, prepnewgrids(), putnod(), putsp(), reclam(), amr_module::rnode, amr_module::store1, amr_module::store2, amr_module::storeaux, and amr_module::timemult.

Referenced by regrid().

8 c
9  use amr_module
10  implicit double precision (a-h,o-z)
11 
12  integer omp_get_thread_num, omp_get_max_threads
13  integer mythread/0/, maxthreads/1/
14  integer newnumgrids(maxlv), listnewgrids(maxnumnewgrids)
15 
16 c
17 c ::::::::::::::::::::::::: GFIXUP ::::::::::::::::::::::::::::::::;
18 c interpolate initial values for the newly created grids.
19 c the start of each level is located in newstl array.
20 c since only levels greater than lbase were examined, start
21 c looking there.
22 c :::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::;
23 c
24 c reclaim old storage (position 8) and list space 15 and 16
25 c before allocating new storage. remember, finest level grids
26 c (if level = mxnest so that error never estimated) don't have
27 c 2 copies of solution values at old and new times.
28 c
29 c
30  call putsp(lbase,lbase,nvar,naux)
31  level = lbase + 1
32  1 if (level .gt. lfine) go to 4
33  call putsp(lbase,level,nvar,naux)
34  mptr = lstart(level)
35  2 if (mptr .eq. 0) go to 3
36  nx = node(ndihi,mptr) - node(ndilo,mptr) + 1
37  ny = node(ndjhi,mptr) - node(ndjlo,mptr) + 1
38  mitot = nx + 2*nghost
39  mjtot = ny + 2*nghost
40  nwords = mitot*mjtot*nvar
41  if (level .lt. mxnest)
42  . call reclam(node(store2, mptr), nwords)
43  node(store2, mptr) = 0
44  mptr = node(levelptr, mptr)
45  go to 2
46  3 level = level + 1
47  go to 1
48 c
49  4 lcheck = lbase + 1
50 
51  time = rnode(timemult, lstart(lbase))
52  5 if (lcheck .gt. mxnest) go to 99
53  hx = hxposs(lcheck)
54  hy = hyposs(lcheck)
55 
56 c
57 c prepare for doing next loop over grids at a given level in parallel
58 c unlike other level loops, these are newly created grids, not yet merged in
59 c so take grids from newstl (NEWSTartOfLevel), not lstart. Dont yet know
60 c how many either.
61  call prepnewgrids(listnewgrids,newnumgrids(lcheck),lcheck)
62 c
63 c interpolate level lcheck
64 c first get space, since cant do that part in parallel
65  do j = 1, newnumgrids(lcheck)
66  mptr = listnewgrids(j)
67  nx = node(ndihi,mptr) - node(ndilo,mptr) + 1
68  ny = node(ndjhi,mptr) - node(ndjlo,mptr) + 1
69  mitot = nx + 2*nghost
70  mjtot = ny + 2*nghost
71  loc = igetsp(mitot * mjtot * nvar)
72  node(store1, mptr) = loc
73  if (naux .gt. 0) then
74  locaux = igetsp(mitot * mjtot * naux)
75  else
76  locaux = 1
77  endif
78  node(storeaux, mptr) = locaux
79  end do
80 
81 c
82 !$OMP PARALLEL DO
83 !$OMP& PRIVATE(j,mptr,nx,ny,mitot,mjtot,corn1,corn2,loc)
84 !$OMP& PRIVATE(locaux,time,mic,mjc,xl,xr,yb,yt,ilo,ihi)
85 !$OMP& PRIVATE(jlo,jhi,sp_over_h,thisSetauxTime)
86 !$OMP& SHARED(newnumgrids,listnewgrids,nghost,node,hx,hy)
87 !$OMP& SHARED(rnode,intratx,intraty,lcheck,nvar,alloc,naux)
88 !$OMP& SCHEDULE(dynamic,1)
89 !$OMP& DEFAULT(none)
90 c
91  do j = 1, newnumgrids(lcheck)
92  mptr = listnewgrids(j)
93 
94 c changed to move setaux out of this loop. instead, copy aux in filval
95 c along with soln.involves changing intcopy to icall and making flag array
96 c can only do this after topo stops moving
97  nx = node(ndihi,mptr) - node(ndilo,mptr) + 1
98  ny = node(ndjhi,mptr) - node(ndjlo,mptr) + 1
99  mitot = nx + 2*nghost
100  mjtot = ny + 2*nghost
101  corn1 = rnode(cornxlo,mptr)
102  corn2 = rnode(cornylo,mptr)
103  loc = node(store1, mptr)
104  if (naux .gt. 0) then
105  locaux = node(storeaux, mptr)
106  else
107  locaux = 1
108  endif
109 
110 c
111 c We now fill in the values for grid mptr using filval. It uses
112 c piecewise linear interpolation to obtain values from the
113 c (lcheck - 1) grid, then overwrites those with whatever (lcheck)
114 c grids are available. We take advantage of the fact that the
115 c (lcheck - 1) grids have already been set, and that the distance
116 c between any point in mptr and a (lcheck - 1) and (lcheck - 2)
117 c interface is at least one (lcheck - 1) cell wide.
118 c
119 
120 c # make a coarsened patch with ghost cells so can use
121 c # grid interpolation routines, but only set "interior".
122 c # extra 2 cells so that can use linear interp. on
123 c # "interior" of coarser patch to fill fine grid.
124  mic = nx/intratx(lcheck-1) + 2
125  mjc = ny/intraty(lcheck-1) + 2
126  xl = rnode(cornxlo,mptr)
127  xr = rnode(cornxhi,mptr)
128  yb = rnode(cornylo,mptr)
129  yt = rnode(cornyhi,mptr)
130  ilo = node(ndilo, mptr)
131  ihi = node(ndihi, mptr)
132  jlo = node(ndjlo, mptr)
133  jhi = node(ndjhi, mptr)
134 
135  call filval(alloc(loc),mitot,mjtot,hx,hy,lcheck,time,
136  1 mic,mjc,
137  2 xl,xr,yb,yt,nvar,
138  3 mptr,ilo,ihi,jlo,jhi,
139  4 alloc(locaux),naux)
140 
141  end do
142 c
143 c done filling new grids at level. move them into lstart from newstl
144 c (so can use as source grids for filling next level). can also
145 c get rid of loc. 7 storage for old level.
146 c
147  80 mptr = lstart(lcheck)
148  85 if (mptr .eq. 0) go to 90
149  nx = node(ndihi,mptr) - node(ndilo,mptr) + 1
150  ny = node(ndjhi,mptr) - node(ndjlo,mptr) + 1
151  mitot = nx + 2*nghost
152  mjtot = ny + 2*nghost
153  call reclam(node(store1,mptr),mitot*mjtot*nvar)
154  if (naux .gt. 0) then
155  call reclam(node(storeaux,mptr),mitot*mjtot*naux)
156  endif
157  mold = mptr
158  mptr = node(levelptr,mptr)
159  call putnod(mold)
160  call freebndrylist(mold)
161  go to 85
162  90 lstart(lcheck) = newstl(lcheck)
163  lcheck = lcheck + 1
164  go to 5
165 c
166  99 lfine = lfnew
167 c
168 c initialize 2nd (old time) storage block for new grids not at top level
169 c
170  levend = min(lfine,mxnest-1)
171  do 110 level = lbase+1, levend
172  mptr = lstart(level)
173  105 if (mptr .eq. 0) go to 110
174  nx = node(ndihi,mptr) - node(ndilo,mptr) + 1
175  ny = node(ndjhi,mptr) - node(ndjlo,mptr) + 1
176  mitot = nx + 2*nghost
177  mjtot = ny + 2*nghost
178  nwords = mitot*mjtot*nvar
179  node(store2,mptr) = igetsp(nwords)
180  mptr = node(levelptr,mptr)
181  go to 105
182  110 continue
183 
184 c
185 c -------------
186 c grid structure now complete again. safe to print, etc. assuming
187 c things initialized to zero in nodget.
188 c -------------
189 c
190  return
integer, parameter timemult
current simulation time on this grid
Definition: amr_module.f90:151
integer, parameter cornxlo
x-coordinate of the left border of this grid
Definition: amr_module.f90:143
function igetsp(nwords)
Allocate contiguous space of length nword in main storage array alloc.
Definition: igetsp.f:9
integer, parameter maxlv
Definition: amr_module.f90:174
subroutine filval(val, mitot, mjtot, dx, dy, level, time, mic, mjc, xleft, xright, ybot, ytop, nvar, mptr, ilo, ihi, jlo, jhi, aux, naux)
Fill grid mptr on level level by copying values from OLD level level grids if available, otherwise by interpolating values from coarser grids.
Definition: filval.f90:6
subroutine reclam(index, nwords)
Definition: reclam.f:5
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
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
real(kind=8), dimension(rsize, maxgr) rnode
Definition: amr_module.f90:193
subroutine putnod(mptr)
Return mptr node to the linked list kept in node array.
Definition: putnod.f:6
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 store1
pointer to the address of memory storing the first copy of solution data on this grid, usually for storing new solution
Definition: amr_module.f90:101
integer, dimension(maxlv) lstart
Definition: amr_module.f90:198
integer, parameter store2
pointer to the address of memory storing the second copy of solution data on this grid...
Definition: amr_module.f90:105
integer, dimension(maxlv) intraty
Definition: amr_module.f90:198
subroutine prepnewgrids(listnewgrids, num, level)
Definition: gfixup.f:201
integer, parameter ndjhi
global j index of upper border of this grid
Definition: amr_module.f90:117
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 mxnest
Definition: amr_module.f90:198
integer, parameter cornxhi
x-coordinate of the right border of this grid
Definition: amr_module.f90:147
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
integer lfine
Definition: amr_module.f90:198
integer, parameter cornyhi
y-coordinate of the upper border of this grid
Definition: amr_module.f90:149
subroutine putsp(lbase, level, nvar, naux)
Reclaim list space in nodes cfluxptr and ffluxptr for all grids at level level
Definition: putsp.f:7
integer, parameter storeaux
pointer to the address of memory storing auxiliary data on this grid
Definition: amr_module.f90:120
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:

◆ prepnewgrids()

subroutine prepnewgrids ( integer, dimension(num)  listnewgrids,
  num,
  level 
)

Definition at line 201 of file gfixup.f.

References amr_module::levelptr, amr_module::newstl, and amr_module::node.

Referenced by gfixup().

201 
202  use amr_module
203  implicit double precision (a-h,o-z)
204  integer listnewgrids(num)
205 
206  mptr = newstl(level)
207  do j = 1, num
208  listnewgrids(j) = mptr
209  mptr = node(levelptr, mptr)
210  end do
211 
212  if (mptr .ne. 0) then
213  write(*,*)" Error in routine setting up grid array "
214  stop
215  endif
216 
217  return
integer, dimension(maxlv) newstl
Definition: amr_module.f90:198
integer, dimension(nsize, maxgr) node
Definition: amr_module.f90:198
integer, parameter levelptr
node number (index) of next grid on the same level
Definition: amr_module.f90:35
The module contains the definition of a "node descriptor" as well as other global variables used duri...
Definition: amr_module.f90:21
Here is the caller graph for this function: