Opens a hdf5 file given by name with mode stated by flag. This can take one of these values: 'R' for read only 'W' for write only (will delete existing file) 'RW' for read-write
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(hdf5_obj), | intent(inout) | :: | this |
A HDF5 object |
||
| character(len=*), | intent(in) | :: | name |
Filename |
||
| character(len=*), | intent(in) | :: | flag |
Access mode |
impure subroutine hdf5_obj_Open(this,name,flag) !> Opens a hdf5 file given by name with mode stated by flag. ! This can take one of these values: ! 'R' for read only ! 'W' for write only (will delete existing file) ! 'RW' for read-write use mpi_f08 implicit none class(hdf5_obj), intent(inout) :: this !! A HDF5 object character(len=*), intent(in) :: name !! Filename character(len=*), intent(in) :: flag !! Access mode ! Work variables type(stringtool_obj) :: stringtool type(sysutils_obj) :: sysutils logical :: file_exists integer :: ierr ! Save filename this%filename = stringtool%RemoveExtension(name)//HDF5_EXTENSION ! Check if file is already open if (this%is_open) call this%parallel%Stop("Unable to open file while another one is open: "//this%filename) associate(plistid => this%plistid, fid => this%fid, parallel => this%parallel) ! MPI stuff call H5Pcreate_f(H5P_FILE_ACCESS_F, plistid, ierr) ! Add MPI IO communicator to the file access property list call H5Pset_fapl_mpio_f(plistid, parallel%comm%g%MPI_VAL, MPI_INFO_NULL%MPI_VAL, ierr) ! Instruct HDF5 to not do metadata operations collectively ! If encountering slowdowns, set to true. ! This will cause MPI rank=0 to read and broadcast to ! others, avoiding a read rush in large MPI jobs. !call H5Pset_all_coll_metadata_ops_f(plistid,.false.,ierr) call H5Pset_all_coll_metadata_ops_f(plistid,.true.,ierr) select case(trim(adjustl(flag))) case ('W','w' ) call H5Fcreate_f(this%filename, H5F_ACC_TRUNC_F, fid, ierr, access_prp=plistid) if (ierr.ne.0) & call parallel%Stop("Unable to create new HDF5 file: "//this%filename) case ('R','r' ) call H5Fopen_f(this%filename, H5F_ACC_RDONLY_F, fid, ierr, access_prp=plistid) if (ierr.ne.0) & call parallel%Stop("Unable to open existing HDF5 file in Read mode: "//this%filename) case ('RW','rw') if (this%parallel%RankIsRoot()) & file_exists = sysutils%FileExists(this%filename) call this%parallel%bcast(file_exists) if (file_exists) then call H5Fopen_f(this%filename, H5F_ACC_RDWR_F, fid, ierr, access_prp=plistid) if (ierr.ne.0) & call parallel%Stop("Unable to open existing HDF5 file in Read/Write mode: "//this%filename) else call H5Fcreate_f(this%filename, H5F_ACC_TRUNC_F, fid, ierr, access_prp=plistid) if (ierr.ne.0) & call parallel%Stop("Unable to create new HDF5 file: "//this%filename) end if case default call parallel%Stop("Invalide acces mode for HDF5 file: "//this%filename) end select end associate ! Toggle switch this%is_open = .true. return end subroutine hdf5_obj_Open