dataset_gridded_hdf5.F90 coverage: 92.31 %func 53.72 %block
1) module Dataset_Gridded_HDF5_class
2)
3) use Dataset_Common_HDF5_class
4)
5) use PFLOTRAN_Constants_module
6)
7) implicit none
8)
9) private
10)
11) #include "petsc/finclude/petscsys.h"
12)
13) type, public, extends(dataset_common_hdf5_type) :: dataset_gridded_hdf5_type
14) PetscBool :: is_cell_centered
15) PetscInt :: data_dim
16) PetscInt :: interpolation_method
17) PetscReal, pointer :: origin(:)
18) PetscReal, pointer :: extent(:)
19) PetscReal, pointer :: discretization(:)
20) end type dataset_gridded_hdf5_type
21)
22) PetscInt, parameter, public :: DIM_NULL = 0
23) PetscInt, parameter, public :: DIM_X = 1
24) PetscInt, parameter, public :: DIM_Y = 2
25) PetscInt, parameter, public :: DIM_Z = 3
26) PetscInt, parameter, public :: DIM_XY = 4
27) PetscInt, parameter, public :: DIM_XZ = 5
28) PetscInt, parameter, public :: DIM_YZ = 6
29) PetscInt, parameter, public :: DIM_XYZ = 7
30)
31) PetscInt, parameter :: default_max_buffer_size = 10
32)
33) public :: DatasetGriddedHDF5Create, &
34) DatasetGriddedHDF5Init, &
35) DatasetGriddedHDF5Cast, &
36) DatasetGriddedHDF5Load, &
37) DatasetGriddedHDF5InterpolateReal, &
38) DatasetGriddedHDF5Print, &
39) DatasetGriddedHDF5Strip, &
40) DatasetGriddedHDF5Destroy
41)
42) contains
43)
44) ! ************************************************************************** !
45)
46) function DatasetGriddedHDF5Create()
47) !
48) ! Creates global dataset class
49) !
50) ! Author: Glenn Hammond
51) ! Date: 05/29/13
52) !
53)
54) implicit none
55)
56) class(dataset_gridded_hdf5_type), pointer :: dataset
57)
58) class(dataset_gridded_hdf5_type), pointer :: DatasetGriddedHDF5Create
59)
60) allocate(dataset)
61) call DatasetGriddedHDF5Init(dataset)
62)
63) DatasetGriddedHDF5Create => dataset
64)
65) end function DatasetGriddedHDF5Create
66)
67) ! ************************************************************************** !
68)
69) subroutine DatasetGriddedHDF5Init(this)
70) !
71) ! Initializes members of global dataset class
72) !
73) ! Author: Glenn Hammond
74) ! Date: 05/29/13
75) !
76)
77) use Dataset_Base_class
78)
79) implicit none
80)
81) class(dataset_gridded_hdf5_type) :: this
82)
83) call DatasetCommonHDF5Init(this)
84) this%is_cell_centered = PETSC_FALSE
85) this%interpolation_method = INTERPOLATION_LINEAR
86) this%data_dim = DIM_NULL
87) nullify(this%origin)
88) nullify(this%extent)
89) nullify(this%discretization)
90)
91) end subroutine DatasetGriddedHDF5Init
92)
93) ! ************************************************************************** !
94)
95) function DatasetGriddedHDF5Cast(this)
96) !
97) ! Casts a dataset_base_type to dataset_gridded_hdf5_type
98) !
99) ! Author: Glenn Hammond
100) ! Date: 08/29/13
101) !
102)
103) use Dataset_Base_class
104)
105) implicit none
106)
107) class(dataset_base_type), pointer :: this
108)
109) class(dataset_gridded_hdf5_type), pointer :: DatasetGriddedHDF5Cast
110)
111) nullify(DatasetGriddedHDF5Cast)
112) select type (this)
113) class is (dataset_gridded_hdf5_type)
114) DatasetGriddedHDF5Cast => this
115) class default
116) ! to catch a class that is not gridded.
117) end select
118)
119) end function DatasetGriddedHDF5Cast
120)
121) ! ************************************************************************** !
122)
123) subroutine DatasetGriddedHDF5Load(this,option)
124) !
125) ! Load new data into dataset buffer
126) !
127) ! Author: Glenn Hammond
128) ! Date: 05/29/13
129) !
130)
131) use Option_module
132) use Time_Storage_module
133) use Dataset_Base_class
134)
135) implicit none
136)
137) class(dataset_gridded_hdf5_type) :: this
138) type(option_type) :: option
139)
140) if (DatasetCommonHDF5Load(this,option)) then
141) #if defined(PETSC_HAVE_HDF5)
142) call DatasetGriddedHDF5ReadData(this,option)
143) #endif
144) ! call this%Reorder(option)
145) call DatasetBaseReorder(this,option)
146) endif
147) call DatasetBaseInterpolateTime(this)
148)
149) end subroutine DatasetGriddedHDF5Load
150)
151) #if defined(PETSC_HAVE_HDF5)
152)
153) ! ************************************************************************** !
154)
155) subroutine DatasetGriddedHDF5ReadData(this,option)
156) !
157) ! Read an hdf5 data into arrays
158) !
159) ! Author: Glenn Hammond
160) ! Date: 10/25/11, 05/29/13
161) !
162)
163) use hdf5
164) use Option_module
165) use Units_module
166) use Logging_module
167) use HDF5_Aux_module
168) use String_module
169)
170) implicit none
171)
172) class(dataset_gridded_hdf5_type) :: this
173) type(option_type) :: option
174)
175) integer(HID_T) :: file_id
176) integer(HID_T) :: file_space_id
177) integer(HID_T) :: memory_space_id
178) integer(HID_T) :: dataset_id
179) integer(HID_T) :: prop_id
180) integer(HID_T) :: grp_id
181) integer(HID_T) :: attribute_id
182) integer(HID_T) :: ndims_h5
183) integer(HID_T) :: atype_id
184) integer(HSIZE_T), allocatable :: dims_h5(:), max_dims_h5(:)
185) integer(HSIZE_T) :: attribute_dim(3)
186) integer(HSIZE_T) :: offset(4), length(4), stride(4)
187) integer(HSIZE_T) :: num_data_values
188) integer(SIZE_T) size_t_int
189) PetscInt :: i, temp_int, temp_array(4)
190) PetscInt :: num_spatial_dims, time_dim, num_times
191) PetscInt :: num_dims_in_h5_file, num_times_in_h5_file
192) PetscMPIInt :: array_rank_mpi, mpi_int
193) PetscBool :: attribute_exists
194) PetscBool :: first_time
195) PetscMPIInt :: hdf5_err
196) PetscErrorCode :: ierr
197) PetscLogDouble :: tstart, tend
198) character(len=MAXWORDLENGTH) :: attribute_name, dataset_name, word
199)
200) call PetscLogEventBegin(logging%event_dataset_gridded_hdf5_read, &
201) ierr);CHKERRQ(ierr)
202)
203) first_time = (this%data_dim == DIM_NULL)
204)
205) !#define TIME_DATASET
206) #ifdef TIME_DATASET
207) call PetscTime(tstart,ierr);CHKERRQ(ierr)
208) #endif
209)
210) #define BROADCAST_DATASET
211) #ifdef BROADCAST_DATASET
212) if (first_time .or. option%myrank == option%io_rank) then
213) #endif
214)
215) ! open the file
216) call h5open_f(hdf5_err)
217) option%io_buffer = 'Opening hdf5 file: ' // trim(this%filename)
218) call printMsg(option)
219)
220) ! set read file access property
221) call h5pcreate_f(H5P_FILE_ACCESS_F,prop_id,hdf5_err)
222) #ifndef SERIAL_HDF5
223) !geh: I don't believe that we ever need this
224) ! if (first_time) then
225) ! call h5pset_fapl_mpio_f(prop_id,option%mycomm,MPI_INFO_NULL,hdf5_err)
226) ! endif
227) #endif
228) call HDF5OpenFileReadOnly(this%filename,file_id,prop_id,option)
229) call h5pclose_f(prop_id,hdf5_err)
230)
231) ! the dataset is actually stored in a group. the group contains
232) ! a "data" dataset and optionally a "time" dataset.
233) option%io_buffer = 'Opening group: ' // trim(this%hdf5_dataset_name)
234) call printMsg(option)
235) call h5gopen_f(file_id,this%hdf5_dataset_name,grp_id,hdf5_err)
236)
237) if (hdf5_err < 0) then
238) option%io_buffer = 'A group named "' // trim(this%hdf5_dataset_name) // &
239) '" not found in HDF5 file "' // trim(this%filename) // '".'
240) call printErrMsg(option)
241) endif
242)
243) ! only want to read on first time through
244) if (this%data_dim == DIM_NULL) then
245) ! read in attributes if they exist
246) attribute_name = "Dimension"
247) call H5aexists_f(grp_id,attribute_name,attribute_exists,hdf5_err)
248) if (attribute_exists) then
249) attribute_dim = 1
250) call h5tcopy_f(H5T_NATIVE_CHARACTER,atype_id,hdf5_err)
251) size_t_int = MAXWORDLENGTH
252) call h5tset_size_f(atype_id,size_t_int,hdf5_err)
253) call h5aopen_f(grp_id,attribute_name,attribute_id,hdf5_err)
254) call h5aread_f(attribute_id,atype_id,word,attribute_dim,hdf5_err)
255) call h5aclose_f(attribute_id,hdf5_err)
256) call h5tclose_f(atype_id,hdf5_err)
257) ! set dimensionality of dataset
258) call DatasetGriddedHDF5SetDimension(this,word)
259) else
260) option%io_buffer = &
261) 'Dimension attribute must be included in hdf5 dataset file.'
262) call printErrMsg(option)
263) endif
264) attribute_name = "Discretization"
265) call H5aexists_f(grp_id,attribute_name,attribute_exists,hdf5_err)
266) if (attribute_exists) then
267) attribute_dim(1) = DatasetGriddedHDF5GetNDimensions(this)
268) allocate(this%discretization(attribute_dim(1)))
269) call h5aopen_f(grp_id,attribute_name,attribute_id,hdf5_err)
270) call h5aread_f(attribute_id,H5T_NATIVE_DOUBLE,this%discretization, &
271) attribute_dim,hdf5_err)
272) call h5aclose_f(attribute_id,hdf5_err)
273) else
274) option%io_buffer = &
275) '"Discretization" attribute must be included in GRIDDED hdf5 ' // &
276) 'dataset file.'
277) call printErrMsg(option)
278) endif
279) attribute_name = "Origin"
280) call H5aexists_f(grp_id,attribute_name,attribute_exists,hdf5_err)
281) if (attribute_exists) then
282) attribute_dim(1) = DatasetGriddedHDF5GetNDimensions(this)
283) allocate(this%origin(attribute_dim(1)))
284) call h5aopen_f(grp_id,attribute_name,attribute_id,hdf5_err)
285) call h5aread_f(attribute_id,H5T_NATIVE_DOUBLE,this%origin, &
286) attribute_dim,hdf5_err)
287) call h5aclose_f(attribute_id,hdf5_err)
288) endif
289) attribute_name = "Cell Centered"
290) call H5aexists_f(grp_id,attribute_name,attribute_exists,hdf5_err)
291) if (attribute_exists) then
292) this%is_cell_centered = PETSC_TRUE
293) endif
294) attribute_name = "Interpolation Method"
295) call H5aexists_f(grp_id,attribute_name,attribute_exists,hdf5_err)
296) if (attribute_exists) then
297) call h5tcopy_f(H5T_NATIVE_CHARACTER,atype_id,hdf5_err)
298) size_t_int = MAXWORDLENGTH
299) call h5tset_size_f(atype_id,size_t_int,hdf5_err)
300) call h5aopen_f(grp_id,attribute_name,attribute_id,hdf5_err)
301) call h5aread_f(attribute_id,atype_id,word,attribute_dim,hdf5_err)
302) call h5aclose_f(attribute_id,hdf5_err)
303) call h5tclose_f(atype_id,hdf5_err)
304) call StringToUpper(word)
305) select case(trim(word))
306) case('STEP')
307) this%interpolation_method = INTERPOLATION_STEP
308) case('LINEAR')
309) this%interpolation_method = INTERPOLATION_LINEAR
310) case default
311) option%io_buffer = '"Interpolation Method" not recognized in ' // &
312) 'Gridded HDF5 Dataset "' // trim(this%name) // '".'
313) call printErrMsg(option)
314) end select
315) endif
316) ! this%max_buffer_size is initially set to UNINITIALIZED_INTEGER to force initializaion
317) ! either here, or in the reading of the dataset block.
318) attribute_name = "Max Buffer Size"
319) call H5aexists_f(grp_id,attribute_name,attribute_exists,hdf5_err)
320) if (attribute_exists .and. this%max_buffer_size < 0) then
321) call h5aopen_f(grp_id,attribute_name,attribute_id,hdf5_err)
322) attribute_dim(1) = 1
323) call h5aread_f(attribute_id,H5T_NATIVE_INTEGER,this%max_buffer_size, &
324) attribute_dim,hdf5_err)
325) call h5aclose_f(attribute_id,hdf5_err)
326) else if (this%max_buffer_size < 0) then
327) this%max_buffer_size = default_max_buffer_size
328) endif
329) endif ! this%data_dim == DIM_NULL
330)
331) #ifdef BROADCAST_DATASET
332) endif
333) #endif
334)
335) num_spatial_dims = DatasetGriddedHDF5GetNDimensions(this)
336)
337) ! num_times and time_dim must be calcualted by all processes; does not
338) ! require communication
339) time_dim = -1
340) num_times = 1
341) if (associated(this%time_storage)) then
342) num_times = this%time_storage%max_time_index
343) time_dim = num_spatial_dims + 1
344) endif
345)
346) #ifdef BROADCAST_DATASET
347) if (first_time .or. option%myrank == option%io_rank) then
348) #endif
349) ! open the "data" dataset
350) dataset_name = 'Data'
351) if (this%realization_dependent) then
352) write(word,'(i9)') option%id
353) dataset_name = trim(dataset_name) // trim(adjustl(word))
354) endif
355) call h5dopen_f(grp_id,dataset_name,dataset_id,hdf5_err)
356) if (hdf5_err < 0) then
357) option%io_buffer = 'A dataset named "Data" not found in HDF5 file "' // &
358) trim(this%filename) // '".'
359) call printErrMsg(option)
360) endif
361) call h5dget_space_f(dataset_id,file_space_id,hdf5_err)
362)
363) ! get dataset dimensions
364) if (.not.associated(this%dims)) then
365) call h5sget_simple_extent_ndims_f(file_space_id,ndims_h5,hdf5_err)
366) allocate(dims_h5(ndims_h5))
367) allocate(max_dims_h5(ndims_h5))
368) call h5sget_simple_extent_dims_f(file_space_id,dims_h5,max_dims_h5,hdf5_err)
369) if (associated(this%time_storage)) then
370) ! if dataset is time dependent, need to remove that time index from
371) ! dimensions and decrement rank (we don't want to include time)
372) this%rank = ndims_h5-1
373) ! the first dimension of dims_h5 is the time dimension
374) num_times_in_h5_file = dims_h5(1)
375) else
376) this%rank = ndims_h5
377) num_times_in_h5_file = 0
378) endif
379) allocate(this%dims(this%rank))
380) ! have to invert dimensions, but do not count the time dimension
381) do i = 1, this%rank
382) this%dims(i) = int(dims_h5(ndims_h5-i+1))
383) enddo
384) deallocate(dims_h5)
385) deallocate(max_dims_h5)
386)
387) allocate(this%extent(num_spatial_dims))
388) do i = 1, num_spatial_dims
389) temp_int = this%dims(i)
390) if (.not.this%is_cell_centered) then
391) temp_int = temp_int - 1
392) endif
393) this%extent(i) = this%origin(i) + this%discretization(i) * temp_int
394) enddo
395)
396) call h5sget_simple_extent_npoints_f(file_space_id,num_data_values,hdf5_err)
397)
398) temp_int = this%dims(1)
399) do i = 2, num_spatial_dims
400) temp_int = temp_int * this%dims(i)
401) enddo
402) if (num_data_values/num_times /= temp_int) then
403) option%io_buffer = &
404) 'Number of values in dataset does not match dimensions.'
405) call printErrMsg(option)
406) endif
407) if (associated(this%time_storage) .and. &
408) num_times_in_h5_file /= num_times) then
409) option%io_buffer = &
410) 'Number of times does not match last dimension of data array.'
411) call printErrMsg(option)
412) endif
413) if (.not.associated(this%rarray)) then
414) allocate(this%rarray(temp_int))
415) endif
416) this%rarray = 0.d0
417) if (associated(this%time_storage) .and. .not.associated(this%rbuffer)) then
418) ! buffered array
419) allocate(this%rbuffer(size(this%rarray)*this%max_buffer_size))
420) this%rbuffer = 0.d0
421) endif
422) endif
423)
424) #ifdef BROADCAST_DATASET
425) endif
426) #endif
427)
428) #ifdef TIME_DATASET
429) call MPI_Barrier(option%mycomm,ierr)
430) call PetscTime(tend,ierr);CHKERRQ(ierr)
431) write(option%io_buffer,'(f6.2," Seconds to set up dataset ",a,".")') &
432) tend-tstart, trim(this%hdf5_dataset_name) // ' (' // &
433) trim(option%group_prefix) // ')'
434) if (option%myrank == option%io_rank) then
435) print *, trim(option%io_buffer)
436) endif
437) call PetscTime(tstart,ierr);CHKERRQ(ierr)
438) #endif
439)
440) call PetscLogEventBegin(logging%event_h5dread_f,ierr);CHKERRQ(ierr)
441)
442) if (associated(this%time_storage)) then
443) num_dims_in_h5_file = this%rank + 1
444) else
445) num_dims_in_h5_file = this%rank
446) endif
447)
448) array_rank_mpi = 1
449) length = 1
450) ! length (or size) must be adjusted according to the size of the
451) ! remaining data in the file
452) this%buffer_nslice = min(this%max_buffer_size,(num_times-this%buffer_slice_offset))
453) if (time_dim > 0) then
454) length(1) = size(this%rarray) * this%buffer_nslice
455) else
456) length(1) = size(this%rarray)
457) endif
458)
459) #ifdef BROADCAST_DATASET
460) if (option%myrank == option%io_rank) then
461) #endif
462)
463) call h5screate_simple_f(array_rank_mpi,length,memory_space_id,hdf5_err,length)
464)
465) length = 1
466) stride = 1
467) offset = 0
468) do i = 1, num_spatial_dims
469) length(i) = this%dims(i)
470) enddo
471) ! cannot read beyond end of the buffer
472) if (time_dim > 0) then
473) length(time_dim) = this%buffer_nslice
474) offset(time_dim) = this%buffer_slice_offset
475) endif
476)
477)
478) !geh: for some reason, we have to invert here. Perhaps because the
479) ! dataset was generated in C???
480) temp_array(1:num_dims_in_h5_file) = int(length(1:num_dims_in_h5_file))
481) do i = 1, num_dims_in_h5_file
482) length(i) = temp_array(num_dims_in_h5_file-i+1)
483) enddo
484) temp_array(1:num_dims_in_h5_file) = int(offset(1:num_dims_in_h5_file))
485) do i = 1, num_dims_in_h5_file
486) offset(i) = temp_array(num_dims_in_h5_file-i+1)
487) enddo
488) ! stride is fine
489)
490) call h5pcreate_f(H5P_DATASET_XFER_F,prop_id,hdf5_err)
491) #ifndef SERIAL_HDF5
492) call h5pset_dxpl_mpio_f(prop_id,H5FD_MPIO_INDEPENDENT_F,hdf5_err)
493) #endif
494) call h5sselect_hyperslab_f(file_space_id, H5S_SELECT_SET_F,offset,length, &
495) hdf5_err,stride,stride)
496) if (associated(this%rbuffer)) then
497) call h5dread_f(dataset_id,H5T_NATIVE_DOUBLE,this%rbuffer,length, &
498) hdf5_err,memory_space_id,file_space_id,prop_id)
499) else
500) call h5dread_f(dataset_id,H5T_NATIVE_DOUBLE,this%rarray,length, &
501) hdf5_err,memory_space_id,file_space_id,prop_id)
502) ! this%rmax = maxval(this%rarray)
503) ! this%rmin = minval(this%rarray)
504) endif
505)
506) call h5pclose_f(prop_id,hdf5_err)
507) if (memory_space_id > -1) call h5sclose_f(memory_space_id,hdf5_err)
508) call h5sclose_f(file_space_id,hdf5_err)
509) call h5dclose_f(dataset_id,hdf5_err)
510)
511) #ifdef BROADCAST_DATASET
512) endif !if (option%myrank == option%io_rank) then
513) if (associated(this%rbuffer)) then
514) mpi_int = size(this%rbuffer)
515) call MPI_Bcast(this%rbuffer,mpi_int,MPI_DOUBLE_PRECISION,option%io_rank, &
516) option%mycomm,ierr)
517) else
518) mpi_int = size(this%rarray)
519) call MPI_Bcast(this%rarray,mpi_int,MPI_DOUBLE_PRECISION,option%io_rank, &
520) option%mycomm,ierr)
521) endif
522) #endif
523)
524) call PetscLogEventEnd(logging%event_h5dread_f,ierr);CHKERRQ(ierr)
525)
526) #ifdef BROADCAST_DATASET
527) if (first_time .or. option%myrank == option%io_rank) then
528) #endif
529) option%io_buffer = 'Closing group: ' // trim(this%hdf5_dataset_name)
530) call printMsg(option)
531) call h5gclose_f(grp_id,hdf5_err)
532) option%io_buffer = 'Closing hdf5 file: ' // trim(this%filename)
533) call printMsg(option)
534) call h5fclose_f(file_id,hdf5_err)
535) call h5close_f(hdf5_err)
536) #ifdef BROADCAST_DATASET
537) endif
538) #endif
539)
540) #ifdef TIME_DATASET
541) call MPI_Barrier(option%mycomm,ierr)
542) call PetscTime(tend,ierr);CHKERRQ(ierr)
543) write(option%io_buffer,'(f6.2," Seconds to read dataset ",a,".")') &
544) tend-tstart, trim(this%hdf5_dataset_name) // ' (' // &
545) trim(option%group_prefix) // ')'
546) if (option%myrank == option%io_rank) then
547) print *, trim(option%io_buffer)
548) endif
549) #endif
550)
551) call PetscLogEventEnd(logging%event_dataset_gridded_hdf5_read, &
552) ierr);CHKERRQ(ierr)
553)
554) end subroutine DatasetGriddedHDF5ReadData
555)
556) #endif
557)
558) ! ************************************************************************** !
559)
560) subroutine DatasetGriddedHDF5SetDimension(this,word)
561) !
562) ! Sets the dimension of the dataset
563) !
564) ! Author: Glenn Hammond
565) ! Date: 10/24/11, 05/29/13
566) !
567)
568) use String_module
569)
570) implicit none
571)
572) class(dataset_gridded_hdf5_type) :: this
573) character(len=MAXWORDLENGTH) :: word
574)
575) call StringToUpper(word)
576) select case(word)
577) case('X')
578) this%data_dim = DIM_X
579) case('Y')
580) this%data_dim = DIM_Y
581) case('Z')
582) this%data_dim = DIM_Z
583) case('XY')
584) this%data_dim = DIM_XY
585) case('XZ')
586) this%data_dim = DIM_XZ
587) case('YZ')
588) this%data_dim = DIM_YZ
589) case('XYZ')
590) this%data_dim = DIM_XYZ
591) end select
592)
593) end subroutine DatasetGriddedHDF5SetDimension
594)
595) ! ************************************************************************** !
596)
597) function DatasetGriddedHDF5GetDimensionString(this)
598) !
599) ! Returns a string describing dimension
600) !
601) ! Author: Glenn Hammond
602) ! Date: 10/24/11, 05/29/13, 10/22/13
603) !
604)
605) implicit none
606)
607) class(dataset_gridded_hdf5_type) :: this
608)
609) character(len=MAXWORDLENGTH) :: DatasetGriddedHDF5GetDimensionString
610)
611) select case(this%data_dim)
612) case(DIM_X)
613) DatasetGriddedHDF5GetDimensionString = 'X'
614) case(DIM_Y)
615) DatasetGriddedHDF5GetDimensionString = 'Y'
616) case(DIM_Z)
617) DatasetGriddedHDF5GetDimensionString = 'Z'
618) case(DIM_XY)
619) DatasetGriddedHDF5GetDimensionString = 'XY'
620) case(DIM_XZ)
621) DatasetGriddedHDF5GetDimensionString = 'XZ'
622) case(DIM_YZ)
623) DatasetGriddedHDF5GetDimensionString = 'YZ'
624) case(DIM_XYZ)
625) DatasetGriddedHDF5GetDimensionString = 'XYZ'
626) end select
627)
628) end function DatasetGriddedHDF5GetDimensionString
629)
630) ! ************************************************************************** !
631)
632) function DatasetGriddedHDF5GetNDimensions(this)
633) !
634) ! Returns the number of dimensions
635) !
636) ! Author: Glenn Hammond
637) ! Date: 10/24/11, 05/29/13
638) !
639)
640) implicit none
641)
642) class(dataset_gridded_hdf5_type) :: this
643)
644) PetscInt :: DatasetGriddedHDF5GetNDimensions
645)
646) select case(this%data_dim)
647) case(DIM_X,DIM_Y,DIM_Z)
648) DatasetGriddedHDF5GetNDimensions = ONE_INTEGER
649) case(DIM_XY,DIM_XZ,DIM_YZ)
650) DatasetGriddedHDF5GetNDimensions = TWO_INTEGER
651) case(DIM_XYZ)
652) DatasetGriddedHDF5GetNDimensions = THREE_INTEGER
653) case default
654) DatasetGriddedHDF5GetNDimensions = ZERO_INTEGER
655) end select
656)
657) end function DatasetGriddedHDF5GetNDimensions
658)
659) ! ************************************************************************** !
660)
661) subroutine DatasetGriddedHDF5InterpolateReal(this,xx,yy,zz,real_value,option)
662) !
663) ! Interpolates data from the dataset
664) !
665) ! Author: Glenn Hammond
666) ! Date: 10/26/11, 05/29/13
667) !
668)
669) use Utility_module, only : InterpolateBilinear
670) use Option_module
671)
672) implicit none
673)
674) class(dataset_gridded_hdf5_type) :: this
675) PetscReal, intent(in) :: xx, yy, zz
676) PetscReal :: real_value
677) type(option_type) :: option
678)
679) PetscInt :: spatial_interpolation_method
680) PetscInt :: i, j, k
681) PetscReal :: x, y, z
682) PetscReal :: x1, x2, y1, y2, z1
683) PetscReal :: v1, v2, v3, v4
684) PetscInt :: index
685) PetscInt :: ii, jj, kk
686) PetscInt :: i_upper, j_upper, k_upper
687) PetscReal :: dx, dy, dz
688) PetscInt :: nx, ny
689) PetscBool :: lerr
690) character(len=MAXWORDLENGTH) :: word
691)
692) call DatasetGriddedHDF5GetIndices(this,xx,yy,zz,i,j,k,x,y,z)
693)
694) ! in the below, i,j,k,xx,yy,zz to not reflect the
695) ! coordinates of the problem domain in 3D. They
696) ! are transfored to the dimensions of the dataset
697) lerr = PETSC_FALSE
698) select case(this%interpolation_method)
699) case(INTERPOLATION_STEP)
700) select case(this%data_dim)
701) case(DIM_X,DIM_Y,DIM_Z)
702) if (this%is_cell_centered) then
703) i_upper = i
704) else
705) i_upper = i+1
706) endif
707) if (i < 1 .or. i_upper > this%dims(1)) then
708) write(word,*) i
709) word = adjustl(word)
710) select case(this%data_dim)
711) case(DIM_X)
712) write(option%io_buffer,*) 'X value (', xx, &
713) ') outside of X bounds (', this%origin(1), this%extent(1), &
714) '), i ='
715) case(DIM_Y)
716) write(option%io_buffer,*) 'Y value (', yy, &
717) ') outside of Y bounds (', this%origin(2), this%extent(2), &
718) '), j ='
719) case(DIM_Z)
720) write(option%io_buffer,*) 'Z value (', zz, &
721) ') outside of Z bounds (', this%origin(3), this%extent(3), &
722) '), k ='
723) end select
724) option%io_buffer = trim(option%io_buffer) // ' ' // &
725) trim(word) // ' for gridded dataset "' // trim(this%name) // '".'
726) call printErrMsgByRank(option)
727) endif
728) index = i
729) if (.not.this%is_cell_centered) then
730) dx = this%discretization(1)
731) x1 = this%origin(1) + (i-1)*dx
732) if ((x-x1) / dx > 0.5) then
733) index = i+1
734) endif
735) endif
736) case(DIM_XY,DIM_XZ,DIM_YZ)
737) if (this%is_cell_centered) then
738) i_upper = i
739) j_upper = j
740) else
741) i_upper = i+1
742) j_upper = j+1
743) endif
744) if (i < 1 .or. i_upper > this%dims(1)) then
745) lerr = PETSC_TRUE
746) write(word,*) i
747) word = adjustl(word)
748) select case(this%data_dim)
749) case(DIM_XY,DIM_XZ)
750) write(option%io_buffer,*) 'X value (', xx, &
751) ') outside of X bounds (', this%origin(1), this%extent(1), &
752) '), i ='
753) case(DIM_YZ)
754) write(option%io_buffer,*) 'Y value (', yy, &
755) ') outside of Y bounds (', this%origin(2), this%extent(2), &
756) '), j ='
757) end select
758) option%io_buffer = trim(option%io_buffer) // ' ' // word
759) call printMsgByRank(option)
760) endif
761) if (j < 1 .or. j_upper > this%dims(2)) then
762) lerr = PETSC_TRUE
763) write(word,*) j
764) word = adjustl(word)
765) select case(this%data_dim)
766) case(DIM_XY)
767) write(option%io_buffer,*) 'Y value (', yy, &
768) ') outside of Y bounds (', this%origin(2), this%extent(2), &
769) '), j ='
770) case(DIM_YZ,DIM_XZ)
771) write(option%io_buffer,*) 'Z value (', zz, &
772) ') outside of Z bounds (', this%origin(3), this%extent(3), &
773) '), k ='
774) end select
775) option%io_buffer = trim(option%io_buffer) // ' ' // word
776) call printMsgByRank(option)
777) endif
778) if (lerr) then
779) word = this%name
780) option%io_buffer = 'gridded dataset "' // trim(word) // &
781) '" out of bounds.'
782) call printErrMsgByRank(option)
783) endif
784) ii = i
785) jj = j
786) nx = this%dims(1)
787) if (.not.this%is_cell_centered) then
788) dx = this%discretization(1)
789) dy = this%discretization(2)
790) x1 = this%origin(1) + (i-1)*dx
791) y1 = this%origin(2) + (j-1)*dy
792) if ((x-x1) / dx > 0.5d0) then
793) ii = i+1
794) endif
795) if ((y-y1) / dy > 0.5d0) then
796) jj = j+1
797) endif
798) endif
799) index = ii + (jj-1)*nx
800) case(DIM_XYZ)
801) if (this%is_cell_centered) then
802) i_upper = i
803) j_upper = j
804) k_upper = k
805) else
806) i_upper = i+1
807) j_upper = j+1
808) k_upper = k+1
809) endif
810) if (i < 1 .or. i_upper > this%dims(1)) then
811) lerr = PETSC_TRUE
812) write(word,*) i
813) word = adjustl(word)
814) write(option%io_buffer,*) 'X value (', xx, &
815) ') outside of X bounds (', this%origin(1), this%extent(1), &
816) '), i = ' // trim(word)
817) call printMsgByRank(option)
818) endif
819) if (j < 1 .or. j_upper > this%dims(2)) then
820) lerr = PETSC_TRUE
821) write(word,*) j
822) word = adjustl(word)
823) write(option%io_buffer,*) 'Y value (', yy, &
824) ') outside of Y bounds (', this%origin(2), this%extent(2), &
825) '), j = ' // trim(word)
826) call printMsgByRank(option)
827) endif
828) if (k < 1 .or. k_upper > this%dims(3)) then
829) lerr = PETSC_TRUE
830) write(word,*) k
831) word = adjustl(word)
832) write(option%io_buffer,*) 'Z value (', zz, &
833) ') outside of Z bounds (', this%origin(3), this%extent(3), &
834) '), k = ' // trim(word)
835) call printMsgByRank(option)
836) endif
837) if (lerr) then
838) word = this%name
839) option%io_buffer = 'gridded dataset "' // trim(word) // &
840) '" out of bounds.'
841) call printErrMsgByRank(option)
842) endif
843) ii = i
844) jj = j
845) kk = k
846) nx = this%dims(1)
847) ny = this%dims(2)
848) if (.not.this%is_cell_centered) then
849) dx = this%discretization(1)
850) dy = this%discretization(2)
851) dz = this%discretization(3)
852) x1 = this%origin(1) + (i-1)*dx
853) y1 = this%origin(2) + (j-1)*dy
854) z1 = this%origin(3) + (k-1)*dz
855) if ((x-x1) / dx > 0.5d0) then
856) ii = i+1
857) endif
858) if ((y-y1) / dy > 0.5d0) then
859) jj = j+1
860) endif
861) if ((z-z1) / dz > 0.5d0) then
862) kk = k+1
863) endif
864) endif
865) index = ii + (jj-1)*nx + (kk-1)*nx*ny
866) end select
867) real_value = this%rarray(index)
868) case(INTERPOLATION_LINEAR)
869) select case(this%data_dim)
870) case(DIM_X,DIM_Y,DIM_Z)
871) if (i < 1 .or. i+1 > this%dims(1)) then
872) write(word,*) i
873) word = adjustl(word)
874) select case(this%data_dim)
875) case(DIM_X)
876) write(option%io_buffer,*) 'X value (', xx, &
877) ') outside of X bounds (', this%origin(1), this%extent(1), &
878) '), i ='
879) case(DIM_Y)
880) write(option%io_buffer,*) 'Y value (', yy, &
881) ') outside of Y bounds (', this%origin(2), this%extent(2), &
882) '), j ='
883) case(DIM_Z)
884) write(option%io_buffer,*) 'Z value (', zz, &
885) ') outside of Z bounds (', this%origin(3), this%extent(3), &
886) '), k ='
887) end select
888) option%io_buffer = trim(option%io_buffer) // ' ' // &
889) trim(word) // ' for gridded dataset "' // trim(this%name) // '".'
890) call printErrMsgByRank(option)
891) endif
892) dx = this%discretization(1)
893) x1 = this%origin(1) + (i-1)*dx
894) if (this%is_cell_centered) x1 = x1 + 0.5d0*dx
895) v1 = this%rarray(i)
896) v2 = this%rarray(i+1)
897) real_value = v1 + (x-x1)/dx*(v2-v1)
898) case(DIM_XY,DIM_XZ,DIM_YZ)
899) if (i < 1 .or. i+1 > this%dims(1)) then
900) write(word,*) i
901) word = adjustl(word)
902) select case(this%data_dim)
903) case(DIM_XY,DIM_XZ)
904) write(option%io_buffer,*) 'X value (', xx, &
905) ') outside of X bounds (', this%origin(1), this%extent(1), &
906) '), i ='
907) case(DIM_YZ)
908) write(option%io_buffer,*) 'Y value (', yy, &
909) ') outside of Y bounds (', this%origin(2), this%extent(2), &
910) '), j ='
911) end select
912) option%io_buffer = trim(option%io_buffer) // ' ' // &
913) trim(word) // ' for gridded dataset "' // trim(this%name) // '".'
914) call printErrMsgByRank(option)
915) endif
916) if (j < 1 .or. j+1 > this%dims(2)) then
917) write(word,*) j
918) word = adjustl(word)
919) select case(this%data_dim)
920) case(DIM_XY)
921) write(option%io_buffer,*) 'Y value (', yy, &
922) ') outside of Y bounds (', this%origin(2), this%extent(2), &
923) '), j ='
924) case(DIM_YZ,DIM_XZ)
925) write(option%io_buffer,*) 'Z value (', zz, &
926) ') outside of Z bounds (', this%origin(3), this%extent(3), &
927) '), k ='
928) end select
929) option%io_buffer = trim(option%io_buffer) // ' ' // &
930) trim(word) // ' for gridded dataset "' // trim(this%name) // '".'
931) call printErrMsgByRank(option)
932) endif
933) dx = this%discretization(1)
934) dy = this%discretization(2)
935) nx = this%dims(1)
936)
937) x1 = this%origin(1) + (i-1)*dx
938) if (this%is_cell_centered) x1 = x1 + 0.5d0*dx
939) x2 = x1 + dx
940)
941) index = i + (j-1)*nx
942) v1 = this%rarray(index)
943) v2 = this%rarray(index+1)
944)
945) y1 = this%origin(2) + (j-1)*dy
946) if (this%is_cell_centered) y1 = y1 + 0.5d0*dy
947) y2 = y1 + dy
948)
949) ! really (j1-1+1)
950) index = i + j*nx
951) v3 = this%rarray(index)
952) v4 = this%rarray(index+1)
953)
954) real_value = InterpolateBilinear(x,y,x1,x2,y1,y2,v1,v2,v3,v4)
955) case(DIM_XYZ)
956) option%io_buffer = 'Trilinear interpolation not yet supported'
957) call printErrMsgByRank(option)
958) end select
959) end select
960)
961) end subroutine DatasetGriddedHDF5InterpolateReal
962)
963) ! ************************************************************************** !
964)
965) subroutine DatasetGriddedHDF5GetIndices(this,xx,yy,zz,i,j,k,x,y,z)
966) !
967) ! Returns bounding indices for point in dataset
968) !
969) ! Author: Glenn Hammond
970) ! Date: 10/26/11, 05/29/13
971) !
972)
973) implicit none
974)
975) class(dataset_gridded_hdf5_type) :: this
976) PetscReal, intent(in) :: xx, yy, zz
977) PetscInt :: i, j, k
978) PetscReal :: x, y, z
979)
980) PetscReal :: discretization_offset
981) PetscInt :: upper_index_offset
982) PetscReal :: tol
983) PetscReal, parameter :: tolerance_scale = 1.d-3
984)
985) select case(this%data_dim)
986) ! since these are 1D array, always use first dimension
987) case(DIM_X)
988) x = xx
989) case(DIM_Y)
990) x = yy
991) case(DIM_XY)
992) x = xx
993) y = yy
994) case(DIM_XYZ)
995) x = xx
996) y = yy
997) z = zz
998) case(DIM_Z)
999) x = zz
1000) case(DIM_XZ)
1001) x = xx
1002) y = zz
1003) case(DIM_YZ)
1004) x = yy
1005) y = zz
1006) end select
1007)
1008) upper_index_offset = 1
1009) if (this%is_cell_centered) then
1010) select case(this%interpolation_method)
1011) case(INTERPOLATION_STEP)
1012) discretization_offset = 1.d0
1013) upper_index_offset = 0
1014) case(INTERPOLATION_LINEAR)
1015) discretization_offset = 0.5d0
1016) end select
1017) i = int((x - this%origin(1))/ &
1018) this%discretization(1) + discretization_offset)
1019) ! i = max(1,min(i,this%dims(1)-1))
1020) if (this%data_dim > DIM_Z) then ! at least 2D
1021) j = int((y - this%origin(2))/ &
1022) this%discretization(2) + discretization_offset)
1023) ! j = max(1,min(j,this%dims(2)-1))
1024) endif
1025) if (this%data_dim > DIM_YZ) then ! at least 3D
1026) k = int((z - this%origin(3))/ &
1027) this%discretization(3) + discretization_offset)
1028) ! k = max(1,min(k,this%dims(3)-1))
1029) endif
1030) else
1031) select case(this%interpolation_method)
1032) case(INTERPOLATION_STEP)
1033) discretization_offset = 1.5d0
1034) case(INTERPOLATION_LINEAR)
1035) discretization_offset = 1.d0
1036) end select
1037) i = int((x - this%origin(1))/ &
1038) this%discretization(1) + discretization_offset)
1039) if (this%data_dim > DIM_Z) then ! at least 2D
1040) j = int((y - this%origin(2))/ &
1041) this%discretization(2) + discretization_offset)
1042) endif
1043) if (this%data_dim > DIM_YZ) then ! at least 3D
1044) k = int((z - this%origin(3))/ &
1045) this%discretization(3) + discretization_offset)
1046) endif
1047) endif
1048)
1049) ! if indices are out of bounds, check if on boundary and reset index
1050) !geh: the tolerance allows one to go outside the bounds
1051) if (i < 1 .or. i+upper_index_offset > this%dims(1)) then
1052) tol = this%discretization(1) * tolerance_scale
1053) if (x >= this%origin(1)-tol .and. x <= this%extent(1)+tol) then
1054) i = min(max(i,1),this%dims(1)-upper_index_offset)
1055) endif
1056) endif
1057) if (this%data_dim > DIM_Z) then ! at least 2D
1058) if (j < 1 .or. j+upper_index_offset > this%dims(2)) then
1059) tol = this%discretization(2) * tolerance_scale
1060) if (y >= this%origin(2)-tol .and. y <= this%extent(2)+tol) then
1061) j = min(max(j,1),this%dims(2)-upper_index_offset)
1062) endif
1063) endif
1064) endif
1065) if (this%data_dim > DIM_YZ) then ! at least 2D
1066) if (k < 1 .or. k+upper_index_offset > this%dims(3)) then
1067) tol = this%discretization(3) * tolerance_scale
1068) if (z >= this%origin(3)-tol .and. z <= this%extent(3)+tol) then
1069) k = min(max(k,1),this%dims(3)-upper_index_offset)
1070) endif
1071) endif
1072) endif
1073)
1074) end subroutine DatasetGriddedHDF5GetIndices
1075)
1076) ! ************************************************************************** !
1077)
1078) subroutine DatasetGriddedHDF5Print(this,option)
1079) !
1080) ! Prints dataset info
1081) !
1082) ! Author: Glenn Hammond
1083) ! Date: 10/22/13
1084) !
1085)
1086) use Option_module
1087)
1088) implicit none
1089)
1090) class(dataset_gridded_hdf5_type), target :: this
1091) type(option_type) :: option
1092)
1093) class(dataset_common_hdf5_type), pointer :: dataset_hdf5
1094)
1095) dataset_hdf5 => this
1096) call DatasetCommonHDF5Print(this,option)
1097)
1098) write(option%fid_out,'(10x,''Grid Dimension: '',a)') &
1099) trim(DatasetGriddedHDF5GetDimensionString(this))
1100) if (this%is_cell_centered) then
1101) write(option%fid_out,'(10x,''Is cell-centered?: yes'')')
1102) else
1103) write(option%fid_out,'(10x,''Is cell-centered?: no'')')
1104) endif
1105) write(option%fid_out,'(10x,''Origin: '',3es12.4)') this%origin(:)
1106) write(option%fid_out,'(10x,''Discretization: '',3es12.4)') &
1107) this%discretization(:)
1108)
1109) end subroutine DatasetGriddedHDF5Print
1110)
1111) ! ************************************************************************** !
1112)
1113) subroutine DatasetGriddedHDF5Strip(this)
1114) !
1115) ! Strips allocated objects within XYZ dataset object
1116) !
1117) ! Author: Glenn Hammond
1118) ! Date: 05/03/13
1119) !
1120)
1121) use Utility_module, only : DeallocateArray
1122)
1123) implicit none
1124)
1125) class(dataset_gridded_hdf5_type) :: this
1126)
1127) call DatasetCommonHDF5Strip(this)
1128)
1129) call DeallocateArray(this%origin)
1130) call DeallocateArray(this%extent)
1131) call DeallocateArray(this%discretization)
1132)
1133) end subroutine DatasetGriddedHDF5Strip
1134)
1135) ! ************************************************************************** !
1136)
1137) subroutine DatasetGriddedHDF5Destroy(this)
1138) !
1139) ! Destroys a dataset
1140) !
1141) ! Author: Glenn Hammond
1142) ! Date: 05/29/13
1143) !
1144)
1145) implicit none
1146)
1147) class(dataset_gridded_hdf5_type), pointer :: this
1148)
1149) if (.not.associated(this)) return
1150)
1151) call DatasetGriddedHDF5Strip(this)
1152)
1153) deallocate(this)
1154) nullify(this)
1155)
1156) end subroutine DatasetGriddedHDF5Destroy
1157)
1158) end module Dataset_Gridded_HDF5_class