source: CIVL/mods/dev.civl.abc/examples/fortran/flash/heat_min/io_readData.F90

main
Last change on this file was aad342c, checked in by Stephen Siegel <siegel@…>, 3 years ago

Performing huge refactor to incorporate ABC, GMC, and SARL into CIVL repo and use Java modules.

git-svn-id: svn://vsl.cis.udel.edu/civl/trunk@5664 fb995dde-84ed-4084-dfe6-e5aef3e2452c

  • Property mode set to 100644
File size: 13.9 KB
Line 
1!!****if* source/IO/IOMain/hdf5/parallel/PM/io_readData
2!!
3!! NAME
4!!
5!! io_readData
6!!
7!!
8!! SYNOPSIS
9!!
10!! io_readData()
11!!
12!!
13!! DESCRIPTION
14!!
15!! This is the reading counterpart to io_writeData. It reads an HDF5
16!! file and distributes it to the processors to restart a simulation.
17!!
18!! All reading is done using parallel HDF5 -- no explicit data movement
19!! is used.
20!!
21!!
22!! ARGUMENTS
23!!
24!!***
25
26!!REORDER(5): unk, facevar[xyz], unkBuf, unkBufGC, face[XYZ]Buf
27
28#ifdef DEBUG_ALL
29#define DEBUG_IO
30#endif
31
32#include "constants.h"
33#include "Flash.h"
34#include "io_flash.h"
35
36subroutine io_readData()
37
38
39
40 use Grid_data, ONLY : gr_globalNumBlocks, gr_nToLeft, &
41 gr_globalOffset, gr_gid
42 use Simulation_interface, ONLY : Simulation_mapStrToInt
43 use Driver_interface, ONLY : Driver_abortFlash
44 use Logfile_interface, ONLY : Logfile_stamp
45 use Grid_interface, ONLY : Grid_putLocalNumBlks, Grid_getBlkIndexLimits, &
46 Grid_receiveInputData
47
48 use IO_data, ONLY : io_globalMe, io_globalComm,&
49 io_baseName, io_checkpointFileNumber, &
50 io_realParmNames, io_realParmValues, io_numRealParms, &
51 io_intParmNames, io_intParmValues, io_numIntParms, &
52 io_logParmNames, io_logParmValues, io_numLogParms, &
53 io_strParmNames, io_strParmValues, io_numStrParms, &
54 io_realScalarNames, io_realScalarValues, io_numRealScalars, &
55 io_intScalarNames, io_intScalarValues, io_numIntScalars, &
56 io_logScalarNames, io_logScalarValues, io_numLogScalars, &
57 io_strScalarNames, io_strScalarValues, io_numStrScalars, &
58 io_logToIntScalarValues, io_logToIntParmValues, io_unklabels, &
59 io_ilo, io_ihi, io_jlo, io_jhi, io_klo, io_khi, &
60 io_outputSplitNum, io_comm, io_chkptFileID, &
61 io_chkGuardCellsInput, io_chkGuardCellsOutput, &
62 io_faceXVarLabels, io_faceYVarLabels, io_faceZVarLabels,&
63 io_realParmNamesPrev, io_realParmValuesPrev, io_numRealParmsPrev, &
64 io_intParmNamesPrev, io_intParmValuesPrev, io_numIntParmsPrev, &
65 io_logParmNamesPrev, io_logParmValuesPrev, io_numLogParmsPrev, &
66 io_strParmNamesPrev, io_strParmValuesPrev, io_numStrParmsPrev, &
67 io_logToIntParmValuesPrev, io_splitNumBlks, &
68 io_unklabelsGlobal, io_meshMe, io_meshNumProcs, tree_data_t
69 use IO_interface, ONLY : IO_getScalar
70
71#ifdef FLASH_GRID_PARAMESH
72 use tree, ONLY : maxblocks_tr, nfaces, nchild, nodetype, &
73 lrefine, bnd_box, coord, bsize, neigh, parent, child
74#ifdef FLASH_GRID_PARAMESH3OR4
75 use tree, ONLY : MFLAGS, which_child, bflags
76 use Grid_data, ONLY : gr_gsurr_blks
77#endif
78#endif
79 use io_typeInterface, ONLY : io_xfer_tree_data
80 use physicaldata, ONLY : unk, facevarx, facevary, facevarz
81
82
83 implicit none
84 type(tree_data_t) :: tree_data
85
86#include "Flash_mpi.h"
87
88
89
90
91
92 integer :: localNumBlocks, ngid
93 integer :: alocalNumBlocks
94
95
96 character (len=4) :: fnumStr
97 character (len=MAX_STRING_LENGTH) :: filename
98 character(len=MAX_STRING_LENGTH), save, allocatable, dimension(:,:) :: strBuff
99
100
101 integer :: blockID, procBlocks, ierr
102
103 integer :: i, lb, j, xx, yy, alnblocks, u
104
105 integer, allocatable :: procnumber(:) !dimension(localNumBlocks)
106
107 integer :: realGlobalNumBlocks
108! real :: unkBufGC(1,GRID_IHI_GC,GRID_JHI_GC,GRID_KHI_GC,MAXBLOCKS)
109 real,allocatable :: unkBufGC(:,:,:,:,:)
110 real,allocatable :: unkBuf(:,:,:,:,:)
111 real,allocatable :: faceXBuf(:,:,:,:,:)
112 real,allocatable :: faceYBuf(:,:,:,:,:)
113 real,allocatable :: faceZBuf(:,:,:,:,:)
114
115 integer :: blkLimits(2,MDIM), blkLimitsGC(2, MDIM)
116
117 integer :: splitOffset, localOffset
118 integer :: doread
119
120 integer :: presentDims
121 integer, parameter :: xferType = IO_READ_XFER, libType = IO_FILE_HDF5
122
123
124 call io_getOutputName(io_checkpointFileNumber, "hdf5", "_chk_", filename, .false.)
125
126
127 call io_h5open_file_for_read(io_chkptFileID, filename, io_comm, io_outputSplitNum)
128
129 if (io_globalMe == MASTER_PE) then
130 allocate (strBuff(2,2))
131 print *, 'file: ', trim(filename), ' opened for restart'
132 write (strBuff(1,1), "(A)") "type"
133 write (strBuff(1,2), "(A)") "checkpoint"
134 write (strBuff(2,1), "(A)") "name"
135 write (strBuff(2,2), "(A)") trim(filename)
136 call Logfile_stamp( strBuff, 2, 2, "[io_readData]")
137 end if
138 if (allocated(strBuff)) deallocate(strBuff)
139
140
141 !!read in the total number of blocks, time, and timestep
142 call io_h5read_header(io_globalMe, &
143 io_chkptFileID, &
144 io_unklabels, &
145 io_outputSplitNum)
146
147
148
149 call io_prepareListsRead()
150
151
152
153 call io_h5read_runtime_parameters(io_chkptFileID, &
154 io_numRealParmsPrev, &
155 io_realParmNamesPrev, &
156 io_realParmValuesPrev, &
157 io_numIntParmsPrev, &
158 io_intParmNamesPrev, &
159 io_intParmValuesPrev, &
160 io_numStrParmsPrev, &
161 io_strParmNamesPrev, &
162 io_strParmValuesPrev, &
163 io_numLogParmsPrev, &
164 io_logParmNamesPrev, &
165 io_logToIntParmValuesPrev)
166
167
168
169 call io_h5read_scalars(io_chkptFileID, &
170 io_numRealScalars, &
171 io_realScalarNames, &
172 io_realScalarValues, &
173 io_numIntScalars, &
174 io_intScalarNames, &
175 io_intScalarValues, &
176 io_numStrScalars, &
177 io_strScalarNames, &
178 io_strScalarValues, &
179 io_numLogScalars, &
180 io_logScalarNames, &
181 io_logToIntScalarValues)
182
183
184 call io_finalizeListsRead()
185
186
187 call IO_getScalar("globalNumBlocks", gr_globalNumBlocks)
188 call io_checkBlockShape(gr_globalNumBlocks)
189 if (io_outputSplitNum /= 1) then
190 call IO_getScalar("splitNumBlocks", io_splitNumBlks)
191 else
192 io_splitNumBlks = gr_globalNumBlocks
193 endif
194
195
196 !---------------------------------------------------------------------------
197 ! compute the number of blocks on each processor -- this will be used to
198 ! get the offset into the file for the parallel read
199 !---------------------------------------------------------------------------
200
201 ! compute the approximate number of blocks per processor
202 alnblocks = int(gr_globalNumBlocks/io_meshNumProcs) + 1
203
204 ! check for error -- if the number of blocks we want to put on each
205 ! processor is greater than maxblocks, then abort
206 if (alnblocks .GT. MAXBLOCKS) then
207
208 print *
209 print *, '********** ERROR in READ_DATA ************'
210 print *
211 print *,' Number of blocks per processor exceeds maxblocks.'
212 print *,' Suggest you reset maxblocks to a larger number or'
213 print *,' run on a larger number of processors.'
214 print *,' globalNumBlocks, io_meshNumProcs = ', gr_globalNumBlocks, io_meshNumProcs
215 print *
216
217 call Driver_abortFlash('[io_readData] ERROR: num blocks per proc exceeds maxblocks')
218
219 end if
220
221 ! figure out the excess blocks
222 yy = (io_meshNumProcs*alnblocks) - gr_globalNumBlocks
223 xx = io_meshNumProcs - yy
224
225 ! loop over all the processor numbers and figure out how many blocks are
226 ! stored to the left of the processor -- this is a little tricky
227
228 gr_nToLeft(0) = 0
229
230 do i = 0, io_meshNumProcs - 2
231 if (i .LT. xx) then
232 procBlocks = alnblocks
233 else
234 procBlocks = alnblocks - 1
235 endif
236
237 if (alnblocks .EQ. 0) then
238 if (i .LT. gr_globalNumBlocks) then
239 procBlocks = 1
240 else
241 procBlocks = 0
242 end if
243 end if
244
245 ! we have the number of blocks on proc i, the number of blocks on i+1 is
246 ! the number of blocks on i + the number of blocks left of i
247 if (i .EQ. 0) then
248 gr_nToLeft(i+1) = procBlocks
249 else
250 gr_nToLeft(i+1) = procBlocks + gr_nToLeft(i)
251 endif
252 enddo
253
254 ! figure out how many blocks are on the current proc.
255 if (io_meshMe < xx) then
256 localNumBlocks = alnblocks
257 else
258 localNumBlocks = alnblocks - 1
259 endif
260
261 if (alnblocks .EQ. 0) then
262 if (io_meshMe < gr_globalNumBlocks) then
263 localNumBlocks = 1
264 else
265 localNumBlocks = 0
266 end if
267 end if
268
269 ! compute the offset into the dataspace in the HDF5 file
270 gr_globalOffset = gr_nToLeft(io_meshMe)
271
272
273
274 call Grid_putLocalNumBlks(localNumBlocks)
275
276 !find our offset into a potentially split file:
277 if(io_outputSplitNum > 1) then
278 call MPI_ALLREDUCE(gr_globalOffset, splitOffset,1,FLASH_INTEGER,&
279 MPI_MIN, io_comm, ierr)
280 localOffset = gr_globalOffset - splitOffset
281
282 else
283 localOffset = gr_globalOffset
284 end if
285
286
287 call io_h5read_present_dims(io_chkptFileID, presentDims)
288
289
290
291 tree_data % bnd_box => bnd_box
292 tree_data % coord => coord
293 tree_data % bsize => bsize
294 tree_data % gid => gr_gid
295 tree_data % nodetype => nodetype
296 tree_data % lrefine => lrefine
297# ifdef FLASH_GRID_PARAMESH3OR4
298 tree_data % bflags => bflags
299 tree_data % which_child => which_child
300# else
301 nullify(tree_data % bflags)
302 nullify(tree_data % which_child)
303# endif
304 nullify(tree_data % procnumber)
305# ifdef FLASH_GRID_PARAMESH4DEV_SURR_BLKS_OPTION
306 tree_data % gsurr_blks => gr_gsurr_blks
307# else
308 nullify(tree_data % gsurr_blks)
309# endif
310
311 call io_h5read_present_dims(io_chkptFileID, presentDims)
312
313 call io_xfer_tree_data(tree_data, &
314 io_chkptFileID, libType, xferType, &
315 localNumBlocks, localOffset, presentDims)
316
317 !Extract data from gr_gid and gr_gsurr_blks arrays
318 call Grid_receiveInputData(localNumBlocks, alnblocks, xx)
319
320
321
322
323 allocate(unkBuf(1, NXB, NYB, NZB, MAXBLOCKS))
324
325 if (io_chkGuardCellsInput) then ! currently, blkLimits is not used otherwise
326 blockID = 1
327 call Grid_getBlkIndexLimits(blockID, blkLimits, blkLimitsGC)
328#ifndef FL_NON_PERMANENT_GUARDCELLS
329 blkLimits = blkLimitsGC
330#endif
331 end if
332
333 do u=1, ubound(io_unklabelsGlobal,1)
334 call Simulation_mapStrToInt(io_unklabelsGlobal(u),i,MAPBLOCK_UNK)
335 doread = 0
336 if(i /= NONEXISTENT) doread = 1
337
338 if(io_chkGuardCellsInput) then
339 allocate(unkBufGC(1,GRID_IHI_GC,GRID_JHI_GC,GRID_KHI_GC,MAXBLOCKS))
340 unkBufGC(1,:GRID_IHI_GC,:GRID_JHI_GC,:GRID_KHI_GC,1:localNumBlocks) = 0.0
341 call io_h5read_unknowns(io_chkptFileID, &
342 GRID_IHI_GC, &
343 GRID_JHI_GC, &
344 GRID_KHI_GC, &
345 unkBufGC, &
346 io_unklabelsGlobal(u), &
347 localNumBlocks, &
348 io_splitNumBlks, &
349 localOffset, &
350 doread)
351
352 if(i /= NONEXISTENT) then
353 unk(i,:,:,:,1:MAXBLOCKS) = &
354 unkBufGC(1,blkLimits(LOW,IAXIS):blkLimits(HIGH,IAXIS), &
355 blkLimits(LOW,JAXIS):blkLimits(HIGH,JAXIS), &
356 blkLimits(LOW,KAXIS):blkLimits(HIGH,KAXIS), 1:MAXBLOCKS)
357 end if
358 deallocate(unkBufGC)
359 else
360 unkBuf(1,:,:,:,1:localNumBlocks) = 0.0
361 call io_h5read_unknowns(io_chkptFileID, &
362 NXB, &
363 NYB, &
364 NZB, &
365 unkBuf, &
366 io_unklabelsGlobal(u), &
367 localNumBlocks, &
368 io_splitNumBlks, &
369 localOffset, &
370 doread)
371
372 if(i /= NONEXISTENT) then
373 unk(i,io_ilo:io_ihi, io_jlo:io_jhi, io_klo:io_khi,1:MAXBLOCKS) = &
374 unkBuf(1,1:NXB,1:NYB,1:NZB,1:MAXBLOCKS)
375 end if
376 end if
377 enddo
378
379 deallocate(unkBuf)
380#ifdef FLASH_GRID_PARAMESH3OR4
381#if(NFACE_VARS>0)
382
383
384 allocate(faceXBuf(1,1:NXB+1,1:NYB,1:NZB,1:MAXBLOCKS))
385 if (NDIM .gt. 1) allocate(faceYBuf(1,1:NXB,1:NYB+1,1:NZB,1:MAXBLOCKS))
386 if (NDIM .gt. 2) allocate(faceZBuf(1,1:NXB,1:NYB,1:NZB+1,1:MAXBLOCKS))
387
388 doread = 1
389 do i = 1,NFACE_VARS
390
391 call io_h5read_unknowns(io_chkptFileID, &
392 NXB+1, &
393 NYB, &
394 NZB, &
395 faceXBuf, &
396 io_faceXVarLabels(i), &
397 localNumBlocks, &
398 io_splitNumBlks, &
399 localOffset, &
400 doread)
401
402 facevarx(i,io_ilo:io_ihi+1,io_jlo:io_jhi,io_klo:io_khi,1:MAXBLOCKS) = &
403 faceXBuf(1,1:NXB+1,1:NYB,1:NZB,1:MAXBLOCKS)
404
405 if (NDIM .gt. 1) then
406
407 call io_h5read_unknowns(io_chkptFileID, &
408 NXB, &
409 NYB+1, &
410 NZB, &
411 faceYBuf, &
412 io_faceYVarLabels(i), &
413 localNumBlocks, &
414 io_splitNumBlks, &
415 localOffset, &
416 doread)
417
418 facevary(i,io_ilo:io_ihi,io_jlo:io_jhi+1,io_klo:io_khi,1:MAXBLOCKS) = &
419 faceYBuf(1,1:NXB,1:NYB+1,1:NZB,1:MAXBLOCKS)
420
421 end if
422
423 if (NDIM .gt. 2) then
424 call io_h5read_unknowns(io_chkptFileID, &
425 NXB, &
426 NYB, &
427 NZB+1, &
428 faceZBuf, &
429 io_faceZVarLabels(i), &
430 localNumBlocks, &
431 io_splitNumBlks, &
432 localOffset, &
433 doread)
434
435 facevarz(i,io_ilo:io_ihi,io_jlo:io_jhi,io_klo:io_khi+1,1:MAXBLOCKS) = &
436 faceZBuf(1,1:NXB,1:NYB,1:NZB+1,1:MAXBLOCKS)
437
438 end if
439
440 end do
441
442 deallocate(faceXBuf)
443 if (NDIM .gt. 1) deallocate(faceYBuf)
444 if (NDIM .gt. 2) deallocate(faceZBuf)
445
446#endif
447#endif
448
449 if (io_globalMe == MASTER_PE) &
450 print *, 'read_data: read ', gr_globalNumBlocks, ' blocks.'
451
452
453 if (io_globalMe == MASTER_PE) then
454 allocate (strBuff(2, 2))
455 write (strBuff(1,1), "(A)") "type"
456 write (strBuff(1,2), "(A)") "checkpoint"
457 write (strBuff(2,1), "(A)") "name"
458 write (strBuff(2,2), "(A)") trim(filename)
459 call Logfile_stamp( strBuff, 2, 2, "[io_readData] file_closed")
460 if (allocated(strBuff)) deallocate(strBuff)
461 end if
462
463
464 call MPI_BARRIER (io_globalComm, ierr)
465 if (io_globalMe == MASTER_PE) &
466 print *, 'io_readData: finished reading input file.'
467
468 return
469end subroutine io_readData
Note: See TracBrowser for help on using the repository browser.