2D AMRCLAW
Functions/Subroutines
restrt.f File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine restrt (nsteps, time, nvar, naux)
 
subroutine initrestoftimers ()
 

Function/Subroutine Documentation

◆ initrestoftimers()

subroutine initrestoftimers ( )

Definition at line 255 of file restrt.f.

References amr_module::timebufnst, amr_module::timeflagger, amr_module::timeflglvl, amr_module::timegrdfit2, and amr_module::timeupdating.

Referenced by restrt().

255  use amr_module
256 c
257 c :::::::::::::::::::::::::::: RESTRT ::::::::::::::::::::::::::::::::
258 c
259 c some timers checkpointed and reinitalized. Set rest to zero to prevent
260 c NaNs. These timers used for developers, not printed for all users, so
261 c not deleted.
262 c
263 c ::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
264  timebufnst = 0.d0
265  timefilval = 0.d0
266  timeflagger = 0.d0
267  timeflglvl = 0.d0
268  timeflglvltot = 0.d0
269  timegrdfit2 = 0.d0
270  timegridfitall = 0.d0
271  timeupdating = 0.d0
272 
273  return
integer timegrdfit2
Definition: amr_module.f90:240
integer timeflagger
Definition: amr_module.f90:242
integer timeupdating
Definition: amr_module.f90:239
integer timebufnst
Definition: amr_module.f90:242
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
Here is the caller graph for this function:

◆ restrt()

subroutine restrt (   nsteps,
  time,
  nvar,
  naux 
)

read back in the check point files written by subr. check.

Definition at line 6 of file restrt.f.

References amr_module::alloc, amr_module::avenumgrids, amr_module::cflmax, amr_module::cfluxptr, amr_module::check_a, amr_module::evol, amr_module::ffluxptr, amr_module::flag_richardson, amr_module::hxposs, amr_module::hyposs, amr_module::icheck, igetsp(), initbndrylist(), initrestoftimers(), amr_module::intratx, amr_module::intraty, amr_module::iorder, amr_module::iregend, amr_module::iregridcount, amr_module::iregst, amr_module::iregsz, amr_module::jregend, amr_module::jregst, amr_module::jregsz, amr_module::kcheck, amr_module::kratio, amr_module::lendim, amr_module::lenf, amr_module::lenmax, amr_module::lentot, amr_module::levelptr, amr_module::lfine, amr_module::lfree, amr_module::listsp, amr_module::liststart, amr_module::lstart, makebndrylist(), makegridlist(), amr_module::matlabu, amr_module::maxlv, amr_module::memsize, amr_module::mstart, amr_module::mxnest, amr_module::ndfree, amr_module::ndfree_bnd, amr_module::ndihi, amr_module::ndilo, amr_module::ndjhi, amr_module::ndjlo, amr_module::newstl, amr_module::nghost, amr_module::node, amr_module::numgrids, amr_module::outunit, amr_module::possk, putnod(), reclam(), restrt_alloc(), amr_module::rnode, amr_module::rstfile, amr_module::rstunit, amr_module::rvol, amr_module::rvoll, amr_module::store1, amr_module::store2, amr_module::storeaux, amr_module::timebound, amr_module::timeboundcpu, amr_module::timeregridding, amr_module::timeregriddingcpu, amr_module::timestepgrid, amr_module::timestepgridcpu, amr_module::timetick, amr_module::timetickcpu, amr_module::timevalout, amr_module::timevaloutcpu, amr_module::tmass0, amr_module::tvoll, and amr_module::tvollcpu.

Referenced by amr2().

6 c
7  use amr_module
8  implicit double precision (a-h,o-z)
9  logical ee
10 
11 
12  logical foundfile
13  dimension intrtx(maxlv),intrty(maxlv),intrtt(maxlv)
14 c
15 c :::::::::::::::::::::::::::: RESTRT ::::::::::::::::::::::::::::::::
17 c
18 c some input variables might have changed, and also the
19 c alloc array could have been written with a smaller size at checkpoint
20 c ::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
21 c
22 c !! Now allow user-specified file name !!
23 c rstfile = 'restart.data'
24 
25  ! If checkpt_style < 0 then alternating between two checkpoint files.
26  ! Set which one to use for the first checkpoint after this restart.
27  ! Set check_a to .true. unless fort.chkaaaaa is the file being read.
28  ! When alternating checkpoint files used, this keeps proper sequence going
29  ! otherwise (checkpt_style > 0) check_a is not used elsewhere.
30  check_a = .not. (rstfile == 'fort.chkaaaaa')
31 
32  write(6,*) 'Attempting to restart computation using '
33  write(6,*) ' checkpoint file: ',trim(rstfile)
34  inquire(file=trim(rstfile),exist=foundfile)
35  if (.not. foundfile) then
36  write(*,*)" Did not find checkpoint file!"
37  stop
38  endif
39  open(rstunit,file=trim(rstfile),status='old',form='unformatted')
40  rewind rstunit
41 
42  read(rstunit) lenmax,lendim,isize
43 
44 c # need to allocate for dynamic memory:
45  call restrt_alloc(isize)
46 
47  read(rstunit) (alloc(i),i=1,lendim)
49  read(rstunit) lfree,lenf
51  1 ibuf,mstart,ndfree,ndfree_bnd,lfine,iorder,mxnold,
52  2 intrtx,intrty,intrtt,iregsz,jregsz,
54  3 numgrids,kcheck1,nsteps,time,
55  3 matlabu
63 
64  close(rstunit)
65 
66  write(outunit,100) nsteps,time
67  write(6,100) nsteps,time
68  100 format(/,' RESTARTING the calculation after ',i5,' steps',
69  1 /,' (time = ',e15.7,')')
70 
71 c
72 c error checking that refinement ratios have not changed
73 c ### new feature: when using variable refinement in time
74 c ### (varRefTime = T) the time ratios are allowed to be different
75 c ### (since they are ignored and calc. on the fly)
76 c ### This is not checked for here, since the same amr2.f is now
77 c used for geoclaw too, and varRefTime is only available in geoclaw.
78 c
79  do i = 1, min(mxnold-1,mxnest-1)
80  if ( (intratx(i) .ne. intrtx(i)) .or.
81  . (intraty(i) .ne. intrty(i)) ) then
82 c . (kratio(i) .ne. intrtt(i) .and. .not. varRefTime) ) then
83  write(outunit,*)
84  . " not allowed to change existing refinement ratios on Restart"
85  write(outunit,*)" Old ratios:"
86  write(*,*) " Old ratios:"
87  write(outunit,903)(intrtx(j),j=1,mxnold-1)
88  write(*,903) (intrtx(j),j=1,mxnold-1)
89  write(outunit,903)(intrty(j),j=1,mxnold-1)
90  write(*,903) (intrty(j),j=1,mxnold-1)
91 c write(outunit,903)(intrtt(j),j=1,mxnold-1)
92 c write(*,903) (intrtt(j),j=1,mxnold-1)
93  903 format(6i3)
94  stop
95  endif
96  end do
97 
98 c if (varRefTime) then ! reset intrat to previously saved ratios, not input ratios
99  if (.true.) then ! reset intrat to previously saved ratios, not input ratios
100  do i = 1, mxnold-1
101  kratio(i) = intrtt(i)
102  end do
103  endif
104 
105 c
106 c adjust free list of storage in case size has changed.
107 c
108  idif = memsize - isize
109  if (idif .gt. 0) then
110  lfree(lenf,1) = isize + 2
111  call reclam(isize+1,idif)
112  else if (idif .lt. 0) then
113  write(outunit,900) isize, memsize
114  write(*,900) isize, memsize
115  900 format(' size of alloc not allowed to shrink with ',/,
116  . ' restart old size ',i7,' current size ',i7)
117  stop
118  endif
119 c
120 c adjust storage in case mxnest has changed - only allow it to increase,
121 c
122  if (mxnest .eq. mxnold) go to 99
123 
124  if (mxnest .lt. mxnold) then
125  if (lfine .lt. mxnest) then
126  go to 99
127  else
128  write(outunit,901) mxnold, mxnest
129  write(*, 901) mxnold, mxnest
130 901 format(' mxnest reduced on restart: ',/,
131  & ' old mxnest ',i4, ' new mxnest ',i4)
132  write(outunit,*)" reclaiming finer levels from",
133  . mxnest+1," to ",mxnold
134  do 95 lev = mxnest,mxnold
135  mptr = lstart(lev)
136  if (lev .gt. mxnest) lstart(lev) = 0
137  85 if (mptr .eq. 0) go to 95
138  if (lev .lt. mxnold) then
139  call reclam(node(cfluxptr,mptr), 5*listsp(lev))
140  node(cfluxptr,mptr) = 0
141  endif
142  nx = node(ndihi,mptr) - node(ndilo,mptr) + 1
143  ny = node(ndjhi,mptr) - node(ndjlo,mptr) + 1
144  ikeep = nx/intrtx(lev-1)
145  jkeep = ny/intrty(lev-1)
146  lenbc = 2*(ikeep+jkeep)
147  if (lev .gt. mxnest) then
148  call reclam
149  . (node(ffluxptr,mptr),2*nvar*lenbc+naux*lenbc)
150  node(ffluxptr,mptr) = 0
151  endif
152  mitot = nx + 2*nghost
153  mjtot = ny + 2*nghost
154  if (lev .gt. mxnest) then ! if level going away take away first storage
155  call reclam(node(store1,mptr),mitot*mjtot*nvar)
156  node(store1,mptr) = 0
157  if (naux .gt. 0) then ! and aux arrays
158  call reclam(node(storeaux,mptr),mitot*mjtot*naux)
159  node(storeaux,mptr) = 0
160  endif
161  endif
162  if (lev .ge. mxnest .and. lev .lt. mxnold) then !reclam 2nd level storage too
163  call reclam(node(store2, mptr), mitot*mjtot*nvar)
164  node(store2,mptr) = 0
165  endif
166  mold = mptr
167  mptr = node(levelptr,mptr)
168  if (lev .gt. mxnest) call putnod(mold)
169  go to 85
170  95 end do
171 c stop
172 c reset lfine to new finest level
173  do lev = mxnest,1,-1
174  if (lstart(lev) .gt. 0) then
175  lfine = lev
176  write(*,*)" resetting finest level to ",lfine
177  go to 199
178  endif
179  end do
180  endif
181  endif
182 
183 c if mxnest has increased, add second storage loc to grids at previous mxnest
184  ee = .false.
185  do 10 level = 1, mxnold
186  if (icheck(level) .ge. kcheck) then
187  ee = .true.
188  endif
189 10 continue
190 
191  write(*,*)" increasing max num levels from ",mxnold,
192  . ' to',mxnest
193  write(outunit,*)" increasing max num levels from ",mxnold,
194  . ' to',mxnest
195 
196  if (ee .and. flag_richardson) then
197 c ## if Richardson used, will delay error est. 1 step til have old soln. vals
198  write(*,*)" first Richardson error estimation step"
199  write(*,*)" will estimate mostly spatial error "
200  write(outunit,*)" first Richardson error estimation step"
201  write(outunit,*)" will estimate mostly spatial error "
202  endif
203 
204 c # add second storage location to previous mxnest level
205  mptr = lstart(mxnold)
206 15 if (mptr .eq. 0) go to 25
207  mitot = node(ndihi,mptr)-node(ndilo,mptr)+1+2*nghost
208  mjtot = node(ndjhi,mptr)-node(ndjlo,mptr)+1+2*nghost
209  node(store2,mptr) = igetsp(mitot*mjtot*nvar)
210  mptr = node(levelptr,mptr)
211  go to 15
212 25 continue
213 c
214 c # add new info. to spatial and counting arrays
215  99 level = lfine + 1
216  rrk = dble(kratio(lfine))
217 35 if (level .gt. mxnest) go to 45
218  hxposs(level) = hxposs(level-1) / dble(intratx(level-1))
219  hyposs(level) = hyposs(level-1) / dble(intraty(level-1))
220  possk(level) = possk(level-1) / rrk
221  iregsz(level) = iregsz(level-1) * intratx(level-1)
222  jregsz(level) = jregsz(level-1) * intraty(level-1)
223  rrk = kratio(level)
224  level = level + 1
225  go to 35
226 45 continue
227 c
228 c
229  199 continue
230 
231 c
232 c save array of grids to avoid copying each advanc or update step
233 c lbase is 1 here, to start building from level 1
234 c only for level 1 is listStart set outside of makeGridList
235 c call with lbase 0 to make level 1
236  liststart(1) = 1
237  call makegridlist(0)
238 c
239 c bndry list for faster ghost cell filling
240  call initbndrylist()
241  do level = 1, lfine
242  call makebndrylist(level)
243  end do
244 c
245 c some timers were checkpointed and restarted.
246 c others not currently used, so initialize to 0 to prevent NaNs
247  call initrestoftimers()
248 
249  return
integer timetick
Definition: amr_module.f90:242
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
subroutine restrt_alloc(isize)
integer, dimension(maxlv) kratio
Definition: amr_module.f90:198
integer lendim
Definition: amr_module.f90:247
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
integer ndfree
Definition: amr_module.f90:198
real(kind=8) evol
Definition: amr_module.f90:237
integer timebound
Definition: amr_module.f90:241
integer lenmax
Definition: amr_module.f90:247
subroutine makegridlist(lbase)
Make (modify) array of grid numbers, listOfGrids, (after sorting them in the linked list so they are ...
Definition: nodget.f:109
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) tmass0
Definition: amr_module.f90:267
real(kind=8) cflmax
Definition: amr_module.f90:255
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 memsize
Definition: amr_module.f90:219
integer, dimension(0:maxlv+1) liststart
Definition: amr_module.f90:189
character(len=200) rstfile
Definition: amr_module.f90:311
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(maxlv) iregridcount
Definition: amr_module.f90:238
real(kind=8) timevaloutcpu
Definition: amr_module.f90:245
integer lentot
Definition: amr_module.f90:247
logical flag_richardson
Definition: amr_module.f90:258
subroutine putnod(mptr)
Return mptr node to the linked list kept in node array.
Definition: putnod.f:6
logical check_a
Definition: amr_module.f90:312
real(kind=8), dimension(maxlv) avenumgrids
Definition: amr_module.f90:237
subroutine initrestoftimers()
Definition: restrt.f:255
integer, dimension(maxlv) icheck
Definition: amr_module.f90:198
integer matlabu
Definition: amr_module.f90:283
integer, dimension(maxlv) jregst
Definition: amr_module.f90:198
integer iorder
Definition: amr_module.f90:198
integer, parameter rstunit
Definition: amr_module.f90:292
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
real(kind=8) rvol
Definition: amr_module.f90:237
real(kind=8) timetickcpu
Definition: amr_module.f90:243
real(kind=8) timeboundcpu
Definition: amr_module.f90:244
integer, parameter store2
pointer to the address of memory storing the second copy of solution data on this grid...
Definition: amr_module.f90:105
real(kind=8), dimension(maxlv) tvollcpu
Definition: amr_module.f90:243
integer, dimension(maxlv) intraty
Definition: amr_module.f90:198
real(kind=8) timeregriddingcpu
Definition: amr_module.f90:244
real(kind=8), dimension(maxlv) rvoll
Definition: amr_module.f90:237
integer, parameter outunit
Definition: amr_module.f90:290
real(kind=8), dimension(maxlv) possk
Definition: amr_module.f90:193
integer, dimension(lfdim, 2) lfree
Definition: amr_module.f90:225
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
integer timestepgrid
Definition: amr_module.f90:241
integer, dimension(maxlv) listsp
Definition: amr_module.f90:198
integer lenf
Definition: amr_module.f90:225
real(kind=8) timestepgridcpu
Definition: amr_module.f90:244
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 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, dimension(maxlv) iregend
Definition: amr_module.f90:198
integer timevalout
Definition: amr_module.f90:239
integer mstart
Definition: amr_module.f90:198
subroutine initbndrylist()
Definition: nodget.f:148
integer kcheck
Definition: amr_module.f90:198
integer, parameter cfluxptr
Pointer to an 5 by maxsp array, which has boundary information for this grid.
Definition: amr_module.f90:80
integer lfine
Definition: amr_module.f90:198
integer, parameter storeaux
pointer to the address of memory storing auxiliary data on this grid
Definition: amr_module.f90:120
integer, dimension(maxlv) tvoll
Definition: amr_module.f90:238
real(kind=8), dimension(:), allocatable alloc
Definition: amr_module.f90:218
integer timeregridding
Definition: amr_module.f90:239
integer, parameter ffluxptr
pointer to the address of memory storing fluxes in a layer around the grid, to be used in conservatio...
Definition: amr_module.f90:97
Here is the call graph for this function:
Here is the caller graph for this function: