2D AMRCLAW
restrt_hdf.f
Go to the documentation of this file.
1 c
2 c ---------------------------------------------------------
3 c
4  subroutine restrt(nsteps,time,nvar)
5 c
6  use amr_module
7  implicit double precision (a-h,o-z)
8 
9  logical ee
10 
11 
12  dimension intrtx(maxlv),intrty(maxlv),intrtt(maxlv)
13  character*16 chkname
14 c
15 c # HDF: Declare variables that describe datasets and HDF files.
16 c
17  integer sd_id
18 c
19 c # HDF: Arrays which will hold output arrays.
20 c
21  dimension iqout(15), qout(4)
22 c
23 c # HDF: Declare external HDF functions
24 c
25  integer sfstart, sfend
26  external sfstart, sfend
27 c
28 c # HDF: Set up HDF constants
29 c
30  integer DFACC_READ, DFACC_WRITE, DFACC_CREATE
31  parameter(dfacc_read = 1, dfacc_write = 2, dfacc_create = 4)
32 
33  integer SUCCEED, FAIL
34  parameter(succeed = 0, fail = -1)
35 c
36 c :::::::::::::::::::::::::::: RESTRT ::::::::::::::::::::::::::::::::
37 c read back in the check point files written by subr. check.
38 c
39 c some input variables might have changed, and also the
40 c alloc array could have been written with a smaller size at checkpoint
41 c ::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
42 c
43 c ### make the file name showing the time step
44 c
45  chkname = 'restart.data.hdf'
46 c
47 c # HDF: open hdf restart file.
48 c
49  sd_id = sfstart(chkname,dfacc_read)
50  if (sd_id.eq.fail) THEN
51  WRITE(*,*) 'Failed to open HDF file',
52  & ' (call to sfstart in restrt_hdf.f)'
53  stop
54  end if
55 
56  index = 0
57  call read_integer_vector(sd_id,15,index,'iqout ',iqout)
58  lenmax = iqout(1)
59  lendim = iqout(2)
60  isize = iqout(3)
61  lenf = iqout(4)
62  ibuf = iqout(5)
63  mstart = iqout(6)
64  ndfree = iqout(7)
65  lfine = iqout(8)
66  iorder = iqout(9)
67  mxnold = iqout(10)
68  kcheck1 = iqout(11)
69  nsteps = iqout(12)
70  matlabu = iqout(13)
71  lentot = iqout(14)
72  ndfree_bnd = iqout(15)
73 
74  index = index + 1
75  call read_integer_vector(sd_id,maxlv,index,'icheck',icheck)
76  index = index + 1
77  call read_integer_vector(sd_id,maxlv,index,'lstart',lstart)
78  index = index + 1
79  call read_integer_vector(sd_id,maxlv,index,'newstl',newstl)
80  index = index + 1
81  call read_integer_vector(sd_id,maxlv,index,'listsp',listsp)
82  index = index + 1
83  call read_integer_vector(sd_id,maxlv,index,'intratx',intrtx)
84  index = index + 1
85  call read_integer_vector(sd_id,maxlv,index,'intraty',intrty)
86  index = index + 1
87  call read_integer_vector(sd_id,maxlv,index,'kratio',intrtt)
88  index = index + 1
89  call read_integer_vector(sd_id,maxlv,index,'iregsz',iregsz)
90  index = index + 1
91  call read_integer_vector(sd_id,maxlv,index,'jregsz',jregsz)
92  index = index + 1
93  call read_integer_array(sd_id,lfdim,2,index,'lfree ',lfree)
94  index = index + 1
95  call read_integer_array(sd_id,rsize,maxgr,index,'node ',node)
96 
97 c # HDF: DFNT_FLOAT64
98 
99  index = index + 1
100  call read_double_vector(sd_id,4,index,'qout ',qout)
101  tl = qout(1)
102  time = qout(2)
103  evol = qout(3)
104  rvol = qout(4)
105 
106  index = index + 1
107  call read_double_vector(sd_id,maxlv,index,'hxposs',hxposs)
108  index = index + 1
109  call read_double_vector(sd_id,maxlv,index,'hyposs',hyposs)
110  index = index + 1
111  call read_double_vector(sd_id,maxlv,index,'possk ',possk)
112  index = index + 1
113  call read_double_vector(sd_id,10,index,'rvoll ',rvoll)
114  index = index + 1
115  call read_double_vector(sd_id,lendim,index,'alloc ',alloc)
116 
117  index = index + 1
118  call read_double_array(sd_id,rsize,maxgr,index,'rnode ',rnode)
119 c
120 c # HDF: Close HDF file.
121 c
122  istat = sfend(sd_id)
123  if (istat.eq.fail) then
124  WRITE(*,*) 'Failed to close SDS',
125  & ' (call to sfend in restrt_hdf.f)'
126  stop
127  end if
128 
129  write(outunit,100) nsteps,time
130  write(6,100) nsteps,time
131  100 format(/,' RESTARTING the calculation after ',i5,' steps',
132  1 /,' (time = ',e15.7,')')
133 
134  do i = 1, mxnold-1
135  if ( (intratx(i) .ne. intrtx(i)) .or.
136  . (intraty(i) .ne. intrty(i)) .or.
137  . (kratio(i) .ne. intrtt(i)) ) then
138  write(outunit,*)
139  . " not allowed to change existing refinement ratios on Restart"
140  write(*,*)
141  . " not allowed to change existing refinement ratios on Restart"
142  write(outunit,*)" Old ratios:"
143  write(*,*) " Old ratios:"
144  write(outunit,903)(intrtx(j),j=1,mxnold-1)
145  write(*,903) (intrtx(j),j=1,mxnold-1)
146  write(outunit,903)(intrty(j),j=1,mxnold-1)
147  write(*,903) (intrty(j),j=1,mxnold-1)
148  write(outunit,903)(intrtt(j),j=1,mxnold-1)
149  write(*,903) (intrtt(j),j=1,mxnold-1)
150  903 format(6i3)
151  stop
152  endif
153  end do
154 c
155 c adjust free list of storage in case size has changed.
156 c
157  idif = memsize - isize
158  if (idif .gt. 0) then
159  lfree(lenf,1) = isize + 2
160  call reclam(isize+1,idif)
161  else if (idif .lt. 0) then
162  write(outunit,900) isize, memsize
163  write(*,900) isize, memsize
164  900 format(' size of alloc not allowed to shrink with restart ',/,
165  . ' old size ',i7,' current size ',i7)
166  stop
167  endif
168 c
169 c adjust storage in case mxnest has changed - only allow it to increase,
170 c and only at non-multiples of error estimation on old mxnest.
171 c
172  if (mxnest .eq. mxnold) go to 99
173 
174  if (mxnest .lt. mxnold) then
175  if (lfine .lt. mxnest) then
176  go to 99
177  else
178  write(outunit,901) mxnold, mxnest
179  write(*, 901) mxnold, mxnest
180  901 format(' only allow mxnest to increase: ',/,
181  & ' old mxnest ',i4, ' new mxnest ',i4)
182  stop
183  endif
184  endif
185 
186 c see if simple enough situation to allow changing mxnest
187  ee = .false.
188  do 10 level = 1, mxnold
189  if (icheck(level) .ge. kcheck) then
190  ee = .true.
191  kmust = icheck(level)
192  endif
193  10 continue
194  if (ee) then
195  write(outunit,902) mxnold, mxnest, kmust
196  write(* ,902) mxnold, mxnest, kmust
197  902 format(/,' only allow changes in mxnest (from ',
198  & i4,' to ',i4,')',/,
199  & ' when not time to error estimate: ',/,
200  & ' please run a few more steps before changing ',/,
201  & ' so that # of steps not greater then kcheck',/,
202  & ' or make kcheck > ',i4 )
203  stop
204  else
205 c # add second storage location to previous mxnest level
206  mptr = lstart(mxnold)
207  15 if (mptr .eq. 0) go to 25
208  mitot = node(ndihi,mptr)-node(ndilo,mptr)+1+2*nghost
209  mjtot = node(ndjhi,mptr)-node(ndjlo,mptr)+1+2*nghost
210  node(store2,mptr) = igetsp(mitot*mjtot*nvar)
211  mptr = node(levelptr,mptr)
212  go to 15
213  25 continue
214  endif
215 c
216 c # add new info. to spatial and counting arrays
217  99 level = lfine + 1
218  35 if (level .gt. mxnest) go to 45
219  hxposs(level) = hxposs(level-1) / dble(intratx(level-1))
220  hyposs(level) = hyposs(level-1) / dble(intraty(level-1))
221  possk(level) = possk(level-1) / dble(kratio(level-1))
222  iregsz(level) = iregsz(level-1) * dble(intratx(level-1))
223  jregsz(level) = jregsz(level-1) * dble(intraty(level-1))
224  level = level + 1
225  go to 35
226  45 continue
227 c
228 c
229  return
230  end
231 c=================================================================
232 c # HDF: NOTE THAT ALL INTEGER HDF CHECKPOINTING ROUTINES FOR
233 c # BOTH WRITING AND READING ARE INCLUDED HERE, ALTHOUGH
234 c # THE INTEGER WRITING ROUTINES ARE CALLED IN check_hdf.f.
235 c # THIS IS DONE TO PREVENT COMPILER ERRORS SINCE THE SAME
236 c # HDF ROUTINES ARE USED TO READ BOTH INTEGER AND DOUBLE
237 c # PRECISION ARRAYS.
238 c=================================================================
239  subroutine read_integer_vector(sd_id,idims,index,qname,iout)
240  implicit double precision (a-h,o-z)
241  character*6 qname
242 c
243 c # HDF: Declare variables that describe datasets and HDF files.
244 c
245  integer sd_id, sds_id
246  dimension iout(idims), istart(1), istride(1), iedges(1)
247 c
248 c # HDF: Declare external HDF functions
249 c
250  integer sfcreate, sfrdata, sfselect, sfendacc
251  external sfcreate, sfrdata, sfselect, sfendacc
252 c
253 c # HDF: Set up HDF constants
254 c
255  integer SUCCEED, FAIL
256  parameter(succeed = 0, fail = -1)
257 c
258 c # HDF: Select dataset in HDF file.
259 c
260  sds_id = sfselect(sd_id,index)
261  if (sds_id.eq.fail) THEN
262  WRITE(*,*) 'Failed to select data set for variable ', qname,
263  & ' in restart HDF file'
264  WRITE(*,*) '(call to sfselect in restrt_hdf.f)'
265  stop
266  end if
267 c
268 c # HDF: Set up dimensions for integer vector.
269 c
270  istart(1) = 0
271  iedges(1) = idims
272  istride(1) = 1
273 c
274 c # HDF: read integer vector from hdf file.
275 c
276  istat = sfrdata(sds_id,istart,istride,iedges,iout)
277  if (istat.eq.fail) THEN
278  WRITE(*,*) 'Failed to read variable ', qname,
279  & ' from restart HDF file'
280  WRITE(*,*) '(call to sfrdata in restrt_hdf.f)'
281  stop
282  end if
283 c
284 c # HDF: End access to integer vector in HDF file.
285 c
286  istat = sfendacc(sds_id)
287  if (istat.eq.fail) THEN
288  WRITE(*,*) 'Failed to end access to variable ', qname,
289  & ' in restart HDF file'
290  WRITE(*,*) '(call to sfendacc in restrt_hdf.f)'
291  stop
292  end if
293 
294  return
295  end
296 
297 c==========================================================
298  subroutine read_integer_array(sd_id,idim1,idim2,index,qname,iout)
299  implicit double precision (a-h,o-z)
300  character*6 qname
301 c
302 c # HDF: Declare variables that describe datasets and HDF files.
303 c
304  integer sd_id, sds_id
305  dimension istart(2), istride(2), iedges(2)
306  dimension iout(idim1,idim2)
307 c
308 c # HDF: Declare external HDF functions
309 c
310  integer sfcreate, sfrdata, sfselect, sfendacc
311  external sfcreate, sfrdata, sfselect, sfendacc
312 c
313 c # HDF: Set up HDF constants
314 c
315  integer SUCCEED, FAIL
316  parameter(succeed = 0, fail = -1)
317 c
318 c # HDF: Select dataset in HDF file.
319 c
320  sds_id = sfselect(sd_id,index)
321  if (sds_id.eq.fail) THEN
322  WRITE(*,*) 'Failed to select data set for variable ', qname,
323  & ' in restart HDF file'
324  WRITE(*,*) '(call to sfselect in restrt_hdf.f)'
325  stop
326  end if
327 c
328 c # HDF: Set up dimensions for integer array.
329 c
330  istart(1) = 0
331  istart(2) = 0
332  iedges(1) = idim1
333  iedges(2) = idim2
334  istride(1) = 1
335  istride(2) = 1
336 c
337 c # HDF: read integer array from hdf file.
338 c
339  istat = sfrdata(sds_id,istart,istride,iedges,iout)
340  if (istat.eq.fail) THEN
341  WRITE(*,*) 'Failed to read variable ', qname,
342  & ' from restart HDF file'
343  WRITE(*,*) '(call to sfrdata in restrt_hdf.f)'
344  stop
345  end if
346 c
347 c # HDF: End access to integer array in HDF file.
348 c
349  istat = sfendacc(sds_id)
350  if (istat.eq.fail) THEN
351  WRITE(*,*) 'Failed to end access to variable ', qname,
352  & ' in restart HDF file'
353  WRITE(*,*) '(call to sfendacc in restrt_hdf.f)'
354  stop
355  end if
356 
357  return
358  end
359 
360 c==========================================================
361  subroutine dump_integer_vector(sd_id,idims,qname,iout)
362  implicit double precision (a-h,o-z)
363  character*6 qname
364 c
365 c # HDF: Declare variables that describe datasets and HDF files.
366 c
367  integer sd_id, sds_id
368  dimension iout(idims), istart(1), istride(1), iedges(1), idim(1)
369 c
370 c # HDF: Declare external HDF functions
371 c
372  integer sfcreate, sfwdata, sfscompress, sfendacc
373  external sfcreate, sfwdata, sfscompress, sfendacc
374 c
375 c # HDF: Set up HDF constants
376 c
377  integer DFNT_FLOAT64, DFNT_INT32
378  parameter(dfnt_float64 = 6, dfnt_int32 = 24)
379 
380  integer SUCCEED, FAIL
381  parameter(succeed = 0, fail = -1)
382 c
383 c # HDF: Set up compression constants for HDF file.
384 c
385  integer COMP_CODE_DEFLATE, DEFLATE_LEVEL
386  parameter(comp_code_deflate = 4, deflate_level = 6)
387 c
388 c # HDF: Create integer vector in HDF file.
389 c # NOTE THAT THE RANK MUST BE ONE.
390 c
391  irank = 1
392  idim(1) = idims
393  sds_id = sfcreate(sd_id,qname,dfnt_int32,irank,idim)
394  if (sds_id.eq.fail) THEN
395  WRITE(*,*) 'Failed to create variable ', qname,
396  & ' in restart HDF file'
397  WRITE(*,*) '(call to sfcreate in check_hdf.f)'
398  stop
399  end if
400 c
401 c # HDF: Set up dimensions for integer vector.
402 c
403  istart(1) = 0
404  iedges(1) = idims
405  istride(1) = 1
406 c
407 c # HDF: set compression mode and write data to hdf file.
408 c
409  istat=sfscompress(sds_id,comp_code_deflate,deflate_level)
410  istat = sfwdata(sds_id,istart,istride,iedges,iout)
411  if (istat.eq.fail) THEN
412  WRITE(*,*) 'Failed to write variable ', qname,
413  & ' in restart HDF file'
414  WRITE(*,*) '(call to sfwdata in check_hdf.f)'
415  stop
416  end if
417 c
418 c # HDF: End access to integer vector in HDF file.
419 c
420  istat = sfendacc(sds_id)
421  if (istat.eq.fail) THEN
422  WRITE(*,*) 'Failed to end access to variable ', qname,
423  & ' in restart HDF file'
424  WRITE(*,*) '(call to sfendacc in check_hdf.f)'
425  stop
426  end if
427 
428  return
429  end
430 
431 c==========================================================
432  subroutine dump_integer_array(sd_id,idim1,idim2,qname,iout)
433  implicit double precision (a-h,o-z)
434  character*6 qname
435 c
436 c # HDF: Declare variables that describe datasets and HDF files.
437 c
438  integer sd_id, sds_id
439  dimension idims(2), istart(2), istride(2), iedges(2)
440  dimension iout(idim1,idim2)
441 c
442 c # HDF: Declare external HDF functions
443 c
444  integer sfcreate, sfwdata, sfscompress, sfendacc
445  external sfcreate, sfwdata, sfscompress, sfendacc
446 c
447 c # HDF: Set up HDF constants
448 c
449  integer DFNT_FLOAT64, DFNT_INT32
450  parameter(dfnt_float64 = 6, dfnt_int32 = 24)
451 
452  integer SUCCEED, FAIL
453  parameter(succeed = 0, fail = -1)
454 c
455 c # HDF: Set up compression constants for HDF file.
456 c
457  integer COMP_CODE_DEFLATE, DEFLATE_LEVEL
458  parameter(comp_code_deflate = 4, deflate_level = 6)
459 c
460 c # HDF: Create integer array in HDF file.
461 c # NOTE THAT THE RANK MUST BE TWO.
462 c
463  irank = 2
464  idims(1) = idim1
465  idims(2) = idim2
466  sds_id = sfcreate(sd_id,qname,dfnt_int32,irank,idims)
467  if (sds_id.eq.fail) THEN
468  WRITE(*,*) 'Failed to create variable ', qname,
469  & ' in restart HDF file'
470  WRITE(*,*) '(call to sfcreate in check_hdf.f)'
471  stop
472  end if
473 c
474 c # HDF: Set up dimensions for integer array.
475 c
476  istart(1) = 0
477  istart(2) = 0
478  iedges(1) = idims(1)
479  iedges(2) = idims(2)
480  istride(1) = 1
481  istride(2) = 1
482 c
483 c # HDF: set compression mode and write data to hdf file.
484 c
485  istat=sfscompress(sds_id,comp_code_deflate,deflate_level)
486  istat = sfwdata(sds_id,istart,istride,iedges,iout)
487  if (istat.eq.fail) THEN
488  WRITE(*,*) 'Failed to write variable ', qname,
489  & ' in restart HDF file'
490  WRITE(*,*) '(call to sfwdata in check_hdf.f)'
491  stop
492  end if
493 c
494 c # HDF: End access to integer array in HDF file.
495 c
496  istat = sfendacc(sds_id)
497  if (istat.eq.fail) THEN
498  WRITE(*,*) 'Failed to end access to variable ', qname,
499  & ' in restart HDF file'
500  WRITE(*,*) '(call to sfendacc in check_hdf.f)'
501  stop
502  end if
503 
504  return
505  end
506 
507 
integer, dimension(maxlv) kratio
Definition: amr_module.f90:198
integer, parameter rsize
Definition: amr_module.f90:30
integer lendim
Definition: amr_module.f90:247
subroutine read_double_array(sd_id, idim1, idim2, index, qname, out)
Definition: check_hdf.f:326
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
integer, parameter maxgr
Definition: amr_module.f90:173
real(kind=8) evol
Definition: amr_module.f90:237
integer lenmax
Definition: amr_module.f90:247
subroutine restrt(nsteps, time, nvar, naux)
Definition: restrt.f:6
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, parameter lfdim
Definition: amr_module.f90:224
integer memsize
Definition: amr_module.f90:219
real(kind=8), dimension(rsize, maxgr) rnode
Definition: amr_module.f90:193
integer ndfree_bnd
Definition: amr_module.f90:198
integer lentot
Definition: amr_module.f90:247
subroutine dump_integer_array(sd_id, idim1, idim2, qname, iout)
Definition: restrt_hdf.f:433
integer, dimension(maxlv) icheck
Definition: amr_module.f90:198
integer matlabu
Definition: amr_module.f90:283
subroutine read_integer_vector(sd_id, idims, index, qname, iout)
Definition: restrt_hdf.f:240
integer iorder
Definition: amr_module.f90:198
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
real(kind=8) rvol
Definition: amr_module.f90:237
subroutine read_double_vector(sd_id, idims, index, qname, out)
Definition: check_hdf.f:267
integer, parameter store2
pointer to the address of memory storing the second copy of solution data on this grid...
Definition: amr_module.f90:105
subroutine read_integer_array(sd_id, idim1, idim2, index, qname, iout)
Definition: restrt_hdf.f:299
integer, dimension(maxlv) intraty
Definition: amr_module.f90:198
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, parameter ndjhi
global j index of upper border of this grid
Definition: amr_module.f90:117
integer, dimension(maxlv) listsp
Definition: amr_module.f90:198
integer lenf
Definition: amr_module.f90:225
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
subroutine dump_integer_vector(sd_id, idims, qname, iout)
Definition: restrt_hdf.f:362
integer mstart
Definition: amr_module.f90:198
integer kcheck
Definition: amr_module.f90:198
integer lfine
Definition: amr_module.f90:198
real(kind=8), dimension(:), allocatable alloc
Definition: amr_module.f90:218