hdf5_obj_Open Subroutine

private impure subroutine hdf5_obj_Open(this, name, flag)

Uses

  • proc~~hdf5_obj_open~~UsesGraph proc~hdf5_obj_open hdf5_obj%hdf5_obj_Open mpi_f08 mpi_f08 proc~hdf5_obj_open->mpi_f08

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 Bound

hdf5_obj

Arguments

Type IntentOptional Attributes Name
class(hdf5_obj), intent(inout) :: this

A HDF5 object

character(len=*), intent(in) :: name

Filename

character(len=*), intent(in) :: flag

Access mode


Calls

proc~~hdf5_obj_open~~CallsGraph proc~hdf5_obj_open hdf5_obj%hdf5_obj_Open h5fcreate_f h5fcreate_f proc~hdf5_obj_open->h5fcreate_f h5fopen_f h5fopen_f proc~hdf5_obj_open->h5fopen_f h5pcreate_f h5pcreate_f proc~hdf5_obj_open->h5pcreate_f h5pset_all_coll_metadata_ops_f h5pset_all_coll_metadata_ops_f proc~hdf5_obj_open->h5pset_all_coll_metadata_ops_f h5pset_fapl_mpio_f h5pset_fapl_mpio_f proc~hdf5_obj_open->h5pset_fapl_mpio_f none~bcast parallel_obj%Bcast proc~hdf5_obj_open->none~bcast proc~parallel_obj_rankisroot parallel_obj%parallel_obj_RankIsRoot proc~hdf5_obj_open->proc~parallel_obj_rankisroot proc~stringtool_obj_removeextension stringtool_obj%stringtool_obj_RemoveExtension proc~hdf5_obj_open->proc~stringtool_obj_removeextension proc~sysutils_obj_fileexists sysutils_obj%sysutils_obj_FileExists proc~hdf5_obj_open->proc~sysutils_obj_fileexists proc~parallel_obj_bcast_char parallel_obj%parallel_obj_Bcast_char none~bcast->proc~parallel_obj_bcast_char proc~parallel_obj_bcast_int_0d parallel_obj%parallel_obj_Bcast_int_0d none~bcast->proc~parallel_obj_bcast_int_0d proc~parallel_obj_bcast_int_1d parallel_obj%parallel_obj_Bcast_int_1d none~bcast->proc~parallel_obj_bcast_int_1d proc~parallel_obj_bcast_int_2d parallel_obj%parallel_obj_Bcast_int_2d none~bcast->proc~parallel_obj_bcast_int_2d proc~parallel_obj_bcast_int_3d parallel_obj%parallel_obj_Bcast_int_3d none~bcast->proc~parallel_obj_bcast_int_3d proc~parallel_obj_bcast_log_0d parallel_obj%parallel_obj_Bcast_log_0d none~bcast->proc~parallel_obj_bcast_log_0d proc~parallel_obj_bcast_log_1d parallel_obj%parallel_obj_Bcast_log_1d none~bcast->proc~parallel_obj_bcast_log_1d proc~parallel_obj_bcast_real_0d parallel_obj%parallel_obj_Bcast_real_0d none~bcast->proc~parallel_obj_bcast_real_0d proc~parallel_obj_bcast_real_1d parallel_obj%parallel_obj_Bcast_real_1d none~bcast->proc~parallel_obj_bcast_real_1d proc~parallel_obj_bcast_real_2d parallel_obj%parallel_obj_Bcast_real_2d none~bcast->proc~parallel_obj_bcast_real_2d proc~parallel_obj_bcast_real_3d parallel_obj%parallel_obj_Bcast_real_3d none~bcast->proc~parallel_obj_bcast_real_3d mpi_bcast mpi_bcast proc~parallel_obj_bcast_char->mpi_bcast proc~parallel_obj_bcast_int_0d->mpi_bcast proc~parallel_obj_bcast_int_1d->mpi_bcast proc~parallel_obj_bcast_int_2d->mpi_bcast proc~parallel_obj_bcast_int_3d->mpi_bcast proc~parallel_obj_bcast_log_0d->mpi_bcast proc~parallel_obj_bcast_log_1d->mpi_bcast proc~parallel_obj_bcast_real_0d->mpi_bcast proc~parallel_obj_bcast_real_1d->mpi_bcast proc~parallel_obj_bcast_real_2d->mpi_bcast proc~parallel_obj_bcast_real_3d->mpi_bcast

Source Code

    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