Writes Eulerian/3D data to a HDF5 file.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(hdf5_obj), | intent(inout) | :: | this |
A HDF5 object |
||
| character(len=*), | intent(in) | :: | groupname |
Groupname |
||
| character(len=*), | intent(in) | :: | name |
Variable name |
||
| class(*), | intent(in), | target | :: | array(:,:,:) |
3-D data array |
|
| integer, | intent(in) | :: | lo(3) |
Low bounds |
||
| integer, | intent(in) | :: | hi(3) |
High bounds |
impure subroutine hdf5_obj_Write3D(this,groupname,name,array,lo,hi) !> Writes Eulerian/3D data to a HDF5 file. implicit none class(hdf5_obj), intent(inout) :: this !! A HDF5 object character(len=*), intent(in) :: groupname !! Groupname character(len=*), intent(in) :: name !! Variable name class(*), intent(in), & target :: array(:,:,:) !! 3-D data array integer, intent(in) :: lo(3) !! Low bounds integer, intent(in) :: hi(3) !! High bounds ! Work variables integer :: ierr integer(HID_T) :: sidc,sidg integer(HID_T) :: did integer(HID_T) :: plist integer(HID_T) :: dtype integer(HSIZE_T) :: dims(3) integer(HSIZE_T) :: offset(3) integer :: ndim integer :: maxhi(3) integer :: minlo(3) integer(HID_T) :: obj integer :: npoint type(c_ptr) :: fptr = c_null_ptr ! Get number of grid points in each direction call this%parallel%max(hi,maxhi) call this%parallel%min(lo,minlo) ! Check for invalid arrays npoint = sum(hi-lo+1) if (npoint.lt.0) then dims = int(0,kind=HSIZE_T) else dims = int(hi-lo+1,kind=HSIZE_T) end if ! Open group, if not already open call this%OpenGroup(groupname) ! Find parent object obj = this%GetGroupObject(groupname) ! writing 3d arrays ndim = 3 ! Create the global data space call H5Screate_simple_f(ndim,int(maxhi-minlo+1,kind=HSIZE_T),sidg,ierr) ! Create the local data space call H5Screate_simple_f(ndim,dims,sidc,ierr) ! Select a subset of the global data space offset = int(lo - minlo,kind=HSIZE_T) call H5Sselect_hyperslab_f(sidg,H5S_SELECT_SET_F,offset,dims,ierr ) ! Select all elements from the local data space offset = int(0,kind=HSIZE_T) call H5Sselect_hyperslab_f(sidc,H5S_SELECT_SET_F,offset,dims,ierr ) ! Create data space select type(array) type is (real(leapDP)) dtype = H5T_NATIVE_DOUBLE type is (real(leapSP)) dtype = H5T_NATIVE_REAL type is (integer(leapI4)) dtype = H5T_NATIVE_INTEGER type is (integer(leapI8)) dtype = H5T_STD_I64LE end select call H5Dcreate_f(obj, trim(adjustl(name)), dtype, sidg, did, ierr) ! Create property list call H5Pcreate_f(H5P_DATASET_XFER_F, plist, ierr ) call H5Pset_dxpl_mpio_f(plist,H5FD_MPIO_COLLECTIVE_F,ierr ) if (size(array).ne.0) then select type(array) type is (real(leapDP)) fptr = c_loc(array(1,1,1)) type is (real(leapSP)) fptr = c_loc(array(1,1,1)) type is (integer(leapI4)) fptr = c_loc(array(1,1,1)) type is (integer(leapI8)) fptr = c_loc(array(1,1,1)) end select end if ! Write data set call H5Dwrite_f(did,dtype,fptr,ierr,mem_space_id=sidc,file_space_id=sidg,xfer_prp=plist) ! Deallocate pointer fptr = c_null_ptr ! Close data set call H5Dclose_f(did,ierr) ! Close data spaces call H5Sclose_f(sidc,ierr) call H5Sclose_f(sidg,ierr) ! Close property list call H5Pclose_f(plist,ierr) ! Close group call this%CloseGroup(groupname) return end subroutine hdf5_obj_Write3D