dataset_global_hdf5.F90 coverage: 75.00 %func 65.09 %block
1) module Dataset_Global_HDF5_class
2)
3) use Dataset_Common_HDF5_class
4) use DM_Kludge_module
5)
6) use PFLOTRAN_Constants_module
7)
8) implicit none
9)
10) private
11)
12) #include "petsc/finclude/petscsys.h"
13)
14) type, public, extends(dataset_common_hdf5_type) :: dataset_global_hdf5_type
15) PetscInt :: local_size ! local number of entries on this process
16) PetscInt :: global_size ! global number of entries
17) type(dm_ptr_type), pointer :: dm_wrapper ! pointer
18) end type dataset_global_hdf5_type
19)
20) PetscInt, parameter :: default_max_buffer_size = 10
21)
22) public :: DatasetGlobalHDF5Create, &
23) DatasetGlobalHDF5Init, &
24) DatasetGlobalHDF5Cast, &
25) DatasetGlobalHDF5Load, &
26) DatasetGlobalHDF5Print, &
27) DatasetGlobalHDF5Destroy
28)
29) contains
30)
31) ! ************************************************************************** !
32)
33) function DatasetGlobalHDF5Create()
34) !
35) ! Creates global dataset class
36) !
37) ! Author: Glenn Hammond
38) ! Date: 05/03/13
39) !
40)
41) implicit none
42)
43) class(dataset_global_hdf5_type), pointer :: dataset
44)
45) class(dataset_global_hdf5_type), pointer :: DatasetGlobalHDF5Create
46)
47) allocate(dataset)
48) call DatasetGlobalHDF5Init(dataset)
49)
50) DatasetGlobalHDF5Create => dataset
51)
52) end function DatasetGlobalHDF5Create
53)
54) ! ************************************************************************** !
55)
56) subroutine DatasetGlobalHDF5Init(this)
57) !
58) ! Initializes members of global dataset class
59) !
60) ! Author: Glenn Hammond
61) ! Date: 05/03/13
62) !
63)
64) implicit none
65)
66) class(dataset_global_hdf5_type) :: this
67)
68) call DatasetCommonHDF5Init(this)
69) this%local_size = 0
70) this%global_size = 0
71) nullify(this%dm_wrapper)
72)
73) end subroutine DatasetGlobalHDF5Init
74)
75) ! ************************************************************************** !
76)
77) function DatasetGlobalHDF5Cast(this)
78) !
79) ! Casts a dataset_base_type to database_global_type
80) !
81) ! Author: Glenn Hammond
82) ! Date: 05/03/13
83) !
84)
85) use Dataset_Base_class
86)
87) implicit none
88)
89) class(dataset_base_type), pointer :: this
90)
91) class(dataset_global_hdf5_type), pointer :: DatasetGlobalHDF5Cast
92)
93) nullify(DatasetGlobalHDF5Cast)
94) select type (this)
95) class is (dataset_global_hdf5_type)
96) DatasetGlobalHDF5Cast => this
97) end select
98)
99) end function DatasetGlobalHDF5Cast
100)
101) ! ************************************************************************** !
102)
103) subroutine DatasetGlobalHDF5Load(this,option)
104) !
105) ! Load new data into dataset buffer
106) !
107) ! Author: Glenn Hammond
108) ! Date: 05/03/13
109) !
110)
111) #if defined(PETSC_HAVE_HDF5)
112) use hdf5, only : H5T_NATIVE_DOUBLE
113) #endif
114) use Option_module
115) use Time_Storage_module
116) use Dataset_Base_class
117)
118) implicit none
119)
120) class(dataset_global_hdf5_type) :: this
121) type(option_type) :: option
122)
123) if (.not.associated(this%dm_wrapper)) then
124) option%io_buffer = 'dm_wrapper not associated in Global Dataset: ' // &
125) trim(this%name)
126) call printErrMsg(option)
127) endif
128)
129) if (DatasetCommonHDF5Load(this,option)) then
130) if (.not.associated(this%rarray)) then
131) if (this%local_size == 0) then
132) option%io_buffer = 'Local size of Global Dataset has not been set.'
133) call printErrMsg(option)
134) endif
135) allocate(this%rarray(this%local_size))
136) this%rarray = 0.d0
137) endif
138) if (.not.associated(this%rbuffer)) then ! not initialized
139) if (this%max_buffer_size < 0) then
140) this%max_buffer_size = default_max_buffer_size
141) endif
142) this%buffer_nslice = min(this%max_buffer_size, &
143) this%time_storage%max_time_index)
144) allocate(this%rbuffer(this%local_size*this%buffer_nslice))
145) this%rbuffer = 0.d0
146) endif
147) #if defined(PETSC_HAVE_HDF5)
148) call DatasetGlobalHDF5ReadData(this,option,H5T_NATIVE_DOUBLE)
149) #endif
150) ! no need to reorder since it is 1D in the h5 file.
151) endif
152) call DatasetBaseInterpolateTime(this)
153)
154) end subroutine DatasetGlobalHDF5Load
155)
156) #if defined(PETSC_HAVE_HDF5)
157)
158) ! ************************************************************************** !
159)
160) subroutine DatasetGlobalHDF5ReadData(this,option,data_type)
161) !
162) ! Read an hdf5 array into a Petsc Vec
163) !
164) ! Author: Glenn Hammond
165) ! Date: 01/12/08
166) !
167)
168) use hdf5
169) use Logging_module
170) use Option_module
171) use HDF5_Aux_module
172)
173) implicit none
174)
175) #include "petsc/finclude/petscvec.h"
176) #include "petsc/finclude/petscvec.h90"
177) #include "petsc/finclude/petscdm.h"
178) #include "petsc/finclude/petscdm.h90"
179) #include "petsc/finclude/petscdmda.h"
180) #include "petsc/finclude/petscviewer.h"
181)
182) ! Default HDF5 Mechanism
183)
184) class(dataset_global_hdf5_type) :: this
185) type(option_type) :: option
186) integer(HID_T) :: data_type
187)
188) integer(HID_T) :: file_id
189) integer(HID_T) :: file_space_id
190) integer(HID_T) :: memory_space_id
191) integer(HID_T) :: data_set_id
192) integer(HID_T) :: prop_id
193) integer(HSIZE_T) :: dims(3), max_dims(3)
194) integer(HSIZE_T) :: offset(3), length(3), stride(3)
195) PetscMPIInt :: ndims
196) PetscMPIInt :: rank_mpi
197) Vec :: natural_vec
198) Vec :: global_vec
199) PetscReal, pointer :: vec_ptr(:)
200) character(len=MAXSTRINGLENGTH) :: string
201) character(len=MAXWORDLENGTH) :: word
202) PetscInt :: i, istart
203) PetscInt :: buffer_size, buffer_rank2_size, file_rank1_size, file_rank2_size
204) integer, allocatable :: integer_buffer_i4(:)
205) PetscErrorCode :: ierr
206) PetscMPIInt :: hdf5_err
207)
208) PetscViewer :: viewer
209)
210) call PetscLogEventBegin(logging%event_read_array_hdf5,ierr);CHKERRQ(ierr)
211)
212) if (this%dm_wrapper%dm /= 0) then
213) call DMCreateGlobalVector(this%dm_wrapper%dm,global_vec, &
214) ierr);CHKERRQ(ierr)
215) call DMDACreateNaturalVector(this%dm_wrapper%dm,natural_vec, &
216) ierr);CHKERRQ(ierr)
217) else
218) option%io_buffer = 'ugdm not yet supported in DatasetGlobalHDF5ReadData()'
219) call printErrMsg(option)
220) endif
221)
222) call VecZeroEntries(natural_vec,ierr);CHKERRQ(ierr)
223) call VecZeroEntries(global_vec,ierr);CHKERRQ(ierr)
224)
225) ! open the file
226) call h5open_f(hdf5_err)
227) option%io_buffer = 'Opening hdf5 file: ' // trim(this%filename)
228) call printMsg(option)
229)
230) ! set read file access property
231) call h5pcreate_f(H5P_FILE_ACCESS_F,prop_id,hdf5_err)
232) #ifndef SERIAL_HDF5
233) call h5pset_fapl_mpio_f(prop_id,option%mycomm,MPI_INFO_NULL,hdf5_err)
234) #endif
235) call HDF5OpenFileReadOnly(this%filename,file_id,prop_id,option)
236) call h5pclose_f(prop_id,hdf5_err)
237)
238) string = trim(this%hdf5_dataset_name) // '/Data'
239) if (this%realization_dependent) then
240) write(word,'(i9)') option%id
241) string = trim(string) // trim(adjustl(word))
242) endif
243) option%io_buffer = 'Opening data set: ' // trim(string)
244) call printMsg(option)
245) call h5dopen_f(file_id,string,data_set_id,hdf5_err)
246) call h5dget_space_f(data_set_id,file_space_id,hdf5_err)
247)
248) call h5sget_simple_extent_ndims_f(file_space_id,ndims,hdf5_err)
249) call h5sget_simple_extent_dims_f(file_space_id,dims,max_dims,hdf5_err)
250)
251) buffer_size = size(this%rbuffer)
252) buffer_rank2_size = buffer_size / this%local_size
253) file_rank1_size = dims(1)
254) if (ndims > 1) then
255) file_rank2_size = dims(2)
256) else
257) if (option%mycommsize > 1) then
258) option%io_buffer = 'Dataset "' // trim(this%hdf5_dataset_name) // &
259) '" in file "' // trim (this%filename) // &
260) '" must be a 2D dataset (time,cell) if PFLOTRAN is run in parallel.'
261) call printErrMsg(option)
262) endif
263) file_rank2_size = 1
264) endif
265)
266) allocate(this%dims(1))
267) this%dims = size(this%rarray)
268) this%rank = 1
269)
270) if (mod(buffer_size,this%local_size) /= 0) then
271) write(option%io_buffer, &
272) '(a," buffer dimension (",i9,") is not a multiple of local domain",&
273) &" dimension (",i9,").")') trim(this%hdf5_dataset_name), &
274) size(this%rbuffer,1), this%local_size
275) call printErrMsg(option)
276) endif
277)
278) if (mod(file_rank1_size,this%global_size) /= 0) then
279) write(option%io_buffer, &
280) '(a," data space dimension (",i9,") is not a multiple of domain",&
281) &" dimension (",i9,").")') trim(this%hdf5_dataset_name), &
282) file_rank1_size, this%global_size
283) call printErrMsg(option)
284) endif
285)
286) istart = 0
287) call MPI_Exscan(this%local_size,istart,ONE_INTEGER_MPI,MPIU_INTEGER,MPI_SUM, &
288) option%mycomm,ierr)
289)
290) call h5pcreate_f(H5P_DATASET_XFER_F,prop_id,hdf5_err)
291) #ifndef SERIAL_HDF5
292) call h5pset_dxpl_mpio_f(prop_id,H5FD_MPIO_INDEPENDENT_F,hdf5_err)
293) #endif
294)
295) ! must initialize here to avoid error below when closing memory space
296) memory_space_id = -1
297)
298) ! offset is zero-based
299) offset = 0
300) length = 0
301) stride = 1
302)
303) offset(1) = istart ! istart is offset in the first dimension
304) length(1) = file_rank1_size
305) if (ndims > 1) then
306) offset(2) = this%buffer_slice_offset
307) length(1) = this%local_size
308) length(2) = min(buffer_rank2_size,file_rank2_size-this%buffer_slice_offset)
309) else
310) offset(1) = this%buffer_slice_offset*this%local_size
311) length(1) = min(buffer_size, &
312) file_rank1_size-this%buffer_slice_offset*this%local_size)
313) endif
314) call h5sselect_hyperslab_f(file_space_id, H5S_SELECT_SET_F,offset, &
315) length,hdf5_err,stride,stride)
316)
317) dims = 0
318) rank_mpi = 1
319) if (ndims > 1) then
320) dims(1) = min(buffer_size, &
321) this%local_size*(file_rank2_size-this%buffer_slice_offset))
322) else
323) dims(1) = min(buffer_size, &
324) file_rank1_size-this%buffer_slice_offset*this%local_size)
325) endif
326) call h5screate_simple_f(rank_mpi,dims,memory_space_id,hdf5_err,dims)
327)
328) ! initialize to UNINITIALIZED_INTEGER to catch errors
329) this%rbuffer = UNINITIALIZED_DOUBLE
330)
331) if (data_type == H5T_NATIVE_DOUBLE) then
332) call PetscLogEventBegin(logging%event_h5dread_f,ierr);CHKERRQ(ierr)
333) call h5dread_f(data_set_id,H5T_NATIVE_DOUBLE,this%rbuffer,dims, &
334) hdf5_err,memory_space_id,file_space_id,prop_id)
335) call PetscLogEventEnd(logging%event_h5dread_f,ierr);CHKERRQ(ierr)
336) else if (data_type == H5T_NATIVE_INTEGER) then
337) allocate(integer_buffer_i4(size(this%rbuffer)))
338) call PetscLogEventBegin(logging%event_h5dread_f,ierr);CHKERRQ(ierr)
339) call h5dread_f(data_set_id,H5T_NATIVE_INTEGER,integer_buffer_i4,dims, &
340) hdf5_err,memory_space_id,file_space_id,prop_id)
341) call PetscLogEventEnd(logging%event_h5dread_f,ierr);CHKERRQ(ierr)
342) do i = 1, min(buffer_size, &
343) this%local_size*(file_rank2_size-this%buffer_slice_offset))
344) this%rbuffer(i) = real(integer_buffer_i4(i))
345) enddo
346) deallocate(integer_buffer_i4)
347) endif
348)
349) call h5pclose_f(prop_id,hdf5_err)
350) if (memory_space_id > -1) call h5sclose_f(memory_space_id,hdf5_err)
351) call h5sclose_f(file_space_id,hdf5_err)
352) string = trim(this%hdf5_dataset_name) // '/Data'
353) option%io_buffer = 'Closing data set: ' // trim(string)
354) call printMsg(option)
355) call h5dclose_f(data_set_id,hdf5_err)
356) option%io_buffer = 'Closing hdf5 file: ' // trim(this%filename)
357) call printMsg(option)
358) call h5fclose_f(file_id,hdf5_err)
359) call h5close_f(hdf5_err)
360)
361) istart = 0
362) do i = 1, min(buffer_rank2_size,file_rank2_size-this%buffer_slice_offset)
363) call VecGetArrayF90(natural_vec,vec_ptr,ierr);CHKERRQ(ierr)
364) vec_ptr(:) = this%rbuffer(istart+1:istart+this%local_size)
365) call VecRestoreArrayF90(natural_vec,vec_ptr,ierr);CHKERRQ(ierr)
366)
367) !geh: for debugging purposes
368) !write(string,*) i
369) !string = trim(adjustl(this%hdf5_dataset_name)) // '_' // &
370) ! trim(adjustl(string)) // '_natural.txt'
371) !call PetscViewerASCIIOpen(option%mycomm,string,viewer,ierr)
372) !call VecView(natural_vec,viewer,ierr)
373) !call PetscViewerDestroy(viewer,ierr)
374)
375) call DMDANaturalToGlobalBegin(this%dm_wrapper%dm,natural_vec, &
376) INSERT_VALUES,global_vec,ierr);CHKERRQ(ierr)
377) call DMDANaturalToGlobalEnd(this%dm_wrapper%dm,natural_vec, &
378) INSERT_VALUES,global_vec,ierr);CHKERRQ(ierr)
379) call VecGetArrayF90(global_vec,vec_ptr,ierr);CHKERRQ(ierr)
380) this%rbuffer(istart+1:istart+this%local_size) = vec_ptr(:)
381) call VecRestoreArrayF90(global_vec,vec_ptr,ierr);CHKERRQ(ierr)
382)
383) !geh: for debugging purposes
384) !write(string,*) i
385) !string = trim(adjustl(this%hdf5_dataset_name)) // '_' // &
386) ! trim(adjustl(string)) // '_global.txt'
387) !call PetscViewerASCIIOpen(option%mycomm,string,viewer,ierr)
388) !call VecView(global_vec,viewer,ierr)
389) !call PetscViewerDestroy(viewer,ierr)
390)
391) istart = istart + this%local_size
392) enddo
393)
394) call VecDestroy(natural_vec,ierr);CHKERRQ(ierr)
395) call VecDestroy(global_vec,ierr);CHKERRQ(ierr)
396)
397) call PetscLogEventEnd(logging%event_read_array_hdf5,ierr);CHKERRQ(ierr)
398) ! End of Default HDF5 Mechanism
399)
400) end subroutine DatasetGlobalHDF5ReadData
401) #endif
402)
403) ! ************************************************************************** !
404)
405) subroutine DatasetGlobalHDF5Print(this,option)
406) !
407) ! Prints dataset info
408) !
409) ! Author: Glenn Hammond
410) ! Date: 10/22/13
411) !
412)
413) use Option_module
414)
415) implicit none
416)
417) class(dataset_global_hdf5_type), target :: this
418) type(option_type) :: option
419)
420) class(dataset_common_hdf5_type), pointer :: dataset_hdf5
421)
422) dataset_hdf5 => this
423) call DatasetCommonHDF5Print(this,option)
424)
425) ! no need to print local_size as it varies by process
426) write(option%fid_out,'(10x,''Global Size: '',i2)') this%global_size
427)
428) end subroutine DatasetGlobalHDF5Print
429)
430) ! ************************************************************************** !
431)
432) subroutine DatasetGlobalHDF5Strip(this)
433) !
434) ! Strips allocated objects within Global dataset object
435) !
436) ! Author: Glenn Hammond
437) ! Date: 05/03/13
438) !
439)
440) use Utility_module, only : DeallocateArray
441)
442) implicit none
443)
444) class(dataset_global_hdf5_type) :: this
445)
446) call DatasetCommonHDF5Strip(this)
447) nullify(this%dm_wrapper) ! do not deallocate, as this is solely a pointer
448)
449) end subroutine DatasetGlobalHDF5Strip
450)
451) ! ************************************************************************** !
452)
453) subroutine DatasetGlobalHDF5Destroy(this)
454) !
455) ! Destroys a dataset
456) !
457) ! Author: Glenn Hammond
458) ! Date: 01/12/11
459) !
460)
461) implicit none
462)
463) class(dataset_global_hdf5_type), pointer :: this
464)
465) if (.not.associated(this)) return
466)
467) call DatasetGlobalHDF5Strip(this)
468)
469) deallocate(this)
470) nullify(this)
471)
472) end subroutine DatasetGlobalHDF5Destroy
473)
474) end module Dataset_Global_HDF5_class