2D AMRCLAW
tick.f
Go to the documentation of this file.
1 c
2 c -------------------------------------------------------------
3 c
4  subroutine tick(nvar,cut,nstart,vtime,time,naux,start_time,
5  & rest,dt_max)
6 c
7  use amr_module
10 
11  implicit double precision (a-h,o-z)
12 c include "call.i"
13 
14  logical vtime, dumpout/.false./, dumpchk/.false./
15  logical rest, dump_final
16  dimension dtnew(maxlv), ntogo(maxlv), tlevel(maxlv)
17  integer clock_start, clock_finish, clock_rate
18  integer tick_clock_start, tick_clock_finish, tick_clock_rate
19 
20 
21 c
22 c :::::::::::::::::::::::::::: TICK :::::::::::::::::::::::::::::
28 
29 c parameters:
30 c nstop = # of coarse grid time steps to be taken
31 c iout = output interval every 'iout' coarse time steps
32 c (if 0, not used - set to inf.)
33 c vtime = true for variable timestep, calculated each coarse step
34 c
35 c integration strategy is to advance a fine grid until it catches
36 c up to the coarse grid. this strategy is applied recursively.
37 c coarse grid goes first.
38 c
39 c nsteps: used to count how number steps left for a level to be
40 c integrated before it catches up with the next coarser level.
41 c ncycle: counts number of coarse grid steps = # cycles.
42 c
43 c icheck: counts the number of steps (incrementing by 1
44 c each step) to keep track of when that level should
45 c have its error estimated and finer levels should be regridded.
46 c ::::::::::::::::::::::::::::::::::::;::::::::::::::::::::::::::
47 c
48  call system_clock(tick_clock_start,tick_clock_rate)
49  call cpu_time(tick_cpu_start)
50 
51 
52  ncycle = nstart
53  call setbestsrc() ! need at very start of run, including restart
54  if (iout .eq. 0) then
55 c # output_style 1 or 2
56  iout = iinfinity
57  nextout = 0
58  if (nout .gt. 0) then
59  nextout = 1
60  if (nstart .gt. 0) then
61 c # restart: make sure output times start after restart time
62  do ii = 1, nout
63  if (tout(ii) .gt. time) then
64  nextout = ii
65  go to 2
66  endif
67  end do
68  2 continue
69  endif
70  endif
71  endif
72 
73  nextchk = 1
74  if ((nstart .gt. 0) .and. (abs(checkpt_style).eq.2)) then
75 c if this is a restart, make sure chkpt times start after restart time
76  do ii = 1, nchkpt
77  if (tchk(ii) .gt. time) then
78  nextchk = ii
79  go to 3
80  endif
81  enddo
82  3 continue
83  endif
84 
85  tlevel(1) = time
86 
87  do 5 i = 2, mxnest
88  tlevel(i) = tlevel(1)
89  5 continue
90 c
91 c ------ start of coarse grid integration loop. ------------------
92 c
93  20 if (ncycle .ge. nstop .or. time .ge. tfinal) goto 999
94 
95  if (nout .gt. 0) then
96  if (nextout .le. nout) then
97  outtime = tout(nextout)
98  else
99  outtime = rinfinity
100  endif
101  else
102  outtime = tfinal
103  endif
104 
105  if (nextchk .le. nchkpt) then
106  chktime = tchk(nextchk)
107  else
108  chktime = rinfinity
109  endif
110 
111  dumpout = .false. !# may be reset below
112 
113  if (time.lt.outtime .and. time+1.001*possk(1) .ge. outtime) then
114 c ## adjust time step to hit outtime exactly, and make output
115 c # apr 2010 mjb: modified so allow slightly larger timestep to
116 c # hit output time exactly, instead of taking minuscule timestep
117 c # should still be stable since increase dt in only 3rd digit.
118  oldposs = possk(1)
119  possk(1) = outtime - time
120 c write(*,*)" old possk is ", possk(1)
121  diffdt = oldposs - possk(1) ! if positive new step is smaller
122 
123 
124  if (.false.) then
125  write(*,122) diffdt,outtime ! notify of change
126  122 format(" Adjusting timestep by ",e10.3,
127  . " to hit output time of ",e13.6)
128 c write(*,*)" new possk is ", possk(1)
129  if (diffdt .lt. 0.) then ! new step is slightly larger
130  pctincrease = -100.*diffdt/oldposs ! minus sign to make whole expr. positive
131  write(*,123) pctincrease
132  123 format(" New step is ",f8.2," % larger.",
133  . " Should still be stable")
134  endif
135  endif
136 
137 
138  do i = 2, mxnest
139  possk(i) = possk(i-1) / kratio(i-1)
140  enddo
141  if (nout .gt. 0) then
142  nextout = nextout + 1
143  dumpout = .true.
144  endif
145  endif
146 
147 
148  if (time.lt.chktime .and. time + possk(1) .ge. chktime) then
149 c ## adjust time step to hit chktime exactly, and do checkpointing
150  possk(1) = chktime - time
151  do 13 i = 2, mxnest
152  13 possk(i) = possk(i-1) / kratio(i-1)
153  nextchk = nextchk + 1
154  dumpchk = .true.
155  else
156  dumpchk = .false.
157  endif
158 
159 c
160  level = 1
161  ntogo(level) = 1
162  do i = 1, maxlv
163  dtnew(i) = rinfinity
164  enddo
165 c
166 c ------------- regridding time? ---------
167 c
168 c check if either
169 c (i) this level should have its error estimated before being advanced
170 c (ii) this level needs to provide boundary values for either of
171 c next 2 finer levels to have their error estimated.
172 c this only affects two grid levels higher, occurs because
173 c previous time step needs boundary vals for giant step.
174 c no error estimation on finest possible grid level
175 c
176  60 continue
177  if (icheck(level) .ge. kcheck) then
178  lbase = level
179  else if (level+1 .ge. mxnest) then
180  go to 90
181  else if (icheck(level+1) .ge. kcheck) then
182  lbase = level+1
183  else if (level+2 .ge. mxnest) then
184  go to 90
185  else if (icheck(level+2) .ge. kcheck) then
186  lbase = level+2
187  else
188  go to 90
189  endif
190  if (lbase .eq. mxnest .or. lbase .gt. lfine) go to 70
191 c
192 c regrid level 'lbase+1' up to finest level.
193 c level 'lbase' stays fixed.
194 c
195  if (rprint) write(outunit,101) lbase
196 101 format(8h level ,i5,32h stays fixed during regridding )
197 
198  call system_clock(clock_start,clock_rate)
199  call cpu_time(cpu_start)
200  call regrid(nvar,lbase,cut,naux,start_time)
201  call system_clock(clock_finish,clock_rate)
202  call cpu_time(cpu_finish)
203  timeregridding = timeregridding + clock_finish - clock_start
204  timeregriddingcpu=timeregriddingcpu+cpu_finish-cpu_start
205 
206  call setbestsrc() ! need at every grid change
207 c call outtre(lstart(lbase+1),.false.,nvar,naux)
208 c note negative time to signal regridding output in plots
209 c call valout(lbase,lfine,-tlevel(lbase),nvar,naux)
210 c
211 c maybe finest level in existence has changed. reset counters.
212 c
213  if (rprint .and. lbase .lt. lfine) then
214  call outtre(lstart(lbase+1),.false.,nvar,naux)
215  endif
216  70 continue
217  do 80 i = lbase, lfine
218  80 icheck(i) = 0
219  do 81 i = lbase+1, lfine
220  81 tlevel(i) = tlevel(lbase)
221 c
222 c MJB: modified to check level where new grids start, which is lbase+1
223  if (verbosity_regrid.ge.lbase+1) then
224  do levnew = lbase+1,lfine
225  write(6,1006) intratx(levnew-1),intraty(levnew-1),
226  & kratio(levnew-1),levnew
227  1006 format(' Refinement ratios... in x:', i3,
228  & ' in y:',i3,' in t:',i3,' for level ',i4)
229  end do
230 
231  endif
232 
233 c ------- done regridding --------------------
234 c
235 c integrate all grids at level 'level'.
236 c
237  90 continue
238 
239 
240  call advanc(level,nvar,dtlevnew,vtime,naux)
241 
242 c # rjl modified 6/17/05 to print out *after* advanc and print cfl
243 c # rjl & mjb changed to cfl_level, 3/17/10
244 
245  timenew = tlevel(level)+possk(level)
246  if (tprint) then
247  write(outunit,100)level,cfl_level,possk(level),timenew
248  endif
249  if (method(4).ge.level) then
250  write(6,100)level,cfl_level,possk(level),timenew
251  endif
252 100 format(' AMRCLAW: level ',i2,' CFL = ',e10.3,
253  & ' dt = ',e11.4, ' final t = ',e13.6)
254 
255 
256 c # to debug individual grid updates...
257 c call valout(level,level,time,nvar,naux)
258 c
259 c done with a level of integration. update counts, decide who next.
260 c
261  ntogo(level) = ntogo(level) - 1
262  dtnew(level) = dmin1(dtnew(level),dtlevnew)
263  tlevel(level) = tlevel(level) + possk(level)
264  icheck(level) = icheck(level) + 1
265 c
266  if (level .lt. lfine) then
267  level = level + 1
268 c # check if should adjust finer grid time step to start wtih
269  if (((possk(level-1) - dtnew(level-1))/dtnew(level-1)) .gt.
270  . .05) then
271  dttemp = dtnew(level-1)/kratio(level-1)
272  ntogo(level) = (tlevel(level-1)-tlevel(level))/dttemp+.9
273  else
274  ntogo(level) = kratio(level-1)
275  endif
276  possk(level) = possk(level-1)/ntogo(level)
277  go to 60
278  endif
279 c
280  105 if (level .eq. 1) go to 110
281  if (ntogo(level) .gt. 0) then
282 c same level goes again. check for ok time step
283  106 if ((possk(level)-dtnew(level))/dtnew(level)
284  . .gt. .05) then
285 c adjust time steps for this and finer levels
286  ntogo(level) = ntogo(level) + 1
287  possk(level) = (tlevel(level-1)-tlevel(level))/
288  . ntogo(level)
289  go to 106
290  endif
291  go to 60
292  else
293  level = level - 1
294  call system_clock(clock_start,clock_rate)
295  call update(level,nvar,naux)
296  call system_clock(clock_finish,clock_rate)
297  timeupdating=timeupdating+clock_finish-clock_start
298  endif
299  go to 105
300 c
301 c --------------one complete coarse grid integration cycle done. -----
302 c
303 c time for output? done with the whole thing?
304 c
305  110 continue
306  time = time + possk(1)
307  ncycle = ncycle + 1
308  call conck(1,nvar,naux,time,rest)
309 
310  if ( vtime) then
311 c
312 c find new dt for next cycle (passed back from integration routine).
313  do 115 i = 2, lfine
314  ii = lfine+1-i
315  dtnew(ii) = min(dtnew(ii),dtnew(ii+1)*kratio(ii))
316  115 continue
317 c make sure not to exceed largest permissible dt
318  dtnew(1) = min(dtnew(1),dt_max)
319  possk(1) = dtnew(1) ! propagate new timestep to hierarchy
320  do 120 i = 2, mxnest
321  120 possk(i) = possk(i-1) / kratio(i-1)
322 
323  endif
324 
325  if ((abs(checkpt_style).eq.3 .and.
326  & mod(ncycle,checkpt_interval).eq.0) .or. dumpchk) then
327  call check(ncycle,time,nvar,naux)
328  dumpchk = .true.
329  endif
330 
331  if ((mod(ncycle,iout).eq.0) .or. dumpout) then
332  call valout(1,lfine,time,nvar,naux)
333  if (printout) call outtre(mstart,.true.,nvar,naux)
334  endif
335 
336  go to 20
337 c
338 999 continue
339 
340 c
341 c # computation is complete to final time or requested number of steps
342 c
343  if (ncycle .ge. nstop .and. tfinal .lt. rinfinity) then
344 c # warn the user that calculation finished prematurely
345  write(outunit,102) nstop
346  write(6,102) nstop
347  102 format('*** Computation halted after nv(1) = ',i8,
348  & ' steps on coarse grid')
349  endif
350 c
351 c # final output (unless we just did it above)
352 c
353  dump_final = ((iout.lt.iinfinity) .and. (mod(ncycle,iout).ne.0))
354  if (.not. dumpout) then
355  if (nout > 0) then
356  dump_final = (tout(nout).eq.tfinal)
357  endif
358  endif
359 
360  if (dump_final) then
361  call valout(1,lfine,time,nvar,naux)
362  if (printout) call outtre(mstart,.true.,nvar,naux)
363  endif
364 
365 c ## tick timing moved here so can be saved in checkpoint file
366 c
367  call system_clock(tick_clock_finish,tick_clock_rate)
368  call cpu_time(tick_cpu_finish)
369  timetick = timetick + tick_clock_finish - tick_clock_start
370  timetickcpu = timetickcpu + tick_cpu_finish - tick_cpu_start
371 
372 
373 c # checkpoint everything for possible future restart
374 c # (unless we just did it based on dumpchk)
375 c
376  if (checkpt_style .ne. 0) then ! want a chckpt
377  ! check if just did it so dont do it twice
378  if (.not. dumpchk) call check(ncycle,time,nvar,naux)
379  else ! no chkpt wanted, so need to print gauges separately
380  if (num_gauges .gt. 0) then
381  do ii = 1, num_gauges
383  end do
384  endif
385  endif
386 
387 
388 
389  write(6,*) "Done integrating to time ",time
390  return
391  end
integer timetick
Definition: amr_module.f90:242
integer, dimension(7) method
Definition: amr_module.f90:253
integer, dimension(maxlv) kratio
Definition: amr_module.f90:198
real(kind=8), parameter rinfinity
Definition: amr_module.f90:169
subroutine regrid(nvar, lbase, cut, naux, start_time)
Definition: regrid.f:5
subroutine print_gauges_and_reset_nextloc(gauge_num)
real(kind=8) tfinal
Definition: amr_module.f90:272
integer, parameter maxlv
Definition: amr_module.f90:174
real(kind=8), dimension(:), allocatable tout
Definition: amr_module.f90:271
logical rprint
Definition: amr_module.f90:297
subroutine conck(level, nvar, naux, time, rest)
Conservation check for specified level.
Definition: conck.f:7
integer checkpt_interval
Definition: amr_module.f90:280
subroutine setbestsrc()
subroutine valout(lst, lend, time, nvar, naux)
Definition: valout.f:5
subroutine check(nsteps, time, nvar, naux)
Definition: check.f:5
integer timeupdating
Definition: amr_module.f90:239
integer iout
Definition: amr_module.f90:270
real(kind=8) cfl_level
Definition: amr_module.f90:255
integer, dimension(maxlv) icheck
Definition: amr_module.f90:198
subroutine update(level, nvar, naux)
Synchronize between all grids on level level and grids on level level+1.
Definition: update.f:17
subroutine advanc(level, nvar, dtlevnew, vtime, naux)
Integrate all grids at the input level by one step of its delta(t)
Definition: advanc.f:11
integer, dimension(maxlv) lstart
Definition: amr_module.f90:198
real(kind=8) timetickcpu
Definition: amr_module.f90:243
integer, parameter iinfinity
Definition: amr_module.f90:170
logical printout
Definition: amr_module.f90:264
integer, dimension(maxlv) intraty
Definition: amr_module.f90:198
real(kind=8), dimension(:), allocatable tchk
Definition: amr_module.f90:281
real(kind=8) timeregriddingcpu
Definition: amr_module.f90:244
integer nstop
Definition: amr_module.f90:270
integer nchkpt
Definition: amr_module.f90:280
integer, parameter outunit
Definition: amr_module.f90:290
real(kind=8), dimension(maxlv) possk
Definition: amr_module.f90:193
integer num_gauges
logical tprint
Definition: amr_module.f90:297
integer checkpt_style
Definition: amr_module.f90:280
integer, dimension(maxlv) intratx
Definition: amr_module.f90:198
subroutine outtre(mlev, outgrd, nvar, naux)
Output a subtree of the grids.
Definition: outtre.f:11
subroutine tick(nvar, cut, nstart, vtime, time, naux, start_time, rest, dt_max)
Definition: tick.f:6
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 nout
Definition: amr_module.f90:270
integer mstart
Definition: amr_module.f90:198
integer kcheck
Definition: amr_module.f90:198
integer lfine
Definition: amr_module.f90:198
integer timeregridding
Definition: amr_module.f90:239