!!******************************************************************************
!!
!!  This file is part of the AMUN source code, a program to perform
!!  Newtonian or relativistic magnetohydrodynamical simulations on uniform or
!!  adaptive mesh.
!!
!!  Copyright (C) 2008-2019 Grzegorz Kowal <grzegorz@amuncode.org>
!!
!!  This program is free software: you can redistribute it and/or modify
!!  it under the terms of the GNU General Public License as published by
!!  the Free Software Foundation, either version 3 of the License, or
!!  (at your option) any later version.
!!
!!  This program is distributed in the hope that it will be useful,
!!  but WITHOUT ANY WARRANTY; without even the implied warranty of
!!  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
!!  GNU General Public License for more details.
!!
!!  You should have received a copy of the GNU General Public License
!!  along with this program.  If not, see <http://www.gnu.org/licenses/>.
!!
!!******************************************************************************
!!
!! module: IO
!!
!!  This module handles data storage and job restart from restart files.
!!
!!
!!******************************************************************************
!
module io

! import external subroutines
!
  use blocks, only : pointer_meta
#ifdef HDF5
  use hdf5  , only : hid_t
#endif /* HDF5 */
  use timers, only : set_timer, start_timer, stop_timer

! module variables are not implicit by default
!
  implicit none

! subroutine interfaces
!
  interface read_snapshot_parameter
#ifdef HDF5
    module procedure read_snapshot_parameter_string_h5
    module procedure read_snapshot_parameter_integer_h5
    module procedure read_snapshot_parameter_integer_vector_h5
    module procedure read_snapshot_parameter_double_h5
#endif /* HDF5 */
  end interface
  interface write_attribute
#ifdef HDF5
    module procedure write_scalar_attribute_string_h5
    module procedure write_scalar_attribute_integer_h5
    module procedure write_scalar_attribute_double_h5
    module procedure write_vector_attribute_integer_h5
    module procedure write_vector_attribute_double_h5
#endif /* HDF5 */
  end interface
  interface read_attribute
#ifdef HDF5
    module procedure read_scalar_attribute_integer_h5
    module procedure read_scalar_attribute_double_h5
    module procedure read_vector_attribute_integer_h5
    module procedure read_vector_attribute_double_h5
#endif /* HDF5 */
  end interface
  interface write_array
#ifdef HDF5
    module procedure write_1d_array_integer_h5
    module procedure write_2d_array_integer_h5
    module procedure write_3d_array_integer_h5
    module procedure write_4d_array_integer_h5
    module procedure write_5d_array_integer_h5
    module procedure write_1d_array_double_h5
    module procedure write_2d_array_double_h5
    module procedure write_3d_array_double_h5
    module procedure write_4d_array_double_h5
    module procedure write_5d_array_double_h5
#endif /* HDF5 */
  end interface
  interface read_array
#ifdef HDF5
    module procedure read_1d_array_integer_h5
    module procedure read_2d_array_integer_h5
    module procedure read_3d_array_integer_h5
    module procedure read_4d_array_integer_h5
    module procedure read_5d_array_integer_h5
    module procedure read_1d_array_double_h5
    module procedure read_2d_array_double_h5
    module procedure read_3d_array_double_h5
    module procedure read_4d_array_double_h5
    module procedure read_5d_array_double_h5
#endif /* HDF5 */
  end interface

! timer indices
!
  integer            , save :: iio
#ifdef PROFILE
  integer            , save :: ioi, iow, ios
#endif /* PROFILE */

! MODULE PARAMETERS:
! =================
!
!   respath - the directory from which the restart snapshots should be read;
!   ftype   - the type of snapshots to write:
!        'p' -> all primitive variables (default);
!        'c' -> all conserved variables;
!   nrest   - for job restarting, this is the number of restart snapshot;
!   irest   - the local counter for the restart snapshots;
!   isnap   - the local counter for the regular snapshots;
!   ishift  - the shift of the snapshot counter for restarting job with
!             different snapshot interval;
!   hrest   - the execution time interval for restart snapshot storing
!             (in hours); the minimum allowed value is 3 minutes;
!   hsnap   - the problem time interval for regular snapshot storing;
!   tsnap   - the next snapshot time;
!
  character(len=255), save :: respath     = "./"
  character         , save :: ftype       = "p"
  character(len=64) , save :: ftype_name  = "primitive"
  integer           , save :: nrest       = -1
  integer(kind=4)   , save :: irest       =  1
  integer(kind=4)   , save :: isnap       =  0
  integer(kind=4)   , save :: ishift      =  0
  real(kind=8)      , save :: hrest       =  6.0d+00
  real(kind=8)      , save :: hsnap       =  1.0d+00
  real(kind=8)      , save :: tsnap       =  0.0d+00

! flag indicating to store snapshots at exact intervals
!
  logical           , save :: precise_snapshots = .false.

! flags to determine the way of data writing
!
  logical           , save :: with_ghosts = .true.

! a flag to determine if XDMF files should be generated
!
  logical           , save :: with_xdmf   = .false.

#ifdef HDF5
! compression type
!
  integer      , parameter :: H5Z_DEFLATE = 1, H5Z_ZSTANDARD = 32015

! compression type (0 for no compressions, 1 for deflate, 32015 for zstandard)
!
  integer           , save :: compression = 0

! compression level
!
  integer           , save :: clevel      = 0

! HDF5 property object identifier
!
  integer(hid_t)    , save :: pid
#endif /* HDF5 */

! local variables to store the number of processors
!
  integer(kind=4)   , save :: nfiles = 1

! array of pointer used during job restart
!
  type(pointer_meta), dimension(:), allocatable, save :: block_array

! by default everything is private
!
  private

! declare public subroutines
!
  public :: initialize_io, finalize_io, print_io
  public :: restart_snapshot_number, restart_from_snapshot
  public :: read_snapshot_parameter
  public :: read_restart_snapshot, write_restart_snapshot, write_snapshot
  public :: next_tout, precise_snapshots

!- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
!
  contains
!
!===============================================================================
!!
!!***  PUBLIC SUBROUTINES  *****************************************************
!!
!===============================================================================
!
!===============================================================================
!
! subroutine INITIALIZE_IO:
! ------------------------
!
!   Subroutine initializes module IO by setting its parameters.
!
!   Arguments:
!
!     verbose - flag determining if the subroutine should be verbose;
!     iret    - return flag of the procedure execution status;
!
!===============================================================================
!
  subroutine initialize_io(verbose, iret)

! import external procedures
!
#ifdef HDF5
    use hdf5           , only : hsize_t
    use hdf5           , only : H5P_DATASET_CREATE_F, H5Z_FLAG_OPTIONAL_F
    use hdf5           , only : h5open_f, h5zfilter_avail_f, h5pcreate_f
    use hdf5           , only : h5pset_deflate_f, h5pset_filter_f
#endif /* HDF5 */
    use iso_fortran_env, only : error_unit
    use parameters     , only : get_parameter

! local variables are not implicit by default
!
    implicit none

! subroutine arguments
!
    logical, intent(in)    :: verbose
    integer, intent(inout) :: iret

! local variables
!
    character(len=255)    :: precise = "off"
    character(len=255)    :: ghosts  = "on"
    character(len=255)    :: xdmf    = "off"
#ifdef HDF5
    logical               :: status = .false.
    integer(hsize_t)      :: cd_nelmts = 1
    integer, dimension(1) :: cd_values = 3
#endif /* HDF5 */

! local parameters
!
    character(len=*), parameter :: loc = 'IO::initialize_io()'
!
!-------------------------------------------------------------------------------
!
! set timer descriptions
!
    call set_timer('SNAPSHOTS I/O'        , iio)
#ifdef PROFILE
    call set_timer('io:: initialization'  , ioi)
    call set_timer('io:: snapshot writing', iow)
    call set_timer('io:: snapshot reading', ios)

! start accounting time for module initialization/finalization
!
    call start_timer(ioi)
#endif /* PROFILE */

! get module parameters
!
    call get_parameter("restart_path"     , respath)
    call get_parameter("restart_number"   , nrest  )
    call get_parameter("restart_interval" , hrest  )
    call get_parameter("snapshot_type"    , ftype  )
    call get_parameter("snapshot_interval", hsnap  )
    call get_parameter("precise_snapshots", precise)
    call get_parameter("include_ghosts"   , ghosts )
    call get_parameter("generate_xdmf"    , xdmf   )

! check the snapshot type
!
    select case(ftype)
    case('c')
      ftype_name = 'conservative variables'
    case('p')
      ftype_name = 'primitive variables'
    case default
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Snapshot type " // trim(ftype) // " is not suppoerted!"
      iret = 1
    end select

! check ghost cell storing flag
!
    select case(trim(precise))
    case ("off", "OFF", "n", "N", "false", "FALSE", "no", "NO")
      precise_snapshots = .false.
    case default
      precise_snapshots = .true.
    end select

! check ghost cell storing flag
!
    select case(trim(ghosts))
    case ("off", "OFF", "n", "N", "false", "FALSE", "no", "NO")
      with_ghosts = .false.
    case default
      with_ghosts = .true.
    end select

! check flag for generating XDMF files
!
    select case(trim(xdmf))
    case ("off", "OFF", "n", "N", "false", "FALSE", "no", "NO")
      with_xdmf = .false.
    case default
      with_xdmf = .true.
    end select

#ifdef HDF5
! initialize the FORTRAN interface
!
    call h5open_f(iret)

! in the case of error, print a message and quit the subroutine
!
    if (iret < 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot initialize the HDF5 Fortran interface!"
      return
    end if

! prepare the property object for compression
!
    call h5pcreate_f(H5P_DATASET_CREATE_F, pid, iret)

! check if the object has been created properly, if not quit
!
    if (iret < 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create the compression property for datasets!"
      return
    end if

! detect available compressions
!
    status = .false.
    if (.not. status) then
      call h5zfilter_avail_f(H5Z_ZSTANDARD, status, iret)
      if (status) compression = H5Z_ZSTANDARD
    end if
    if (.not. status) then
      call h5zfilter_avail_f(H5Z_DEFLATE, status, iret)
      if (status) compression = H5Z_DEFLATE
    end if

! get compression_level
!
    call get_parameter("compression_level", clevel)

! initialize proper compressor
!
    select case(compression)
    case(H5Z_ZSTANDARD)
      if (clevel < 1 .or. clevel > 20) clevel = 3
      cd_values(:) = clevel
      call h5pset_filter_f(pid, H5Z_ZSTANDARD, H5Z_FLAG_OPTIONAL_F, cd_nelmts  &
                         , cd_values, iret)
    case(H5Z_DEFLATE)
      if (clevel < 1 .or. clevel  > 9) clevel = 6
      call h5pset_deflate_f(pid, clevel, iret)
    case default
    end select
#endif /* HDF5 */

#ifdef PROFILE
! stop accounting time for module initialization/finalization
!
    call stop_timer(ioi)
#endif /* PROFILE */

!-------------------------------------------------------------------------------
!
  end subroutine initialize_io
!
!===============================================================================
!
! subroutine FINALIZE_IO:
! ----------------------
!
!   Subroutine releases memory used by the module.
!
!   Arguments:
!
!     iret    - an integer flag for error return value;
!
!===============================================================================
!
  subroutine finalize_io(iret)

! import external procedures
!
#ifdef HDF5
    use hdf5           , only : h5pclose_f, h5close_f
    use iso_fortran_env, only : error_unit
#endif /* HDF5 */

! local variables are not implicit by default
!
    implicit none

! subroutine arguments
!
    integer, intent(inout) :: iret

! local parameters
!
    character(len=*), parameter :: loc = 'IO::finalize_io()'
!
!-------------------------------------------------------------------------------
!
#ifdef PROFILE
! start accounting time for module initialization/finalization
!
    call start_timer(ioi)
#endif /* PROFILE */

#ifdef HDF5
! close the property object for compression
!
    call h5pclose_f(pid, iret)

! check if the object has been closed properly
!
    if (iret < 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot close the compression property for datasets!"
      return
    end if

! close the FORTRAN interface
!
    call h5close_f(iret)

! check if the interface has been closed successfuly
!
    if (iret > 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot close the HDF5 Fortran interface!"
      return
    end if
#endif /* HDF5 */

#ifdef PROFILE
! stop accounting time for module initialization/finalization
!
    call stop_timer(ioi)
#endif /* PROFILE */

!-------------------------------------------------------------------------------
!
  end subroutine finalize_io
!
!===============================================================================
!
! subroutine PRINT_IO:
! -------------------
!
!   Subroutine prints IO parameters.
!
!   Arguments:
!
!     verbose - flag determining if the subroutine should be verbose;
!
!===============================================================================
!
  subroutine print_io(verbose)

! import external procedures and variables
!
    use helpers, only : print_section, print_parameter

! local variables are not implicit by default
!
    implicit none

! subroutine arguments
!
    logical, intent(in) :: verbose

! local variables
!
    character(len=80) :: sfmt, msg
    integer           :: dd, hh, mm, ss
!
!-------------------------------------------------------------------------------
!
    if (verbose) then
      call print_section(verbose, "Snapshots")
      if (precise_snapshots) then
        call print_parameter(verbose, "precise snapshot intervals", "on" )
      else
        call print_parameter(verbose, "precise snapshot intervals", "off")
      end if
      call print_parameter(verbose, "snapshot type", ftype_name)
      if (with_ghosts) then
        call print_parameter(verbose, "with ghosts cells", "on" )
      else
        call print_parameter(verbose, "with ghosts cells", "off")
      end if
#ifdef HDF5
      select case(compression)
      case(H5Z_ZSTANDARD)
        call print_parameter(verbose, "HDF5 compression" , "zstd"   )
        call print_parameter(verbose, "compression level", clevel   )
      case(H5Z_DEFLATE)
        call print_parameter(verbose, "HDF5 compression" , "deflate")
        call print_parameter(verbose, "compression level", clevel   )
      case default
        call print_parameter(verbose, "HDF5 compression" , "none"   )
      end select
#endif /* HDF5 */
      if (with_xdmf) then
        call print_parameter(verbose, "generate XDMF files", "on" )
      else
        call print_parameter(verbose, "generate XDMF files", "off")
      end if
        call print_parameter(verbose, "snapshot interval"  , hsnap)
      if (hrest > 0.0d+00) then
        dd = int(hrest / 2.4d+01)
        hh = int(mod(hrest, 2.4d+01))
        mm = int(mod(6.0d+01 * hrest, 6.0d+01))
        ss = int(mod(3.6d+03 * hrest, 6.0d+01))
        sfmt = "(i2.2,'d',i2.2,'h',i2.2,'m',i2.2,'s')"
        write(msg,sfmt) dd, hh, mm, ss
        call print_parameter(verbose, "restart interval"     , msg    )
      end if
      if (restart_from_snapshot()) then
        call print_parameter(verbose, "restart from path"    , respath)
        call print_parameter(verbose, "restart from snapshot", nrest  )
      end if
    end if

!-------------------------------------------------------------------------------
!
  end subroutine print_io
!
!===============================================================================
!
! function RESTART_SNAPSHOT_NUMBER:
! --------------------------------
!
!   Subroutine returns the number of restart snapshot.
!
!
!===============================================================================
!
  integer function restart_snapshot_number()

! local variables are not implicit by default
!
    implicit none
!
!-------------------------------------------------------------------------------
!
    restart_snapshot_number = nrest

!-------------------------------------------------------------------------------
!
  end function restart_snapshot_number
!
!===============================================================================
!
! function RESTART_FROM_SNAPSHOT:
! ------------------------------
!
!   Subroutine returns true if the current job is the restarted one.
!
!
!===============================================================================
!
  logical function restart_from_snapshot()

! local variables are not implicit by default
!
    implicit none
!
!-------------------------------------------------------------------------------
!
    restart_from_snapshot = (nrest > 0)

!-------------------------------------------------------------------------------
!
  end function restart_from_snapshot
!
!===============================================================================
!
! subroutine READ_RESTART_SNAPSHOT:
! --------------------------------
!
!   Subroutine reads restart snapshot files in order to resume the job.
!   This is a wrapper calling specific format subroutine.
!
!   Arguments:
!
!     iret - the return flag to inform if subroutine succeeded or failed;
!
!===============================================================================
!
  subroutine read_restart_snapshot(iret)

! import external variables
!
    use evolution      , only : time

! local variables are not implicit by default
!
    implicit none

! input and output arguments
!
    integer, intent(out) :: iret
!
!-------------------------------------------------------------------------------
!
#ifdef PROFILE
! start accounting time for the data reading
!
    call start_timer(ios)
#endif /* PROFILE */

! reset the return flag
!
    iret = 0

! start accounting time for I/O
!
    call start_timer(iio)

#ifdef HDF5
! read HDF5 restart file and rebuild the meta and data block structures
!
    call read_restart_snapshot_h5(iret)
#endif /* HDF5 */

! stop accounting time for I/O
!
    call stop_timer(iio)

! calculate the shift of the snapshot counter, and the next snapshot time
!
    ishift = int(time / hsnap) - isnap + 1
    tsnap  = (ishift + isnap) * hsnap

#ifdef PROFILE
! stop accounting time for the data reading
!
    call stop_timer(ios)
#endif /* PROFILE */

!-------------------------------------------------------------------------------
!
  end subroutine read_restart_snapshot
!
!===============================================================================
!
! subroutine WRITE_RESTART_SNAPSHOTS:
! ----------------------------------
!
!   Subroutine stores current restart snapshot files.  This is a wrapper
!   calling specific format subroutine.
!
!   Arguments:
!
!     thrs - the current execution time in hours;
!     nrun - the run number;
!     iret - the return flag;
!
!===============================================================================
!
  subroutine write_restart_snapshot(thrs, nrun, iret)

! local variables are not implicit by default
!
    implicit none

! input and output arguments
!
    real(kind=8), intent(in)  :: thrs
    integer     , intent(in)  :: nrun
    integer     , intent(out) :: iret
!
!-------------------------------------------------------------------------------
!
! reset the return flag
!
    iret = 0

! check if conditions for storing the restart snapshot have been met
!
    if (hrest < 5.0d-02 .or. thrs < irest * hrest) return

#ifdef PROFILE
! start accounting time for the data writing
!
    call start_timer(iow)
#endif /* PROFILE */

! start accounting time for I/O
!
    call start_timer(iio)

#ifdef HDF5
! store restart file
!
    call write_restart_snapshot_h5(nrun, iret)
#endif /* HDF5 */

! stop accounting time for I/O
!
    call stop_timer(iio)

! increase the restart snapshot counter
!
    irest = irest + 1

#ifdef PROFILE
! stop accounting time for the data writing
!
    call stop_timer(iow)
#endif /* PROFILE */

!-------------------------------------------------------------------------------
!
  end subroutine write_restart_snapshot
!
!===============================================================================
!
! subroutine WRITE_SNAPSHOT:
! -------------------------
!
!   Subroutine stores block data in snapshots.  Block variables are grouped
!   together and stored in big 4D arrays separately.  This is a wrapper for
!   specific format storing.
!
!
!===============================================================================
!
  subroutine write_snapshot()

! import external variables
!
    use evolution      , only : time
    use mpitools       , only : master

! local variables are not implicit by default
!
    implicit none
!
!-------------------------------------------------------------------------------
!
! check if conditions for storing the regular snapshot have been met
!
    if (hsnap <= 0.0e+00 .or. time < tsnap) return

#ifdef PROFILE
! start accounting time for the data writing
!
    call start_timer(iow)
#endif /* PROFILE */

! start accounting time for I/O
!
    call start_timer(iio)

#ifdef HDF5
! store variable snapshot file
!
    call write_snapshot_h5()
    if (with_xdmf) then
      call write_snapshot_xdmf()
      if (master) call write_snapshot_xdmf_master()
    end if
#endif /* HDF5 */

! stop accounting time for I/O
!
    call stop_timer(iio)

! increase the snapshot counter and calculate the next snapshot time
!
    isnap = isnap + 1
    tsnap = (ishift + isnap) * hsnap

#ifdef PROFILE
! stop accounting time for the data writing
!
    call stop_timer(iow)
#endif /* PROFILE */

!-------------------------------------------------------------------------------
!
  end subroutine write_snapshot
!
!===============================================================================
!
! function NEXT_TOUT:
! ------------------
!
!   Function returns the next data snapshot time.
!
!
!===============================================================================
!
  real(kind=8) function next_tout()

! local variables are not implicit by default
!
    implicit none
!
!-------------------------------------------------------------------------------
!
    if (hsnap > 0.0d+00) then
      next_tout = tsnap
    else
      next_tout = huge(hsnap)
    end if

!-------------------------------------------------------------------------------
!
  end function next_tout
!
!===============================================================================
!!
!!***  PRIVATE SUBROUTINES  ****************************************************
!!
!===============================================================================
!
#ifdef HDF5
!
!===============================================================================
!
! subroutine READ_SNAPSHOT_PARAMETER_STRING_H5:
! --------------------------------------------
!
!   Subroutine reads a string parameter from the restart snapshot.
!
!   Arguments:
!
!     pname   - the parameter name;
!     pvalue  - the parameter value;
!     iret    - the success flag (the success is 0, failure otherwise);
!
!===============================================================================
!
  subroutine read_snapshot_parameter_string_h5(pname, pvalue, iret)

! import external procedures
!
    use hdf5           , only : hid_t
    use hdf5           , only : H5F_ACC_RDONLY_F
    use hdf5           , only : H5T_NATIVE_CHARACTER
    use hdf5           , only : hid_t, hsize_t, size_t
    use hdf5           , only : h5fopen_f, h5fclose_f
    use hdf5           , only : h5gopen_f, h5gclose_f
    use hdf5           , only : h5aexists_by_name_f
    use hdf5           , only : h5tcopy_f, h5tset_size_f
    use hdf5           , only : h5aopen_by_name_f, h5aread_f, h5aclose_f
    use iso_fortran_env, only : error_unit
    use mpitools       , only : nproc
    use parameters     , only : get_parameter

! local variables are not implicit by default
!
    implicit none

! subroutine arguments
!
    character(len=*), intent(in)    :: pname
    character(len=*), intent(inout) :: pvalue
    integer         , intent(inout) :: iret

! local variables
!
    logical            :: info
    character(len=255) :: rname
    integer            :: np
    integer(hid_t)     :: fid, gid, tid, aid
    integer(size_t)    :: aln
    integer(hsize_t)   :: am(1) = 1

! local parameters
!
    character(len=*), parameter :: loc =                                       &
                                     'IO::read_snapshot_parameter_string_h5()'
!
!-------------------------------------------------------------------------------
!
! reset the success flag
!
    iret = 0

! generate the filename of the restart snapshot
!
    info  = .false.
    np    = nproc + 1
    do while (.not. info .and. np >= 0)
      np = np - 1
      write (rname, "(a,'r',i6.6,'_',i5.5,'.h5')") trim(respath), nrest, np
      inquire(file = rname, exist = info)
    end do

! procees if file exists
!
    if (info) then
      call h5fopen_f(rname, H5F_ACC_RDONLY_F, fid, iret)
      if (iret >= 0) then
        call h5gopen_f(fid, 'attributes', gid, iret)
        if (iret >= 0) then
          call h5aexists_by_name_f(gid, '.', trim(pname), info, iret)
          if (info .and. iret >= 0) then
            call h5aopen_by_name_f(gid, '.', trim(pname), aid, iret)
            if (iret >= 0) then
              aln = len(pvalue)
              call h5tcopy_f(H5T_NATIVE_CHARACTER, tid, iret)
              call h5tset_size_f(tid, aln, iret)
              call h5aread_f(aid, tid, pvalue, am(:), iret)
              call h5aclose_f(aid, iret)
            end if
          end if
          call h5gclose_f(gid, iret)
        end if
        call h5fclose_f(fid, iret)
      end if
    else
      write(error_unit,"('[', a, ']: ', a)") trim(loc)                         &
                       , "Snapshot " // trim(rname) // " file does not exist!"
      iret = 1
    end if

!-------------------------------------------------------------------------------
!
  end subroutine read_snapshot_parameter_string_h5
!
!===============================================================================
!
! subroutine READ_SNAPSHOT_PARAMETER_INTEGER_H5:
! ---------------------------------------------
!
!   Subroutine reads an integer parameter from the restart snapshot.
!
!   Arguments:
!
!     pname   - the parameter name;
!     pvalue  - the parameter value;
!     iret    - the success flag (the success is 0, failure otherwise);
!
!===============================================================================
!
  subroutine read_snapshot_parameter_integer_h5(pname, pvalue, iret)

! import external procedures
!
    use hdf5           , only : hid_t
    use hdf5           , only : H5F_ACC_RDONLY_F
    use hdf5           , only : H5T_NATIVE_INTEGER
    use hdf5           , only : hid_t, hsize_t, size_t
    use hdf5           , only : h5fopen_f, h5fclose_f
    use hdf5           , only : h5gopen_f, h5gclose_f
    use hdf5           , only : h5aexists_by_name_f
    use hdf5           , only : h5aopen_by_name_f, h5aread_f, h5aclose_f
    use iso_fortran_env, only : error_unit
    use mpitools       , only : nproc
    use parameters     , only : get_parameter

! local variables are not implicit by default
!
    implicit none

! subroutine arguments
!
    character(len=*), intent(in)    :: pname
    integer         , intent(inout) :: pvalue
    integer         , intent(inout) :: iret

! local variables
!
    logical            :: info
    character(len=255) :: rname
    integer            :: np
    integer(hid_t)     :: fid, gid, aid
    integer(size_t)    :: aln
    integer(hsize_t)   :: am(1) = 1

! local parameters
!
    character(len=*), parameter :: loc =                                       &
                                     'IO::read_snapshot_parameter_integer_h5()'
!
!-------------------------------------------------------------------------------
!
! reset the success flag
!
    iret = 0

! generate the filename of the restart snapshot
!
    info  = .false.
    np    = nproc + 1
    do while (.not. info .and. np >= 0)
      np = np - 1
      write (rname, "(a,'r',i6.6,'_',i5.5,'.h5')") trim(respath), nrest, np
      inquire(file = rname, exist = info)
    end do

! procees if file exists
!
    if (info) then
      call h5fopen_f(rname, H5F_ACC_RDONLY_F, fid, iret)
      if (iret >= 0) then
        call h5gopen_f(fid, 'attributes', gid, iret)
        if (iret >= 0) then
          call h5aexists_by_name_f(gid, '.', trim(pname), info, iret)
          if (info .and. iret >= 0) then
            call h5aopen_by_name_f(gid, '.', trim(pname), aid, iret)
            if (iret >= 0) then
              call h5aread_f(aid, H5T_NATIVE_INTEGER, pvalue, am(:), iret)
              call h5aclose_f(aid, iret)
            end if
          end if
          call h5gclose_f(gid, iret)
        end if
        call h5fclose_f(fid, iret)
      end if
    else
      write(error_unit,"('[', a, ']: ', a)") trim(loc)                         &
                       , "Snapshot " // trim(rname) // " file does not exist!"
      iret = 1
    end if

!-------------------------------------------------------------------------------
!
  end subroutine read_snapshot_parameter_integer_h5
!
!===============================================================================
!
! subroutine READ_SNAPSHOT_PARAMETER_INTEGER_VECTOR_H5:
! ----------------------------------------------------
!
!   Subroutine reads an integer vector parameter from the restart snapshot.
!
!   Arguments:
!
!     pname   - the parameter name;
!     pvalue  - the parameter value;
!     iret    - the success flag (the success is 0, failure otherwise);
!
!===============================================================================
!
  subroutine read_snapshot_parameter_integer_vector_h5(pname, pvalue, iret)

! import external procedures
!
    use hdf5           , only : hid_t
    use hdf5           , only : H5F_ACC_RDONLY_F
    use hdf5           , only : H5T_NATIVE_INTEGER
    use hdf5           , only : hid_t, hsize_t, size_t
    use hdf5           , only : h5fopen_f, h5fclose_f
    use hdf5           , only : h5gopen_f, h5gclose_f
    use hdf5           , only : h5aexists_by_name_f
    use hdf5           , only : h5aopen_by_name_f, h5aread_f, h5aclose_f
    use iso_fortran_env, only : error_unit
    use mpitools       , only : nproc
    use parameters     , only : get_parameter

! local variables are not implicit by default
!
    implicit none

! subroutine arguments
!
    character(len=*)     , intent(in)    :: pname
    integer, dimension(:), intent(inout) :: pvalue
    integer              , intent(inout) :: iret

! local variables
!
    logical            :: info
    character(len=255) :: rname
    integer            :: np
    integer(hid_t)     :: fid, gid, aid
    integer(size_t)    :: aln
    integer(hsize_t)   :: am(1) = 1

! local parameters
!
    character(len=*), parameter :: loc =                                       &
                             'IO::read_snapshot_parameter_integer_vector_h5()'
!
!-------------------------------------------------------------------------------
!
! reset the success flag
!
    iret = 0

! generate the filename of the restart snapshot
!
    info  = .false.
    np    = nproc + 1
    do while (.not. info .and. np >= 0)
      np = np - 1
      write (rname, "(a,'r',i6.6,'_',i5.5,'.h5')") trim(respath), nrest, np
      inquire(file = rname, exist = info)
    end do

! procees if file exists
!
    if (info) then
      call h5fopen_f(rname, H5F_ACC_RDONLY_F, fid, iret)
      if (iret >= 0) then
        call h5gopen_f(fid, 'attributes', gid, iret)
        if (iret >= 0) then
          call h5aexists_by_name_f(gid, '.', trim(pname), info, iret)
          if (info .and. iret >= 0) then
            call h5aopen_by_name_f(gid, '.', trim(pname), aid, iret)
            if (iret >= 0) then
              am(1) = size(pvalue)
              call h5aread_f(aid, H5T_NATIVE_INTEGER, pvalue, am(:), iret)
              call h5aclose_f(aid, iret)
            end if
          end if
          call h5gclose_f(gid, iret)
        end if
        call h5fclose_f(fid, iret)
      end if
    else
      write(error_unit,"('[', a, ']: ', a)") trim(loc)                         &
                       , "Snapshot " // trim(rname) // " file does not exist!"
      iret = 1
    end if

!-------------------------------------------------------------------------------
!
  end subroutine read_snapshot_parameter_integer_vector_h5
!
!===============================================================================
!
! subroutine READ_SNAPSHOT_PARAMETER_DOUBLE_H5:
! --------------------------------------------
!
!   Subroutine reads a double precision real parameter from the restart
!   snapshot.
!
!   Arguments:
!
!     pname   - the parameter name;
!     pvalue  - the parameter value;
!     iret    - the success flag (the success is 0, failure otherwise);
!
!===============================================================================
!
  subroutine read_snapshot_parameter_double_h5(pname, pvalue, iret)

! import external procedures
!
    use hdf5           , only : hid_t
    use hdf5           , only : H5F_ACC_RDONLY_F
    use hdf5           , only : H5T_NATIVE_DOUBLE
    use hdf5           , only : hid_t, hsize_t, size_t
    use hdf5           , only : h5fopen_f, h5fclose_f
    use hdf5           , only : h5gopen_f, h5gclose_f
    use hdf5           , only : h5aexists_by_name_f
    use hdf5           , only : h5aopen_by_name_f, h5aread_f, h5aclose_f
    use iso_fortran_env, only : error_unit
    use mpitools       , only : nproc
    use parameters     , only : get_parameter

! local variables are not implicit by default
!
    implicit none

! subroutine arguments
!
    character(len=*), intent(in)    :: pname
    real(kind=8)    , intent(inout) :: pvalue
    integer         , intent(inout) :: iret

! local variables
!
    logical            :: info
    character(len=255) :: rname
    integer            :: np
    integer(hid_t)     :: fid, gid, aid
    integer(size_t)    :: aln
    integer(hsize_t)   :: am(1) = 1

! local parameters
!
    character(len=*), parameter :: loc =                                       &
                                     'IO::read_snapshot_parameter_double_h5()'
!
!-------------------------------------------------------------------------------
!
! reset the success flag
!
    iret = 0

! generate the filename of the restart snapshot
!
    info  = .false.
    np    = nproc + 1
    do while (.not. info .and. np >= 0)
      np = np - 1
      write (rname, "(a,'r',i6.6,'_',i5.5,'.h5')") trim(respath), nrest, np
      inquire(file = rname, exist = info)
    end do

! procees if file exists
!
    if (info) then
      call h5fopen_f(rname, H5F_ACC_RDONLY_F, fid, iret)
      if (iret >= 0) then
        call h5gopen_f(fid, 'attributes', gid, iret)
        if (iret >= 0) then
          call h5aexists_by_name_f(gid, '.', trim(pname), info, iret)
          if (info .and. iret >= 0) then
            call h5aopen_by_name_f(gid, '.', trim(pname), aid, iret)
            if (iret >= 0) then
              call h5aread_f(aid, H5T_NATIVE_DOUBLE, pvalue, am(:), iret)
              call h5aclose_f(aid, iret)
            end if
          end if
          call h5gclose_f(gid, iret)
        end if
        call h5fclose_f(fid, iret)
      end if
    else
      write(error_unit,"('[', a, ']: ', a)") trim(loc)                         &
                       , "Snapshot " // trim(rname) // " file does not exist!"
      iret = 1
    end if

!-------------------------------------------------------------------------------
!
  end subroutine read_snapshot_parameter_double_h5
!
!===============================================================================
!
! subroutine READ_RESTART_SNAPSHOT_H5:
! -----------------------------------
!
!   Subroutine reads restart snapshot, i.e. parameters, meta and data blocks
!   stored in the HDF5 format restart files and reconstructs the data structure
!   in order to resume a terminated job.
!
!   Arguments:
!
!     iret - the return flag to inform if subroutine succeeded or failed;
!
!===============================================================================
!
  subroutine read_restart_snapshot_h5(iret)

! import external procedures and variables
!
    use blocks         , only : change_blocks_process
    use hdf5           , only : hid_t
    use hdf5           , only : H5F_ACC_RDONLY_F
    use hdf5           , only : h5fis_hdf5_f, h5fopen_f, h5fclose_f
    use iso_fortran_env, only : error_unit
#ifdef MPI
    use mesh           , only : redistribute_blocks
#endif /* MPI */
    use mpitools       , only : nprocs, npmax, nproc

! local variables are not implicit by default
!
    implicit none

! input and output arguments
!
    integer, intent(out) :: iret

! local variables
!
    character(len=255) :: fl, msg
    integer(hid_t)     :: fid
    integer            :: err, lfile
    logical            :: info

! local parameters
!
    character(len=*), parameter :: loc = 'IO::read_restart_snapshot_h5()'
!
!-------------------------------------------------------------------------------
!
! initialize success flag
!
    iret = 0

!! 1. RESTORE PARAMETERS AND META BLOCKS FROM THE FIRST FILE
!!
! prepare the filename using the current process number; in case the file does
! not exist decrease it until the file corresponding to lower process number
! is found;
!
    info  = .false.
    lfile = nproc + 1
    do while (.not. info .and. lfile > 0)
      lfile = lfile - 1
      write (fl, "(a,'r',i6.6,'_',i5.5,'.h5')") trim(respath), nrest, lfile
      inquire(file = fl, exist = info)
    end do

! quit, if file does not exist
!
    if (.not. info) then

      iret = 121
      msg  = "File " // trim(fl) // " does not exist!"

    else ! file does exist

! check if this file is in the HDF5 format
!
      call h5fis_hdf5_f(fl, info, err)

! quit, if the format verification failed or file is not in HDF5 format
!
      if (err < 0 .or. .not. info) then

        iret = 122
        if (err < 0)    msg = "Cannot check the file format!"
        if (.not. info) msg = "File " // trim(fl) // " is not an HDF5 file!"

      else ! file is HDF5

! open the HDF5 file
!
        call h5fopen_f(fl, H5F_ACC_RDONLY_F, fid, err)

! quit, if file could not be opened
!
        if (err < 0) then

          iret = 123
          msg = "Cannot open file: " // trim(fl)

        else ! file is opened

! read global attributes
!
          call read_attributes_h5(fid)

! read meta blocks and recreate the meta block hierarchy
!
          call read_metablocks_h5(fid)

! close the file
!
          call h5fclose_f(fid, err)

! quit, if file could not be closed
!
          if (err > 0) then
            iret = 124
            msg  = "Cannot close file: " // trim(fl)
          end if
        end if ! file is opened
      end if ! file is HDF5
    end if ! file does exist

!! 1. RESTORE DATA BLOCKS
!!
! separate data blocks reading into two cases, when the number of processors is
! larger or equal to the number of files, and when we have less processors than
! files
!
! first, read data blocks to processes which have corresponding restart file
!
    if (nproc < nfiles) then

! prepare the filename
!
      write (fl, "(a,'r',i6.6,'_',i5.5,'.h5')") trim(respath), nrest, nproc

! check if the HDF5 file exists
!
      inquire(file = fl, exist = info)

! quit, if file does not exist
!
      if (.not. info) then

        iret = 121
        msg  = "File " // trim(fl) // " does not exist!"

      else ! file does exist

! check if this file is in the HDF5 format
!
        call h5fis_hdf5_f(fl, info, err)

! quit, if the format verification failed or file is not in HDF5 format
!
        if (err < 0 .or. .not. info) then

          iret = 122
          if (err < 0)    msg = "Cannot check the file format!"
          if (.not. info) msg = "File " // trim(fl) // " is not an HDF5 file!"

        else ! file is HDF5

! open the HDF5 file
!
          call h5fopen_f(fl, H5F_ACC_RDONLY_F, fid, err)

! quit, if file could not be opened
!
          if (err < 0) then

            iret = 123
            msg = "Cannot open file: " // trim(fl)

          else ! file is opened

! read data blocks
!
            call read_datablocks_h5(fid)

! close the file
!
            call h5fclose_f(fid, err)

! quit, if file could not be closed
!
            if (err > 0) then
              iret = 124
              msg  = "Cannot close file: " // trim(fl)
            end if
          end if ! file is opened
        end if ! file is HDF5
      end if ! file exists
    end if ! nproc < nfiles

! if there are more files than processes, read the remaining files by
! the last process and redistribute blocks after each processed file,
! otherwise only redistribute blocks
!
    if (nprocs < nfiles) then

! iterate over remaining files and read one by one to the last block
!
      do lfile = nprocs, nfiles - 1

! switch meta blocks from the read file to belong to the reading process
!
        call change_blocks_process(lfile, npmax)

! read the remaining files by the last process only
!
        if (nproc == npmax) then

! prepare the filename
!
          write (fl, "(a,'r',i6.6,'_',i5.5,'.h5')") trim(respath), nrest, lfile

! check if the HDF5 file exists
!
          inquire(file = fl, exist = info)

! quit, if file does not exist
!
          if (.not. info) then

            iret = 121
            msg  = "File " // trim(fl) // " does not exist!"

          else ! file does exist

! check if this file is in the HDF5 format
!
            call h5fis_hdf5_f(fl, info, err)

! quit, if the format verification failed or file is not in HDF5 format
!
            if (err < 0 .or. .not. info) then

              iret = 122
              if (err < 0)    msg = "Cannot check the file format!"
              if (.not. info) msg = "File " // trim(fl) //                     &
                                                       " is not an HDF5 file!"

            else ! file is HDF5

! open the HDF5 file
!
              call h5fopen_f(fl, H5F_ACC_RDONLY_F, fid, err)

! quit, if file could not be opened
!
              if (err < 0) then

                iret = 123
                msg = "Cannot open file: " // trim(fl)

              else ! file is opened

! read data blocks
!
                call read_datablocks_h5(fid)

! close the file
!
                call h5fclose_f(fid, err)

! quit, if file could not be closed
!
                if (err > 0) then
                  iret = 124
                  msg  = "Cannot close file: " // trim(fl)
                end if
              end if ! file is opened
            end if ! file is HDF5
          end if ! file exists
        end if ! nproc == npmax

#ifdef MPI
! redistribute blocks between processors
!
        call redistribute_blocks()
#endif /* MPI */

      end do ! lfile = nprocs, nfiles - 1

    else ! nprocs < nfiles

#ifdef MPI
! redistribute blocks between processors
!
      call redistribute_blocks()
#endif /* MPI */

    end if ! nprocs < nfiles

! deallocate the array of block pointers
!
    if (allocated(block_array)) deallocate(block_array)

! if there was any problem, print the message
!
    if (iret > 0) write(error_unit,"('[',a,']: ',a)") trim(loc), trim(msg)

!-------------------------------------------------------------------------------
!
  end subroutine read_restart_snapshot_h5
!
!===============================================================================
!
! subroutine WRITE_RESTART_SNAPSHOT_H5:
! ------------------------------------
!
!   Subroutine writes restart snapshot, i.e. parameters, meta and data blocks
!   to the HDF5 format restart files in order to resume a terminated job later.
!
!   Arguments:
!
!     nrun - the snapshot number;
!     iret - the return flag to inform if subroutine succeeded or failed;
!
!===============================================================================
!
  subroutine write_restart_snapshot_h5(nrun, iret)

! import external procedures and variables
!
    use hdf5           , only : hid_t
    use hdf5           , only : H5F_ACC_TRUNC_F, H5F_SCOPE_GLOBAL_F
    use hdf5           , only : h5fcreate_f, h5fflush_f, h5fclose_f
    use iso_fortran_env, only : error_unit
    use mpitools       , only : nproc

! local variables are not implicit by default
!
    implicit none

! input and output arguments
!
    integer, intent(in)  :: nrun
    integer, intent(out) :: iret

! local variables
!
    character(len=64) :: fl
    integer(hid_t)    :: fid
    integer           :: err

! local parameters
!
    character(len=*), parameter :: loc = 'IO::write_restart_snapshot_h5()'
!
!-------------------------------------------------------------------------------
!
! prepare the restart snapshot filename
!
    write (fl, "('r',i6.6,'_',i5.5,'.h5')") nrun, nproc

! create the new HDF5 file to store the snapshot
!
    call h5fcreate_f(fl, H5F_ACC_TRUNC_F, fid, err)

! if the file could not be created, print message and quit
!
    if (err < 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create file: " // trim(fl)
      iret = 201
      return
    end if

! write the global attributes
!
    call write_attributes_h5(fid)

! write all metablocks which represent the internal structure of domain
!
    call write_metablocks_h5(fid)

! write all datablocks which represent the all variables
!
    call write_datablocks_h5(fid)

! flush the file
!
    call h5fflush_f(fid, H5F_SCOPE_GLOBAL_F, err)

! close the file
!
    call h5fclose_f(fid, err)

! if the file could not be closed print message and quit
!
    if (err > 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot close file: " // trim(fl)
      iret = 203
      return
    end if

!-------------------------------------------------------------------------------
!
  end subroutine write_restart_snapshot_h5
!
!===============================================================================
!
! subroutine WRITE_SNAPSHOT_H5:
! ----------------------------
!
!   Subroutine writes the current simulation snapshot, i.e. parameters,
!   coordinates and variables to the HDF5 format files for further processing.
!
!
!===============================================================================
!
  subroutine write_snapshot_h5()

! import external procedures and variables
!
    use hdf5           , only : hid_t
    use hdf5           , only : H5F_ACC_TRUNC_F, H5F_SCOPE_GLOBAL_F
    use hdf5           , only : h5fcreate_f, h5fflush_f, h5fclose_f
    use iso_fortran_env, only : error_unit
    use mpitools       , only : nproc

! local variables are not implicit by default
!
    implicit none

! local variables
!
    character(len=64) :: fl
    integer(hid_t)    :: fid
    integer           :: err

! local parameters
!
    character(len=*), parameter :: loc = 'IO::write_snapshot_h5()'
!
!-------------------------------------------------------------------------------
!
! prepare the restart snapshot filename
!
    write (fl, "(a1,i6.6,'_',i5.5,'.h5')") ftype, isnap, nproc

! create the new HDF5 file to store the snapshot
!
    call h5fcreate_f(fl, H5F_ACC_TRUNC_F, fid, err)

! if the file could not be created, print message and quit
!
    if (err < 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create file: " // trim(fl)
      return
    end if

! write the global attributes
!
    call write_attributes_h5(fid)

! write the coordinates (data block bounds, refinement levels, etc.)
!
    call write_coordinates_h5(fid)

! depending on the selected type of output file write the right groups
!
    select case(ftype)

    case('c')

! write the variables stored in data blocks (leafs)
!
      call write_conservative_variables_h5(fid)

    case('p')

! write the variables stored in data blocks (leafs)
!
      call write_primitive_variables_h5(fid)

    case default

! print information about unsupported file format and quit
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "File type is not suppoerted!"
      call h5fclose_f(fid, err)
      return

    end select

! flush the file
!
    call h5fflush_f(fid, H5F_SCOPE_GLOBAL_F, err)

! close the file
!
    call h5fclose_f(fid, err)

! if the file could not be closed print message and quit
!
    if (err > 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot close file: " // trim(fl)
      return
    end if

!-------------------------------------------------------------------------------
!
  end subroutine write_snapshot_h5
!
!===============================================================================
!
! subroutine WRITE_ATTRIBUTES_H5:
! ------------------------------
!
!   Subroutine stores global attributes in the HDF5 file provided by an
!   identifier.
!
!   Arguments:
!
!     fid - the HDF5 file identifier;
!
!===============================================================================
!
  subroutine write_attributes_h5(fid)

! import external procedures and variables
!
    use blocks         , only : get_mblocks, get_dblocks, get_nleafs
    use blocks         , only : get_last_id
    use coordinates    , only : minlev, maxlev
    use coordinates    , only : ncells, nghosts
    use coordinates    , only : ddims => domain_base_dims
    use coordinates    , only : xmin, xmax, ymin, ymax, zmin, zmax
    use coordinates    , only : periodic
    use equations      , only : eqsys, eos, gamma, csnd
    use evolution      , only : step, time, dt, dtn
    use hdf5           , only : hid_t
    use hdf5           , only : h5gcreate_f, h5gclose_f
    use iso_fortran_env, only : error_unit
    use mpitools       , only : nprocs, nproc
    use problems       , only : problem_name
    use random         , only : nseeds, get_seeds

! local variables are not implicit by default
!
    implicit none

! subroutine arguments
!
    integer(hid_t), intent(in) :: fid

! local variables
!
    integer(hid_t)                :: gid
    integer                       :: err

! local vectors
!
    integer, dimension(3)         :: bdims = 1
    integer, dimension(3)         :: per

! local allocatable arrays
!
    integer(kind=4), dimension(:), allocatable :: seeds

! local parameters
!
    character(len=*), parameter :: loc = 'IO::write_attributes_h5()'
!
!-------------------------------------------------------------------------------
!
! store the code name in order to determine the format of data
!
    call write_attribute(fid, 'code'   , 'AMUN')
    call write_attribute(fid, 'version', 'v1.0')

! create a group to store the global attributes
!
    call h5gcreate_f(fid, 'attributes', gid, err)

! check if the group has been created successfuly
!
    if (err < 0) then

! print error about the problem with creating the group
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create the group!"

! return from the subroutine
!
      return

    end if

! convert periodic(:) to an integer vector
!
    per(:) = merge(1, 0, periodic(:))

! store string attributes
!
    call write_attribute(gid, 'problem', problem_name )
    call write_attribute(gid, 'eqsys'  , eqsys        )
    call write_attribute(gid, 'eos'    , eos          )

! store the integer attributes
!
    call write_attribute(gid, 'ndims'   , NDIMS        )
    call write_attribute(gid, 'last_id' , get_last_id())
    call write_attribute(gid, 'mblocks' , get_mblocks())
    call write_attribute(gid, 'dblocks' , get_dblocks())
    call write_attribute(gid, 'nleafs'  , get_nleafs() )
    call write_attribute(gid, 'ncells'  , ncells       )
    call write_attribute(gid, 'nghosts' , nghosts      )
    call write_attribute(gid, 'minlev'  , minlev       )
    call write_attribute(gid, 'maxlev'  , maxlev       )
    call write_attribute(gid, 'nprocs'  , nprocs       )
    call write_attribute(gid, 'nproc'   , nproc        )
    call write_attribute(gid, 'nseeds'  , nseeds       )
    call write_attribute(gid, 'step'    , step         )
    call write_attribute(gid, 'isnap'   , isnap        )
    call write_attribute(gid, 'periodic', per(:)       )

! store the real attributes
!
    call write_attribute(gid, 'xmin', xmin)
    call write_attribute(gid, 'xmax', xmax)
    call write_attribute(gid, 'ymin', ymin)
    call write_attribute(gid, 'ymax', ymax)
    call write_attribute(gid, 'zmin', zmin)
    call write_attribute(gid, 'zmax', zmax)
    call write_attribute(gid, 'time', time)
    call write_attribute(gid, 'dt'  , dt  )
    call write_attribute(gid, 'dtn' , dtn )
    if (eos == 'adi') then
      call write_attribute(gid, 'gamma', gamma)
    end if
    if (eos == 'iso') then
      call write_attribute(gid, 'csnd' , csnd )
    end if

! store the vector attributes
!
    bdims(1:NDIMS) = ncells
    call write_attribute(gid, 'dims'            , bdims)
    call write_attribute(gid, 'domain_base_dims', ddims)

! store random number generator seed values
!
    if (nseeds > 0) then

! allocate space for seeds
!
      allocate(seeds(nseeds))

! get the seed values
!
      call get_seeds(seeds)

! store them in the current group
!
      call write_attribute(gid, 'seeds', seeds(:))

! deallocate seed array
!
      deallocate(seeds)

    end if ! nseeds > 0

! close the group
!
    call h5gclose_f(gid, err)

! check if the group has been closed successfuly
!
    if (err < 0) then

! print error about the problem with closing the group
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot close the group!"

    end if

!-------------------------------------------------------------------------------
!
  end subroutine write_attributes_h5
!
!===============================================================================
!
! subroutine READ_ATTRIBUTES_H5:
! -----------------------------
!
!   Subroutine restores global attributes from an HDF5 file provided by its
!   identifier.
!
!   Arguments:
!
!     fid - the HDF5 file identifier;
!
!===============================================================================
!
  subroutine read_attributes_h5(fid)

! import external procedures and variables
!
    use blocks         , only : block_meta
    use blocks         , only : append_metablock
    use blocks         , only : set_last_id, get_last_id
    use blocks         , only : get_mblocks, get_dblocks, get_nleafs
    use coordinates    , only : ncells, nghosts
    use coordinates    , only : xmin, xmax, ymin, ymax, zmin, zmax
    use evolution      , only : step, time, dt, dtn
    use hdf5           , only : hid_t, hsize_t
    use hdf5           , only : h5gopen_f, h5gclose_f
    use iso_fortran_env, only : error_unit
    use mpitools       , only : nprocs, nproc
    use random         , only : nseeds, set_seeds

! local variables are not implicit by default
!
    implicit none

! subroutine arguments
!
    integer(hid_t), intent(in) :: fid

! local variables
!
    integer(hid_t) :: gid
    integer        :: ierr, l
    integer        :: lndims, lmblocks, lnleafs, llast_id
    integer        :: lncells, lnghost, lnproc, lnseeds

! local pointers
!
    type(block_meta), pointer :: pmeta

! allocatable arrays
!
    integer(kind=4), dimension(:), allocatable :: seeds

! local parameters
!
    character(len=*), parameter :: loc = 'IO::read_attributes_h5()'
!
!-------------------------------------------------------------------------------
!
! open the global attributes group
!
    call h5gopen_f(fid, 'attributes', gid, ierr)

! check if the group has been opened successfuly
!
    if (ierr < 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot open the group!"
      return
    end if

! restore integer attributes
!
    call read_attribute(gid, 'ndims'  , lndims  )
    call read_attribute(gid, 'nprocs' , nfiles  )
    call read_attribute(gid, 'nproc'  , lnproc  )
    call read_attribute(gid, 'mblocks', lmblocks)
    call read_attribute(gid, 'nleafs' , lnleafs )
    call read_attribute(gid, 'last_id', llast_id)
    call read_attribute(gid, 'ncells' , lncells )
    call read_attribute(gid, 'nghosts', lnghost )
    call read_attribute(gid, 'nseeds' , lnseeds )
    call read_attribute(gid, 'step'   , step    )
    call read_attribute(gid, 'isnap'  , isnap   )

! restore double precision attributes
!
    call read_attribute(gid, 'xmin', xmin)
    call read_attribute(gid, 'xmax', xmax)
    call read_attribute(gid, 'ymin', ymin)
    call read_attribute(gid, 'ymax', ymax)
    call read_attribute(gid, 'zmin', zmin)
    call read_attribute(gid, 'zmax', zmax)
    call read_attribute(gid, 'time', time)
    call read_attribute(gid, 'dt'  , dt  )
    call read_attribute(gid, 'dtn' , dtn )

! check the number of dimensions
!
    if (lndims /= NDIMS) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "The number of dimensions does not match!"
      return
    end if

! check the block dimensions
!
    if (lncells /= ncells) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "The block dimensions do not match!"
    end if

! check the number of ghost layers
!
    if (lnghost /= nghosts) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "The number of ghost layers does not match!"
    end if

! allocate space for seeds
!
    allocate(seeds(lnseeds))

! store them in the current group
!
    call read_attribute(gid, 'seeds', seeds(:))

! set the seed values
!
    call set_seeds(lnseeds, seeds(:), nproc /= lnproc)

! deallocate seed array
!
    deallocate(seeds)

! allocate all metablocks
!
    do l = 1, lmblocks
      call append_metablock(pmeta)
    end do

! check if the number of created metablocks is equal to lbmcloks
!
    if (lmblocks /= get_mblocks()) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Number of metablocks does not match!"
    end if

! allocate an array of pointers with the size llast_id
!
    allocate(block_array(llast_id))

! set the last_id
!
    call set_last_id(llast_id)

! close the group
!
    call h5gclose_f(gid, ierr)

! check if the group has been closed successfuly
!
    if (ierr /= 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot close the group!"
    end if

!-------------------------------------------------------------------------------
!
  end subroutine read_attributes_h5
!
!===============================================================================
!
! subroutine WRITE_METABLOCKS_H5:
! ------------------------------
!
!   Subroutine stores all meta blocks with their complete fields in 'metablock'
!   group in a provided file identifier.
!
!   Arguments:
!
!     fid - the HDF5 file identifier;
!
!===============================================================================
!
  subroutine write_metablocks_h5(fid)

! import procedures and variables from other modules
!
    use blocks         , only : block_meta, list_meta
    use blocks         , only : ndims, nchildren, nsides
    use blocks         , only : get_last_id, get_mblocks
    use hdf5           , only : hid_t, hsize_t
    use hdf5           , only : h5gcreate_f, h5gclose_f
    use iso_fortran_env, only : error_unit

! local variables are not implicit by default
!
    implicit none

! subroutine arguments
!
    integer(hid_t), intent(in) :: fid

! local variables
!
    integer(hid_t)                 :: gid
    integer(kind=4)                :: i, j, k, l, n, p
    integer                        :: iret
    integer(hsize_t), dimension(1) :: am, cm
    integer(hsize_t), dimension(2) :: dm, pm
#if NDIMS == 2
    integer(hsize_t), dimension(4) :: nm
#endif /* NDIMS == 2 */
#if NDIMS == 3
    integer(hsize_t), dimension(5) :: nm
#endif /* NDIMS == 3 */

! local allocatable arrays
!
    integer(kind=4), dimension(:)  , allocatable :: idx
    integer(kind=4), dimension(:)  , allocatable :: par, dat
    integer(kind=4), dimension(:)  , allocatable ::  id, cpu, lev, cfg, ref, lea
    real   (kind=8), dimension(:)  , allocatable :: xmn, xmx, ymn, ymx, zmn, zmx
    integer(kind=4), dimension(:,:), allocatable :: chl, pos, cor
#if NDIMS == 2
    integer(kind=4), dimension(:,:,:,:)  , allocatable :: edges
    integer(kind=4), dimension(:,:,:)    , allocatable :: corners
#endif /* NDIMS == 2 */
#if NDIMS == 3
    integer(kind=4), dimension(:,:,:,:,:), allocatable :: faces
    integer(kind=4), dimension(:,:,:,:,:), allocatable :: edges
    integer(kind=4), dimension(:,:,:,:)  , allocatable :: corners
#endif /* NDIMS == 3 */

! local pointers
!
    type(block_meta), pointer :: pmeta

! local parameters
!
    character(len=*), parameter :: loc = 'IO::write_metablocks_h5()'
!
!-------------------------------------------------------------------------------
!
! create the group for metadata
!
    call h5gcreate_f(fid, 'metablocks', gid, iret)

! check if the group has been created successfuly
!
    if (iret >= 0) then

! prepate dimensions
!
      am(1) = get_mblocks()
      cm(1) = get_last_id()
      dm(1) = get_mblocks()
      dm(2) = nchildren
      pm(1) = get_mblocks()
      pm(2) = NDIMS
      nm(1) = get_mblocks()
      nm(2) = nsides
      nm(3) = nsides
#if NDIMS == 2
      nm(4) = ndims
#endif /* NDIMS == 2 */
#if NDIMS == 3
      nm(4) = nsides
      nm(5) = ndims
#endif /* NDIMS == 3 */

! only store data from processes that have any meta blocks
!
      if (am(1) > 0) then

! allocate arrays to store meta block fields
!
        allocate(idx(cm(1)))
        allocate(par(am(1)))
        allocate(dat(am(1)))
        allocate(id (am(1)))
        allocate(cpu(am(1)))
        allocate(lev(am(1)))
        allocate(cfg(am(1)))
        allocate(ref(am(1)))
        allocate(lea(am(1)))
        allocate(xmn(am(1)))
        allocate(xmx(am(1)))
        allocate(ymn(am(1)))
        allocate(ymx(am(1)))
        allocate(zmn(am(1)))
        allocate(zmx(am(1)))
        allocate(chl(dm(1),dm(2)))
        allocate(pos(pm(1),pm(2)))
        allocate(cor(pm(1),pm(2)))
#if NDIMS == 2
        allocate(edges  (nm(1),nm(2),nm(3),nm(4)))
        allocate(corners(nm(1),nm(2),nm(3)))
#endif /* NDIMS == 2 */
#if NDIMS == 3
        allocate(faces  (nm(1),nm(2),nm(3),nm(4),nm(5)))
        allocate(edges  (nm(1),nm(2),nm(3),nm(4),nm(5)))
        allocate(corners(nm(1),nm(2),nm(3),nm(4)))
#endif /* NDIMS == 3 */

! reset stored arrays
!
        idx(:)           = -1
        par(:)           = -1
        dat(:)           = -1
        lea(:)           = -1
        chl(:,:)         = -1
#if NDIMS == 2
        edges(:,:,:,:)   = -1
        corners(:,:,:)   = -1
#endif /* NDIMS == 2 */
#if NDIMS == 3
        faces(:,:,:,:,:) = -1
        edges(:,:,:,:,:) = -1
        corners(:,:,:,:) = -1
#endif /* NDIMS == 3 */

! reset the block counter
!
        l = 0

! associate pmeta with the first block on the meta block list
!
        pmeta => list_meta

! iterate over all meta blocks and fill in the arrays for storage
!
        do while(associated(pmeta))

! increase the block counter
!
          l = l + 1

! store meta block fields
!
          idx(pmeta%id) = l

          if (associated(pmeta%parent)) par(l) = pmeta%parent%id
          if (associated(pmeta%data)  ) dat(l) = 1

          id (l)   = pmeta%id
          cpu(l)   = pmeta%process
          lev(l)   = pmeta%level
          cfg(l)   = pmeta%conf
          ref(l)   = pmeta%refine
          pos(l,:) = pmeta%pos(:)
          cor(l,:) = pmeta%coords(:)

          if (pmeta%leaf) lea(l) = 1

          xmn(l)   = pmeta%xmin
          xmx(l)   = pmeta%xmax
          ymn(l)   = pmeta%ymin
          ymx(l)   = pmeta%ymax
          zmn(l)   = pmeta%zmin
          zmx(l)   = pmeta%zmax

          do p = 1, nchildren
            if (associated(pmeta%child(p)%ptr)) chl(l,p) = pmeta%child(p)%ptr%id
          end do

! store face, edge and corner neighbor pointers
!
#if NDIMS == 2
          do i = 1, nsides
            do j = 1, nsides
              do n = 1, ndims
                if (associated(pmeta%edges(i,j,n)%ptr))                        &
                                    edges(l,i,j,n) = pmeta%edges(i,j,n)%ptr%id
              end do ! ndims
              if (associated(pmeta%corners(i,j)%ptr))                          &
                                    corners(l,i,j) = pmeta%corners(i,j)%ptr%id
            end do ! i = 1, nsides
          end do ! j = 1, nsides
#endif /* NDIMS == 2 */
#if NDIMS == 3
          do i = 1, nsides
            do j = 1, nsides
              do k = 1, nsides
                do n = 1, ndims
                  if (associated(pmeta%faces(i,j,k,n)%ptr))                    &
                                faces(l,i,j,k,n) = pmeta%faces(i,j,k,n)%ptr%id
                  if (associated(pmeta%edges(i,j,k,n)%ptr))                    &
                                edges(l,i,j,k,n) = pmeta%edges(i,j,k,n)%ptr%id
                end do ! ndims
                if (associated(pmeta%corners(i,j,k)%ptr))                      &
                                corners(l,i,j,k) = pmeta%corners(i,j,k)%ptr%id
              end do ! i = 1, nsides
            end do ! j = 1, nsides
          end do ! k = 1, nsides
#endif /* NDIMS == 3 */

! associate pmeta with the next block on the list
!
          pmeta => pmeta%next

        end do ! over all meta blocks

! store meta block data in the HDF5 file
!
        call write_array(gid, 'indices', cm(1)  , idx)
        call write_array(gid, 'parent' , am(1)  , par)
        call write_array(gid, 'data'   , am(1)  , dat)
        call write_array(gid, 'id'     , am(1)  , id )
        call write_array(gid, 'cpu'    , am(1)  , cpu)
        call write_array(gid, 'level'  , am(1)  , lev)
        call write_array(gid, 'config' , am(1)  , cfg)
        call write_array(gid, 'refine' , am(1)  , ref)
        call write_array(gid, 'leaf'   , am(1)  , lea)
        call write_array(gid, 'xmin'   , am(1)  , xmn)
        call write_array(gid, 'xmax'   , am(1)  , xmx)
        call write_array(gid, 'ymin'   , am(1)  , ymn)
        call write_array(gid, 'ymax'   , am(1)  , ymx)
        call write_array(gid, 'zmin'   , am(1)  , zmn)
        call write_array(gid, 'zmax'   , am(1)  , zmx)
        call write_array(gid, 'child'  , dm(:)  , chl(:,:))
        call write_array(gid, 'pos'    , pm(:)  , pos(:,:))
        call write_array(gid, 'coord'  , pm(:)  , cor(:,:))
#if NDIMS == 2
        call write_array(gid, 'edges'  , nm(1:4), edges(:,:,:,:))
        call write_array(gid, 'corners', nm(1:3), corners(:,:,:))
#endif /* NDIMS == 2 */
#if NDIMS == 3
        call write_array(gid, 'faces'  , nm(1:5), faces(:,:,:,:,:))
        call write_array(gid, 'edges'  , nm(1:5), edges(:,:,:,:,:))
        call write_array(gid, 'corners', nm(1:4), corners(:,:,:,:))
#endif /* NDIMS == 3 */

! deallocate allocated arrays
!
        if (allocated(idx))     deallocate(idx)
        if (allocated(par))     deallocate(par)
        if (allocated(dat))     deallocate(dat)
        if (allocated(id) )     deallocate(id)
        if (allocated(cpu))     deallocate(cpu)
        if (allocated(lev))     deallocate(lev)
        if (allocated(cfg))     deallocate(cfg)
        if (allocated(ref))     deallocate(ref)
        if (allocated(lea))     deallocate(lea)
        if (allocated(xmn))     deallocate(xmn)
        if (allocated(xmx))     deallocate(xmx)
        if (allocated(ymn))     deallocate(ymn)
        if (allocated(ymx))     deallocate(ymx)
        if (allocated(zmn))     deallocate(zmn)
        if (allocated(zmx))     deallocate(zmx)
        if (allocated(chl))     deallocate(chl)
        if (allocated(cor))     deallocate(cor)
#if NDIMS == 3
        if (allocated(faces))   deallocate(faces)
#endif /* NDIMS == 3 */
        if (allocated(edges))   deallocate(edges)
        if (allocated(corners)) deallocate(corners)

      end if ! meta blocks > 0

! close the group
!
      call h5gclose_f(gid, iret)

! check if the group has been closed successfuly
!
      if (iret > 0) then

! print error about the problem with closing the group
!
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot close the group!"

      end if

    else

! print error about the problem with creating the group
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create the group!"

    end if

!-------------------------------------------------------------------------------
!
  end subroutine write_metablocks_h5
!
!===============================================================================
!
! subroutine READ_METABLOCKS_H5:
! -----------------------------
!
!   Subroutine restores all meta blocks with their complete fields from
!   'metablock' group in a provided restart file identifier.
!
!   Arguments:
!
!     fid - the HDF5 file identifier;
!
!===============================================================================
!
  subroutine read_metablocks_h5(fid)

! import procedures and variables from other modules
!
    use blocks         , only : block_meta, list_meta
    use blocks         , only : ndims, nchildren, nsides
    use blocks         , only : get_mblocks
    use blocks         , only : metablock_set_id, metablock_set_process
    use blocks         , only : metablock_set_refinement
    use blocks         , only : metablock_set_configuration
    use blocks         , only : metablock_set_level, metablock_set_position
    use blocks         , only : metablock_set_coordinates, metablock_set_bounds
    use blocks         , only : metablock_set_leaf
    use hdf5           , only : hid_t, hsize_t
    use hdf5           , only : h5gopen_f, h5gclose_f
    use iso_fortran_env, only : error_unit
    use mpitools       , only : nprocs

! local variables are not implicit by default
!
    implicit none

! subroutine arguments
!
    integer(hid_t), intent(in) :: fid

! local variables
!
    integer(hid_t)                 :: gid
    integer(kind=4)                :: i, j, k, l, p, n, ip
    integer                        :: err
    integer(hsize_t), dimension(1) :: am
    integer(hsize_t), dimension(2) :: dm, pm
#if NDIMS == 2
    integer(hsize_t), dimension(4) :: nm
#endif /* NDIMS == 2 */
#if NDIMS == 3
    integer(hsize_t), dimension(5) :: nm
#endif /* NDIMS == 3 */

! local allocatable arrays
!
    integer(kind=4), dimension(:)  , allocatable :: idx
    integer(kind=4), dimension(:)  , allocatable :: par, dat
    integer(kind=4), dimension(:)  , allocatable ::  id, cpu, lev, cfg, ref, lea
    real   (kind=8), dimension(:)  , allocatable :: xmn, xmx, ymn, ymx, zmn, zmx
    integer(kind=4), dimension(:,:), allocatable :: chl, pos, cor
#if NDIMS == 2
    integer(kind=4), dimension(:,:,:,:)  , allocatable :: edges
    integer(kind=4), dimension(:,:,:)    , allocatable :: corners
#endif /* NDIMS == 2 */
#if NDIMS == 3
    integer(kind=4), dimension(:,:,:,:,:), allocatable :: faces
    integer(kind=4), dimension(:,:,:,:,:), allocatable :: edges
    integer(kind=4), dimension(:,:,:,:)  , allocatable :: corners
#endif /* NDIMS == 3 */

! local pointers
!
    type(block_meta), pointer :: pmeta

! local parameters
!
    character(len=*), parameter :: loc = 'IO::read_metablocks_h5()'
!
!-------------------------------------------------------------------------------
!
! open metablock group
!
    call h5gopen_f(fid, 'metablocks', gid, err)

! check if the group has been opened successfuly
!
    if (err >= 0) then

! prepate dimensions
!
      am(1) = get_mblocks()
      dm(1) = get_mblocks()
      dm(2) = nchildren
      pm(1) = get_mblocks()
      pm(2) = NDIMS
      nm(1) = get_mblocks()
      nm(2) = nsides
      nm(3) = nsides
#if NDIMS == 2
      nm(4) = ndims
#endif /* NDIMS == 2 */
#if NDIMS == 3
      nm(4) = nsides
      nm(5) = ndims
#endif /* NDIMS == 3 */

! allocate arrays to restore metablocks data
!
      allocate(id (am(1)))
      allocate(cpu(am(1)))
      allocate(lev(am(1)))
      allocate(par(am(1)))
      allocate(dat(am(1)))
      allocate(cfg(am(1)))
      allocate(ref(am(1)))
      allocate(lea(am(1)))
      allocate(xmn(am(1)))
      allocate(xmx(am(1)))
      allocate(ymn(am(1)))
      allocate(ymx(am(1)))
      allocate(zmn(am(1)))
      allocate(zmx(am(1)))
      allocate(chl(dm(1),dm(2)))
      allocate(pos(pm(1),pm(2)))
      allocate(cor(pm(1),pm(2)))
#if NDIMS == 2
      allocate(edges  (nm(1),nm(2),nm(3),nm(4)))
      allocate(corners(nm(1),nm(2),nm(3)))
#endif /* NDIMS == 2 */
#if NDIMS == 3
      allocate(faces  (nm(1),nm(2),nm(3),nm(4),nm(5)))
      allocate(edges  (nm(1),nm(2),nm(3),nm(4),nm(5)))
      allocate(corners(nm(1),nm(2),nm(3),nm(4)))
#endif /* NDIMS == 3 */

! reset vectors
!
      par(:)           = -1
      dat(:)           = -1
      lea(:)           = -1
      chl(:,:)         = -1
#if NDIMS == 2
      edges(:,:,:,:)   = -1
      corners(:,:,:)   = -1
#endif /* NDIMS == 2 */
#if NDIMS == 3
      faces(:,:,:,:,:) = -1
      edges(:,:,:,:,:) = -1
      corners(:,:,:,:) = -1
#endif /* NDIMS == 3 */

! read metablock fields from the HDF5 file
!
      call read_array(gid, 'id'     , am(:), id (:))
      call read_array(gid, 'cpu'    , am(:), cpu(:))
      call read_array(gid, 'level'  , am(:), lev(:))
      call read_array(gid, 'config' , am(:), cfg(:))
      call read_array(gid, 'refine' , am(:), ref(:))
      call read_array(gid, 'leaf'   , am(:), lea(:))
      call read_array(gid, 'parent' , am(:), par(:))
      call read_array(gid, 'xmin'   , am(:), xmn(:))
      call read_array(gid, 'xmax'   , am(:), xmx(:))
      call read_array(gid, 'ymin'   , am(:), ymn(:))
      call read_array(gid, 'ymax'   , am(:), ymx(:))
      call read_array(gid, 'zmin'   , am(:), zmn(:))
      call read_array(gid, 'zmax'   , am(:), zmx(:))
      call read_array(gid, 'pos'    , pm(:), pos(:,:))
      call read_array(gid, 'coord'  , pm(:), cor(:,:))
      call read_array(gid, 'child'  , dm(:), chl(:,:))
#if NDIMS == 2
      call read_array(gid, 'edges'  , nm(1:4), edges(:,:,:,:))
      call read_array(gid, 'corners', nm(1:3), corners(:,:,:))
#endif /* NDIMS == 2 */
#if NDIMS == 3
      call read_array(gid, 'faces'  , nm(1:5), faces(:,:,:,:,:))
      call read_array(gid, 'edges'  , nm(1:5), edges(:,:,:,:,:))
      call read_array(gid, 'corners', nm(1:4), corners(:,:,:,:))
#endif /* NDIMS == 3 */

! reset the block counter
!
      l = 0

! associate pmeta with the first block on the meta block list
!
      pmeta => list_meta

! iterate over all meta blocks and restore their fields
!
      do while(associated(pmeta))

! increase the block counter
!
        l = l + 1

! restore meta block fields
!
        block_array(id(l))%ptr => pmeta

        call metablock_set_id           (pmeta, id (l))
        call metablock_set_process      (pmeta, cpu(l))
        call metablock_set_refinement   (pmeta, ref(l))
        call metablock_set_configuration(pmeta, cfg(l))
        call metablock_set_level        (pmeta, lev(l))
        call metablock_set_position     (pmeta, pos(l,1), pos(l,2), pos(l,3))
        call metablock_set_coordinates  (pmeta, cor(l,1), cor(l,2), cor(l,3))
        call metablock_set_bounds       (pmeta, xmn(l), xmx(l), ymn(l), ymx(l) &
                                                              , zmn(l), zmx(l))

        if (lea(l) == 1) call metablock_set_leaf(pmeta)

! associate pmeta with the next block on the list
!
        pmeta => pmeta%next

      end do ! over all meta blocks

! reset the block counter
!
      l = 0

! associate pmeta with the first block on the meta block list
!
      pmeta => list_meta

! iterate over all meta blocks and restore their pointers
!
      do while(associated(pmeta))

! increase the block counter
!
        l = l + 1

! restore %parent pointer
!
        if (par(l) > 0) pmeta%parent => block_array(par(l))%ptr

! restore %child pointers
!
        do p = 1, nchildren
          if (chl(l,p) > 0) then
            pmeta%child(p)%ptr => block_array(chl(l,p))%ptr
          end if
        end do ! p = 1, nchildren

! restore %faces, %edges and %corners neighbor pointers
!
#if NDIMS == 2
        do i = 1, nsides
          do j = 1, nsides
            do n = 1, ndims
              ip = edges(l,i,j,n)
              if (ip > 0) pmeta%edges(i,j,n)%ptr => block_array(ip)%ptr
            end do ! n = 1, ndims
            ip = corners(l,i,j)
            if (ip > 0) pmeta%corners(i,j)%ptr => block_array(ip)%ptr
          end do ! i = 1, nsides
        end do ! j = 1, nsides
#endif /* NDIMS == 2 */
#if NDIMS == 3
        do i = 1, nsides
          do j = 1, nsides
            do k = 1, nsides
              do n = 1, ndims
                ip = faces(l,i,j,k,n)
                if (ip > 0) pmeta%faces(i,j,k,n)%ptr => block_array(ip)%ptr
                ip = edges(l,i,j,k,n)
                if (ip > 0) pmeta%edges(i,j,k,n)%ptr => block_array(ip)%ptr
              end do ! n = 1, ndims
              ip = corners(l,i,j,k)
              if (ip > 0) pmeta%corners(i,j,k)%ptr => block_array(ip)%ptr
            end do ! i = 1, nsides
          end do ! j = 1, nsides
        end do ! k = 1, nsides
#endif /* NDIMS == 3 */

! associate pmeta with the next block on the list
!
        pmeta => pmeta%next

      end do ! over all meta blocks

! deallocate allocatable arrays
!
      if (allocated(id) )     deallocate(id )
      if (allocated(par))     deallocate(par)
      if (allocated(dat))     deallocate(dat)
      if (allocated(cpu))     deallocate(cpu)
      if (allocated(lev))     deallocate(lev)
      if (allocated(cfg))     deallocate(cfg)
      if (allocated(ref))     deallocate(ref)
      if (allocated(lea))     deallocate(lea)
      if (allocated(xmn))     deallocate(xmn)
      if (allocated(xmx))     deallocate(xmx)
      if (allocated(ymn))     deallocate(ymn)
      if (allocated(ymx))     deallocate(ymx)
      if (allocated(zmn))     deallocate(zmn)
      if (allocated(zmx))     deallocate(zmx)
      if (allocated(chl))     deallocate(chl)
      if (allocated(cor))     deallocate(cor)
#if NDIMS == 3
      if (allocated(faces))   deallocate(faces)
#endif /* NDIMS == 3 */
      if (allocated(edges))   deallocate(edges)
      if (allocated(corners)) deallocate(corners)

! close the group
!
      call h5gclose_f(gid, err)

! check if the group has been closed successfuly
!
      if (err > 0) then

! print error about the problem with closing the group
!
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot close the metablock group!"

      end if

    else

! print error about the problem with opening the group
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot open the metablock group!"

    end if

!-------------------------------------------------------------------------------
!
  end subroutine read_metablocks_h5
!
!===============================================================================
!
! subroutine WRITE_DATABLOCKS_H5:
! ------------------------------
!
!   Subroutine writes all data block fields in the new group 'datablocks'
!   in the provided handler to the HDF5 file.
!
!   Arguments:
!
!     fid - the HDF5 file identifier;
!
!===============================================================================
!
  subroutine write_datablocks_h5(fid)

! import external procedures and variables
!
    use blocks         , only : ndims
    use blocks         , only : block_meta, block_data, list_data
    use blocks         , only : get_dblocks
    use coordinates    , only : nn => bcells
    use equations      , only : nv
    use hdf5           , only : hid_t, hsize_t
    use hdf5           , only : h5gcreate_f, h5gclose_f
    use iso_fortran_env, only : error_unit

! local variables are not implicit by default
!
    implicit none

! subroutine variables
!
    integer(hid_t), intent(in) :: fid

! local pointers
!
    type(block_meta), pointer  :: pmeta
    type(block_data), pointer  :: pdata

! local variables
!
    character(len=16)          :: bname
    integer(hid_t)             :: gid, bid
    integer(kind=4)            :: l
    integer                    :: err

! local arrays
!
    integer(hsize_t), dimension(4) :: dm = 1

! local parameters
!
    character(len=*), parameter :: loc = 'IO::write_datablocks_h5()'
!
!-------------------------------------------------------------------------------
!
! create a new group for storing data blocks
!
    call h5gcreate_f(fid, 'datablocks', gid, err)

! check if the group has been created successfuly
!
    if (err >= 0) then

! store data blocks only if there is at least one belonging to
! the current process
!
      if (get_dblocks() > 0) then

! prepate the dimensions
!
        dm(1) = nv
        dm(2) = nn
        dm(3) = nn
#if NDIMS == 3
        dm(4) = nn
#endif /* NDIMS == 3 */

! reset the block counter
!
        l = 0

! associate the pointer with the first block in the data block list
!
        pdata => list_data

! iterate over all data blocks and fill in the arrays id, u, and q
!
        do while(associated(pdata))

! increase the block counter
!
          l = l + 1

! create name for the current block
!
          write(bname, "('dblk_', i11.11)") l

! create a group for storing the current data block fields
!
          call h5gcreate_f(gid, bname, bid, err)

! store the corresponding meta block index
!
          call write_attribute(bid, 'meta', pdata%meta%id)

! store the primitive and conservative variables
!
          call write_array(bid, 'pvar' , dm(:), pdata%q (:,:,:,:))
          call write_array(bid, 'cvar0', dm(:), pdata%u0(:,:,:,:))
          call write_array(bid, 'cvar1', dm(:), pdata%u1(:,:,:,:))

! close the block group
!
          call h5gclose_f(bid, err)

! associate the pointer with the next data block on the list
!
          pdata => pdata%next

        end do ! data blocks

      end if ! dblocks > 0

! close the group
!
      call h5gclose_f(gid, err)

! check if the group has been closed successfuly
!
      if (err > 0) then

! print error about the problem with closing the group
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot close the group!"

      end if

    else

! print error about the problem with creating the group
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create the group!"

    end if

!-------------------------------------------------------------------------------
!
  end subroutine write_datablocks_h5
!
!===============================================================================
!
! subroutine READ_DATABLOCKS_H5:
! -----------------------------
!
!   Subroutine reads all data blocks stored in the group 'datablocks' of
!   the provided handler to the HDF5 restart file.
!
!   Arguments:
!
!     fid - the HDF5 file identifier;
!
!===============================================================================
!
  subroutine read_datablocks_h5(fid)

! import external procedures and variables
!
    use blocks         , only : block_meta, block_data, list_data
    use blocks         , only : append_datablock, link_blocks
    use coordinates    , only : nn => bcells
    use equations      , only : nv
    use hdf5           , only : hid_t, hsize_t
    use hdf5           , only : h5gopen_f, h5gclose_f, h5lexists_f
    use iso_fortran_env, only : error_unit

! local variables are not implicit by default
!
    implicit none

! subroutine variables
!
    integer(hid_t), intent(in)     :: fid

! local pointers
!
    type(block_data), pointer      :: pdata

! local variables
!
    logical                        :: flag
    character(len=16)              :: bname
    integer(hid_t)                 :: gid, bid
    integer(kind=4)                :: l, i
    integer                        :: dblocks, ierr

! local arrays
!
    integer(hsize_t), dimension(5) :: dm = 1

! local allocatable arrays
!
    integer(kind=4), dimension(:)        , allocatable :: id
    real(kind=8)   , dimension(:,:,:,:,:), allocatable :: uv, qv

! local parameters
!
    character(len=*), parameter :: loc = 'IO::read_datablocks_h5()'
!
!-------------------------------------------------------------------------------
!
! read the number of data blocks
!
    call h5gopen_f(fid, 'attributes', gid, ierr)
    if (ierr /= 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot open the attribute group!"
      return
    end if
    call read_attribute(gid, 'dblocks', dblocks)
    call h5gclose_f(gid, ierr)
    if (ierr /= 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot close the attribute group!"
      return
    end if

! restore data blocks only if there are any
!
    if (dblocks > 0) then

! open the group 'datablocks'
!
      call h5gopen_f(fid, 'datablocks', gid, ierr)
      if (ierr /= 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot open the data block group!"
        return
      end if

! fill out dimensions dm(:)
!
      dm(1) = dblocks
      dm(2) = nv
      dm(3) = nn
      dm(4) = nn
#if NDIMS == 3
      dm(5) = nn
#endif /* NDIMS == 3 */

! check if the old data format is used, otherwise use the new one
!
      call h5lexists_f(gid, 'meta', flag, ierr)

! restart files are in the old format
!
      if (flag) then

! allocate arrays to read data
!
        allocate(id(dm(1)))
        allocate(uv(dm(1),dm(2),dm(3),dm(4),dm(5)))
        allocate(qv(dm(1),dm(2),dm(3),dm(4),dm(5)))

! read array data from the HDF5 file
!
        call read_array(gid, 'meta', dm(1:1), id(:)        )
        call read_array(gid, 'uvar', dm(1:5), uv(:,:,:,:,:))
        call read_array(gid, 'qvar', dm(1:5), qv(:,:,:,:,:))

! iterate over data blocks, allocate them and fill out their fields
!
        do l = 1, dm(1)

! allocate and append to the end of the list a new datablock
!
          call append_datablock(pdata)

! associate the corresponding meta block with the current data block
!
          call link_blocks(block_array(id(l))%ptr, pdata)

! fill out the array of conservative and primitive variables
!
          pdata%u(:,:,:,:) = uv(l,:,:,:,:)
          pdata%q(:,:,:,:) = qv(l,:,:,:,:)

        end do ! l = 1, dm(1)

! deallocate allocatable arrays
!
        if (allocated(id)) deallocate(id)
        if (allocated(uv)) deallocate(uv)
        if (allocated(qv)) deallocate(qv)

! restart files are in the new format
!
      else ! flag

! iterate over data blocks, allocate them and fill out their fields
!
        do l = 1, dm(1)

! allocate and append to the end of the list a new datablock
!
          call append_datablock(pdata)

! create name for the current block
!
          write(bname, "('dblk_', i11.11)") l

! open the group for the current block fields
!
          call h5gopen_f(gid, bname, bid, ierr)

! get the id of the linked meta block
!
          call read_attribute(bid, 'meta', i)

! associate the corresponding meta block with the current data block
!
          call link_blocks(block_array(i)%ptr, pdata)

! fill out the array of primitive and conservative variables
!
          call read_array(bid, 'pvar' , dm(2:5), pdata%q (:,:,:,:))
          call read_array(bid, 'cvar0', dm(2:5), pdata%u0(:,:,:,:))
          call read_array(bid, 'cvar1', dm(2:5), pdata%u1(:,:,:,:))

! close the current data block group
!
          call h5gclose_f(bid, ierr)

        end do ! l = 1, dm(1)

      end if ! flag

! close the data block group
!
      call h5gclose_f(gid, ierr)
      if (ierr /= 0) then
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot close the data block group!"
        return
      end if

    end if ! dblocks > 0

!-------------------------------------------------------------------------------
!
  end subroutine read_datablocks_h5
!
!===============================================================================
!
! write_coordinates_h5: subroutine writes data block coordinates and other
!                       variables which determine geometrical position of
!                       the blocks
!
! info: this subroutine stores coordinates
!
! arguments:
!   fid - the HDF5 file identifier;
!
!===============================================================================
!
  subroutine write_coordinates_h5(fid)

! references to other modules
!
    use blocks         , only : block_meta, block_data, list_data
    use blocks         , only : nsides
    use blocks         , only : get_dblocks
    use coordinates    , only : maxlev
    use coordinates    , only : adx, ady, adz
    use hdf5           , only : hid_t, hsize_t
    use hdf5           , only : h5gcreate_f, h5gclose_f
    use iso_fortran_env, only : error_unit

! declare variables
!
    implicit none

! input variables
!
    integer(hid_t), intent(in) :: fid

! HDF5 variables
!
    integer(hid_t)    :: gid

! local variables
!
    integer           :: err
    integer(kind=4)   :: l
    integer(hsize_t)  :: am(1), cm(2), rm(2), dm(3)

! local allocatable arrays
!
    integer(kind=4), dimension(:)    , allocatable :: ids, lev, ref
    integer(kind=4), dimension(:,:)  , allocatable :: cor
    real   (kind=8), dimension(:,:,:), allocatable :: bnd

! local pointers
!
    type(block_data), pointer :: pdata

! local parameters
!
    character(len=*), parameter :: loc = 'IO::write_coordinates_h5()'
!
!-------------------------------------------------------------------------------
!
! create a group to store global attributes
!
    call h5gcreate_f(fid, 'coordinates', gid, err)

! check if the group has been created successfuly
!
    if (err .ge. 0) then

! store coordinates only if there are some data blocks on the current processor
!
      if (get_dblocks() .gt. 0) then

! prepare dimensions
!
        am(1) = maxlev
        cm(1) = get_dblocks()
        cm(2) = NDIMS
        rm(1) = maxlev
        rm(2) = NDIMS
        dm(1) = get_dblocks()
        dm(2) = NDIMS
        dm(3) = nsides

! allocate arrays to store coordinates
!
        allocate(ids(cm(1)))
        allocate(lev(cm(1)))
        allocate(ref(cm(1)))
        allocate(cor(cm(1),cm(2)))
        allocate(bnd(dm(1),dm(2),dm(3)))

! iterate over all data blocks and fill in the arrays
!
        l = 1
        pdata => list_data
        do while(associated(pdata))

! fill in the IDs array
!
          ids(l)     = pdata%meta%id

! fill in the level array
!
          lev(l)     = pdata%meta%level

! fill in the refinement flag
!
          ref(l)     = pdata%meta%refine

! fill in the coordinate array
!
          cor(l,:)   = pdata%meta%coords(:)

! fill in the bounds array
!
          bnd(l,1,1) = pdata%meta%xmin
          bnd(l,1,2) = pdata%meta%xmax
          bnd(l,2,1) = pdata%meta%ymin
          bnd(l,2,2) = pdata%meta%ymax
#if NDIMS == 3
          bnd(l,3,1) = pdata%meta%zmin
          bnd(l,3,2) = pdata%meta%zmax
#endif /* NDIMS == 3 */

          l = l + 1
          pdata => pdata%next
        end do

! write the arrays to the HDF5 file
!
        call write_array(gid, 'ids'   , cm(1), ids)
        call write_array(gid, 'levels', cm(1), lev)
        call write_array(gid, 'refine', cm(1), ref)
        call write_array(gid, 'coords', cm(:), cor)
        call write_array(gid, 'bounds', dm(:), bnd)
        call write_array(gid, 'dx'    , am(1), adx(1:maxlev))
        call write_array(gid, 'dy'    , am(1), ady(1:maxlev))
        call write_array(gid, 'dz'    , am(1), adz(1:maxlev))

! deallocate temporary arrays
!
        if (allocated(ids)) deallocate(ids)
        if (allocated(lev)) deallocate(lev)
        if (allocated(ref)) deallocate(ref)
        if (allocated(cor)) deallocate(cor)
        if (allocated(bnd)) deallocate(bnd)

      end if ! dblocks > 0

! close the attribute group
!
      call h5gclose_f(gid, err)

! check if the group has been closed successfuly
!
      if (err .gt. 0) then

! print error about the problem with closing the group
!
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot close the group!"

      end if

    else

! print error about the problem with creating the group
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create the group!"

    end if

!-------------------------------------------------------------------------------
!
  end subroutine write_coordinates_h5
!
!===============================================================================
!
! subroutine WRITE_PRIMITIVE_VARIABLES_H5:
! ---------------------------------------
!
!   Subroutine groups each primitive variable from all data blocks and writes
!   it as an array in the HDF5 dataset connected to the input HDF file
!   identifier.
!
!   Arguments:
!
!     fid - the HDF5 file identifier;
!
!===============================================================================
!
  subroutine write_primitive_variables_h5(fid)

! references to other modules
!
    use blocks         , only : block_data, list_data
    use blocks         , only : get_dblocks
    use coordinates    , only : ni => ncells, nn => bcells
    use coordinates    , only : nb, ne
    use equations      , only : nv, pvars
    use hdf5           , only : hid_t, hsize_t
    use hdf5           , only : h5gcreate_f, h5gclose_f
    use iso_fortran_env, only : error_unit

! local variables are not implicit by default
!
    implicit none

! subroutine arguments
!
    integer(hid_t), intent(in) :: fid

! HDF5 variables
!
    integer(hid_t)             :: gid
    integer(hsize_t)           :: dm(4) = 1

! local variables
!
    integer         :: err
    integer(kind=4) :: i, j, k, l, n
    integer(kind=4) :: il, jl, kl = 1
    integer(kind=4) :: iu, ju, ku = 1

! local allocatable arrays
!
    real(kind=8), dimension(:,:,:,:), allocatable :: qarr

! local pointers
!
    type(block_data), pointer :: pdata

! local parameters
!
    character(len=*), parameter :: loc = 'IO::write_primitive_variables_h5()'
!
!-------------------------------------------------------------------------------
!
! create a group to store variables
!
    call h5gcreate_f(fid, 'variables', gid, err)

! check if the group was created successfuly
!
    if (err >= 0) then

! store variables only if there is at least one data block associated with
! the current process
!
      if (get_dblocks() > 0) then

! prepare dimensions and index limits
!
        dm(1) = get_dblocks()
        if (with_ghosts) then
          dm(2) = nn
          dm(3) = nn
#if NDIMS == 3
          dm(4) = nn
#endif /* NDIMS == 3 */

          il    =  1
          jl    =  1
#if NDIMS == 3
          kl    =  1
#endif /* NDIMS == 3 */
          iu    = nn
          ju    = nn
#if NDIMS == 3
          ku    = nn
#endif /* NDIMS == 3 */
        else
          dm(2) = ni
          dm(3) = ni
#if NDIMS == 3
          dm(4) = ni
#endif /* NDIMS == 3 */

          il    = nb
          jl    = nb
#if NDIMS == 3
          kl    = nb
#endif /* NDIMS == 3 */
          iu    = ne
          ju    = ne
#if NDIMS == 3
          ku    = ne
#endif /* NDIMS == 3 */
        end if

! allocate array to group a variable from all data blocks
!
        allocate(qarr(dm(1),dm(2),dm(3),dm(4)))

! iterate over all variables
!
        do n = 1, nv

! reset the block counter
!
          l = 0

! assosiate the block pointer with the first data block on the list
!
          pdata => list_data

! iterate over all data blocks and copy the variable from each of them to
! the allocate array
!
          do while(associated(pdata))

! increase the data block counter
!
            l = l + 1

! copy the variable from the current data block
!
            qarr(l,1:dm(2),1:dm(3),1:dm(4)) = pdata%q(n,il:iu,jl:ju,kl:ku)

! assign the pointer with the next data block on the list
!
            pdata => pdata%next

          end do ! pdata=>list_data

! write the variable array to the HDF5 file
!
          call write_array(gid, trim(pvars(n)), dm, qarr)

        end do ! n = 1, nv

! deallocate allocatable array
!
        if (allocated(qarr)) deallocate(qarr)

      end if ! dblocks > 0

! close the variable group
!
      call h5gclose_f(gid, err)

! check if the group has been closed successfuly
!
      if (err > 0) then

! print error about the problem with closing the group
!
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot close the group!"

      end if

    else ! error with creating a group

! print error about the problem with creating the group
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create the group!"

    end if

!-------------------------------------------------------------------------------
!
  end subroutine write_primitive_variables_h5
!
!===============================================================================
!
! subroutine WRITE_CONSERVATIVE_VARIABLES_H5:
! ------------------------------------------
!
!   Subroutine groups each conservative variable from all data blocks and writes
!   it as an array in the HDF5 dataset connected to the input HDF file
!   identifier.
!
!   Arguments:
!
!     fid - the HDF5 file identifier;
!
!===============================================================================
!
  subroutine write_conservative_variables_h5(fid)

! references to other modules
!
    use blocks         , only : block_data, list_data
    use blocks         , only : get_dblocks
    use coordinates    , only : ni => ncells, nn => bcells
    use coordinates    , only : nb, ne
    use equations      , only : nv, cvars
    use hdf5           , only : hid_t, hsize_t
    use hdf5           , only : h5gcreate_f, h5gclose_f
    use iso_fortran_env, only : error_unit

! local variables are not implicit by default
!
    implicit none

! subroutine arguments
!
    integer(hid_t), intent(in) :: fid

! HDF5 variables
!
    integer(hid_t)             :: gid
    integer(hsize_t)           :: dm(4) = 1

! local variables
!
    integer         :: err
    integer(kind=4) :: i, j, k, l, n
    integer(kind=4) :: il, jl, kl = 1
    integer(kind=4) :: iu, ju, ku = 1

! local allocatable arrays
!
    real(kind=8), dimension(:,:,:,:), allocatable :: qarr

! local pointers
!
    type(block_data), pointer :: pdata

! local parameters
!
    character(len=*), parameter :: loc = 'IO::write_conservative_variables_h5()'
!
!-------------------------------------------------------------------------------
!
! create a group to store variables
!
    call h5gcreate_f(fid, 'variables', gid, err)

! check if the group was created successfuly
!
    if (err >= 0) then

! store variables only if there is at least one data block associated with
! the current process
!
      if (get_dblocks() > 0) then

! prepare dimensions and index limits
!
        dm(1) = get_dblocks()
        if (with_ghosts) then
          dm(2) = nn
          dm(3) = nn
#if NDIMS == 3
          dm(4) = nn
#endif /* NDIMS == 3 */

          il    =  1
          jl    =  1
#if NDIMS == 3
          kl    =  1
#endif /* NDIMS == 3 */
          iu    = nn
          ju    = nn
#if NDIMS == 3
          ku    = nn
#endif /* NDIMS == 3 */
        else
          dm(2) = ni
          dm(3) = ni
#if NDIMS == 3
          dm(4) = ni
#endif /* NDIMS == 3 */

          il    = nb
          jl    = nb
#if NDIMS == 3
          kl    = nb
#endif /* NDIMS == 3 */
          iu    = ne
          ju    = ne
#if NDIMS == 3
          ku    = ne
#endif /* NDIMS == 3 */
        end if

! allocate array to group a variable from all data blocks
!
        allocate(qarr(dm(1),dm(2),dm(3),dm(4)))

! iterate over all variables
!
        do n = 1, nv

! reset the block counter
!
          l = 0

! assosiate the block pointer with the first data block on the list
!
          pdata => list_data

! iterate over all data blocks and copy the variable from each of them to
! the allocate array
!
          do while(associated(pdata))

! increase the data block counter
!
            l = l + 1

! copy the variable from the current data block
!
            qarr(l,1:dm(2),1:dm(3),1:dm(4)) = pdata%u(n,il:iu,jl:ju,kl:ku)

! assign the pointer with the next data block on the list
!
            pdata => pdata%next

          end do ! pdata=>list_data

! write the variable array to the HDF5 file
!
          call write_array(gid, trim(cvars(n)), dm, qarr)

        end do ! n = 1, nv

! deallocate allocatable array
!
        if (allocated(qarr)) deallocate(qarr)

      end if ! dblocks > 0

! close the variable group
!
      call h5gclose_f(gid, err)

! check if the group has been closed successfuly
!
      if (err > 0) then

! print error about the problem with closing the group
!
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot close the group!"

      end if

    else ! error with creating a group

! print error about the problem with creating the group
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create the group!"

    end if

!-------------------------------------------------------------------------------
!
  end subroutine write_conservative_variables_h5
!
!===============================================================================
!
! WRITE_ATTRIBUTE SUBROUTINES
!
!===============================================================================
!
! subroutine WRITE_SCALAR_ATTRIBUTE_STRING_H5:
! --------------------------------------------
!
!   Subroutine stores a value of the string attribute in the group provided
!   by an identifier and the attribute name.
!
!   Arguments:
!
!     gid    - the group identifier to which the attribute should be linked;
!     aname  - the attribute name;
!     avalue - the attribute value;
!
!===============================================================================
!
  subroutine write_scalar_attribute_string_h5(gid, aname, avalue)

! import procedures and variables from other modules
!
    use hdf5           , only : H5T_NATIVE_CHARACTER
    use hdf5           , only : hid_t, hsize_t
    use hdf5           , only : h5screate_simple_f, h5sclose_f
    use hdf5           , only : h5acreate_f, h5awrite_f, h5aclose_f
    use hdf5           , only : h5tcopy_f, h5tset_size_f, h5tclose_f
    use iso_fortran_env, only : error_unit

! local variables are not implicit by default
!
    implicit none

! subroutine arguments
!
    integer(hid_t)  , intent(in) :: gid
    character(len=*), intent(in) :: aname
    character(len=*), intent(in) :: avalue

! local variables
!
    integer(hid_t)                 :: sid, aid, atype
    integer(hsize_t)               :: alen
    integer(hsize_t), dimension(1) :: am = (/ 1 /)
    integer                        :: ierr

! local parameters
!
    character(len=*), parameter :: loc = 'IO::write_scalar_attribute_string_h5()'
!
!-------------------------------------------------------------------------------
!
! copy the attribute type and set its size
!
    call h5tcopy_f(H5T_NATIVE_CHARACTER, atype, ierr)
    if (ierr /= 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot copy type for attribute :" // trim(aname)
      return
    end if

! get the string length
!
    alen = len(trim(avalue))
    if (alen <= 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "String attribute has wrong length:" // trim(aname)
      return
    end if

! set the attribute type size
!
    call h5tset_size_f(atype, alen, ierr)
    if (ierr /= 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot set the type size for attribute :" // trim(aname)
      return
    end if

! create space for the attribute value
!
    call h5screate_simple_f(1, am, sid, ierr)
    if (ierr /= 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create space for attribute :" // trim(aname)
      return
    end if

! create the attribute in the given group
!
    call h5acreate_f(gid, aname, atype, sid, aid, ierr)
    if (ierr == 0) then

! write the attribute data
!
      call h5awrite_f(aid, atype, trim(avalue), am, ierr)
      if (ierr /= 0) then
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot write the attribute data in :" // trim(aname)
      end if

! close the attribute
!
      call h5aclose_f(aid, ierr)
      if (ierr /= 0) then
        write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                        , "Cannot close attribute :" // trim(aname)
      end if

    else
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create attribute :" // trim(aname)
    end if

! release the space
!
    call h5sclose_f(sid, ierr)
    if (ierr /= 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot close space for attribute :" // trim(aname)
    end if

!-------------------------------------------------------------------------------
!
  end subroutine write_scalar_attribute_string_h5
!
!===============================================================================
!
! subroutine WRITE_SCALAR_ATTRIBUTE_INTEGER_H5:
! --------------------------------------------
!
!   Subroutine stores a value of the integer attribute in the group provided
!   by an identifier and the attribute name.
!
!   Arguments:
!
!     gid    - the group identifier to which the attribute should be linked;
!     aname  - the attribute name;
!     avalue - the attribute value;
!
!===============================================================================
!
  subroutine write_scalar_attribute_integer_h5(gid, aname, avalue)

! import procedures and variables from other modules
!
    use hdf5           , only : H5T_NATIVE_INTEGER
    use hdf5           , only : hid_t, hsize_t
    use hdf5           , only : h5screate_simple_f, h5sclose_f
    use hdf5           , only : h5acreate_f, h5awrite_f, h5aclose_f
    use iso_fortran_env, only : error_unit

! local variables are not implicit by default
!
    implicit none

! subroutine arguments
!
    integer(hid_t)  , intent(in) :: gid
    character(len=*), intent(in) :: aname
    integer(kind=4) , intent(in) :: avalue

! local variables
!
    integer(hid_t)                 :: sid, aid
    integer(hsize_t), dimension(1) :: am = (/ 1 /)
    integer                        :: ierr

! local parameters
!
    character(len=*), parameter :: loc = 'IO::write_scalar_attribute_integer_h5()'
!
!-------------------------------------------------------------------------------
!
! create space for the attribute value
!
    call h5screate_simple_f(1, am, sid, ierr)
    if (ierr /= 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create space for attribute :" // trim(aname)
      return
    end if

! create the attribute in the given group
!
    call h5acreate_f(gid, aname, H5T_NATIVE_INTEGER, sid, aid, ierr)
    if (ierr == 0) then

! write the attribute data
!
      call h5awrite_f(aid, H5T_NATIVE_INTEGER, avalue, am, ierr)
      if (ierr /= 0) then
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot write the attribute data in :" // trim(aname)
      end if

! close the attribute
!
      call h5aclose_f(aid, ierr)
      if (ierr /= 0) then
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot close attribute :" // trim(aname)
      end if

    else
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create attribute :" // trim(aname)
    end if

! release the space
!
    call h5sclose_f(sid, ierr)
    if (ierr /= 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot close space for attribute :" // trim(aname)
    end if

!-------------------------------------------------------------------------------
!
  end subroutine write_scalar_attribute_integer_h5
!
!===============================================================================
!
! subroutine WRITE_SCALAR_ATTRIBUTE_DOUBLE_H5:
! -------------------------------------------
!
!   Subroutine stores a value of the double precision attribute in the group
!   provided by an identifier and the attribute name.
!
!   Arguments:
!
!     gid    - the group identifier to which the attribute should be linked;
!     aname  - the attribute name;
!     avalue - the attribute value;
!
!===============================================================================
!
  subroutine write_scalar_attribute_double_h5(gid, aname, avalue)

! import procedures and variables from other modules
!
    use hdf5           , only : H5T_NATIVE_DOUBLE
    use hdf5           , only : hid_t, hsize_t
    use hdf5           , only : h5screate_simple_f, h5sclose_f
    use hdf5           , only : h5acreate_f, h5awrite_f, h5aclose_f
    use iso_fortran_env, only : error_unit

! local variables are not implicit by default
!
    implicit none

! attribute arguments
!
    integer(hid_t)  , intent(in) :: gid
    character(len=*), intent(in) :: aname
    real(kind=8)    , intent(in) :: avalue

! local variables
!
    integer(hid_t)                 :: sid, aid
    integer(hsize_t), dimension(1) :: am = (/ 1 /)
    integer                        :: ierr

! local parameters
!
    character(len=*), parameter :: loc = 'IO::write_scalar_attribute_double_h5()'
!
!-------------------------------------------------------------------------------
!
! create space for the attribute value
!
    call h5screate_simple_f(1, am, sid, ierr)
    if (ierr /= 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create space for attribute :" // trim(aname)
      return
    end if

! create the attribute in the given group
!
    call h5acreate_f(gid, aname, H5T_NATIVE_DOUBLE, sid, aid, ierr)
    if (ierr == 0) then

! write the attribute data
!
      call h5awrite_f(aid, H5T_NATIVE_DOUBLE, avalue, am, ierr)
      if (ierr /= 0) then
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot write the attribute data in :" // trim(aname)
      end if

! close the attribute
!
      call h5aclose_f(aid, ierr)
      if (ierr /= 0) then
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot close attribute :" // trim(aname)
      end if

    else
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create attribute :" // trim(aname)
    end if

! release the space
!
    call h5sclose_f(sid, ierr)
    if (ierr /= 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot close space for attribute :" // trim(aname)
    end if

!-------------------------------------------------------------------------------
!
  end subroutine write_scalar_attribute_double_h5
!
!===============================================================================
!
! subroutine WRITE_VECTOR_ATTRIBUTE_INTEGER_H5:
! --------------------------------------------
!
!   Subroutine stores a vector of the integer attribute in the group provided
!   by an identifier and the attribute name.
!
!   Arguments:
!
!     gid    - the group identifier to which the attribute should be linked;
!     aname  - the attribute name;
!     avalue - the attribute values;
!
!===============================================================================
!
  subroutine write_vector_attribute_integer_h5(gid, aname, avalue)

! import procedures and variables from other modules
!
    use hdf5           , only : hid_t, hsize_t, H5T_NATIVE_INTEGER
    use hdf5           , only : h5screate_simple_f, h5sclose_f
    use hdf5           , only : h5acreate_f, h5awrite_f, h5aclose_f
    use iso_fortran_env, only : error_unit

! local variables are not implicit by default
!
    implicit none

! attribute arguments
!
    integer(hid_t)                , intent(in) :: gid
    character(len=*)              , intent(in) :: aname
    integer(kind=4) , dimension(:), intent(in) :: avalue

! local variables
!
    integer(hid_t)                 :: sid, aid
    integer(hsize_t), dimension(1) :: am = (/ 1 /)
    integer                        :: ierr

! local parameters
!
    character(len=*), parameter :: loc = 'IO::write_vector_attribute_integer_h5()'
!
!-------------------------------------------------------------------------------
!
! set the proper attribute length
!
    am(1) = size(avalue)

! create space for the attribute value
!
    call h5screate_simple_f(1, am, sid, ierr)
    if (ierr /= 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create space for attribute :" // trim(aname)
      return
    end if

! create the attribute in the given group
!
    call h5acreate_f(gid, aname, H5T_NATIVE_INTEGER, sid, aid, ierr)
    if (ierr == 0) then

! write the attribute data
!
      call h5awrite_f(aid, H5T_NATIVE_INTEGER, avalue, am, ierr)
      if (ierr /= 0) then
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot write the attribute data in :" // trim(aname)
      end if

! close the attribute
!
      call h5aclose_f(aid, ierr)
      if (ierr /= 0) then
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot close attribute :" // trim(aname)
      end if

    else
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create attribute :" // trim(aname)
    end if

! release the space
!
    call h5sclose_f(sid, ierr)
    if (ierr /= 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot close space for attribute :" // trim(aname)
    end if

!-------------------------------------------------------------------------------
!
  end subroutine write_vector_attribute_integer_h5
!
!===============================================================================
!
! subroutine WRITE_VECTOR_ATTRIBUTE_DOUBLE_H5:
! -------------------------------------------
!
!   Subroutine stores a vector of the double precision attribute in the group
!   provided by an identifier and the attribute name.
!
!   Arguments:
!
!     gid    - the group identifier to which the attribute should be linked;
!     aname  - the attribute name;
!     avalue - the attribute values;
!
!===============================================================================
!
  subroutine write_vector_attribute_double_h5(gid, aname, avalue)

! import procedures and variables from other modules
!
    use hdf5           , only : H5T_NATIVE_DOUBLE
    use hdf5           , only : hid_t, hsize_t
    use hdf5           , only : h5screate_simple_f, h5sclose_f
    use hdf5           , only : h5acreate_f, h5awrite_f, h5aclose_f
    use iso_fortran_env, only : error_unit

! local variables are not implicit by default
!
    implicit none

! attribute arguments
!
    integer(hid_t)                , intent(in) :: gid
    character(len=*)              , intent(in) :: aname
    real(kind=8)    , dimension(:), intent(in) :: avalue

! local variables
!
    integer(hid_t)                 :: sid, aid
    integer(hsize_t), dimension(1) :: am = (/ 1 /)
    integer                        :: ierr

! local parameters
!
    character(len=*), parameter :: loc = 'IO::write_vector_attribute_double_h5()'
!
!-------------------------------------------------------------------------------
!
! set the proper attribute length
!
    am(1) = size(avalue)

! create space for the attribute value
!
    call h5screate_simple_f(1, am, sid, ierr)
    if (ierr /= 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create space for attribute :" // trim(aname)
      return
    end if

! create the attribute in the given group
!
    call h5acreate_f(gid, aname, H5T_NATIVE_DOUBLE, sid, aid, ierr)
    if (ierr == 0) then

! write the attribute data
!
      call h5awrite_f(aid, H5T_NATIVE_DOUBLE, avalue, am, ierr)
      if (ierr /= 0) then
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot write the attribute data in :" // trim(aname)
      end if

! close the attribute
!
      call h5aclose_f(aid, ierr)
      if (ierr /= 0) then
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot close attribute :" // trim(aname)
      end if

    else
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create attribute :" // trim(aname)
    end if

! release the space
!
    call h5sclose_f(sid, ierr)
    if (ierr /= 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot close space for attribute :" // trim(aname)
    end if

!-------------------------------------------------------------------------------
!
  end subroutine write_vector_attribute_double_h5

!===============================================================================
!
! READ_ATTRIBUTE SUBROUTINES
!
!===============================================================================
!
! subroutine READ_SCALAR_ATTRIBUTE_INTEGER_H5:
! -------------------------------------------
!
!   Subroutine reads a value of the integer attribute provided by the group
!   identifier to which it is linked and its name.
!
!   Arguments:
!
!     gid    - the group identifier to which the attribute is linked;
!     aname  - the attribute name;
!     avalue - the attribute value;
!
!===============================================================================
!
  subroutine read_scalar_attribute_integer_h5(gid, aname, avalue)

! import procedures and variables from other modules
!
    use hdf5           , only : H5T_NATIVE_INTEGER
    use hdf5           , only : hid_t, hsize_t
    use hdf5           , only : h5aexists_by_name_f
    use hdf5           , only : h5aopen_by_name_f, h5aread_f, h5aclose_f
    use iso_fortran_env, only : error_unit

! local variables are not implicit by default
!
    implicit none

! attribute arguments
!
    integer(hid_t)  , intent(in)    :: gid
    character(len=*), intent(in)    :: aname
    integer(kind=4) , intent(inout) :: avalue

! local variables
!
    logical                         :: exists = .false.
    integer(hid_t)                  :: aid
    integer(hsize_t), dimension(1)  :: am     = (/ 1 /)
    integer                         :: ierr

! local parameters
!
    character(len=*), parameter :: loc = 'IO::read_scalar_attribute_integer_h5()'
!
!-------------------------------------------------------------------------------
!
! check if the attribute exists in the group provided by gid
!
    call h5aexists_by_name_f(gid, '.', aname, exists, ierr)
    if (ierr /= 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot check if attribute exists :" // trim(aname)
      return
    end if
    if (.not. exists) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Attribute does not exist :" // trim(aname)
      return
    end if

! open the attribute
!
    call h5aopen_by_name_f(gid, '.', aname, aid, ierr)
    if (ierr /= 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot open attribute :" // trim(aname)
      return
    end if

! read attribute value
!
    call h5aread_f(aid, H5T_NATIVE_INTEGER, avalue, am(:), ierr)
    if (ierr /= 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot read attribute :" // trim(aname)
    end if

! close the attribute
!
    call h5aclose_f(aid, ierr)
    if (ierr /= 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot close attribute :" // trim(aname)
      return
    end if

!-------------------------------------------------------------------------------
!
  end subroutine read_scalar_attribute_integer_h5
!
!===============================================================================
!
! subroutine READ_SCALAR_ATTRIBUTE_DOUBLE_H5:
! ------------------------------------------
!
!   Subroutine reads a value of the double precision attribute provided by
!   the group identifier to which it is linked and its name.
!
!   Arguments:
!
!     gid    - the group identifier to which the attribute is linked;
!     aname  - the attribute name;
!     avalue - the attribute value;
!
!===============================================================================
!
  subroutine read_scalar_attribute_double_h5(gid, aname, avalue)

! import procedures and variables from other modules
!
    use hdf5           , only : H5T_NATIVE_DOUBLE
    use hdf5           , only : hid_t, hsize_t
    use hdf5           , only : h5aexists_by_name_f
    use hdf5           , only : h5aopen_by_name_f, h5aread_f, h5aclose_f
    use iso_fortran_env, only : error_unit

! local variables are not implicit by default
!
    implicit none

! attribute arguments
!
    integer(hid_t)  , intent(in)    :: gid
    character(len=*), intent(in)    :: aname
    real(kind=8)    , intent(inout) :: avalue

! local variables
!
    logical                         :: exists = .false.
    integer(hid_t)                  :: aid
    integer(hsize_t), dimension(1)  :: am     = (/ 1 /)
    integer                         :: ierr

! local parameters
!
    character(len=*), parameter :: loc = 'IO::read_scalar_attribute_double_h5()'
!
!-------------------------------------------------------------------------------
!
! check if the attribute exists in the group provided by gid
!
    call h5aexists_by_name_f(gid, '.', aname, exists, ierr)
    if (ierr /= 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot check if attribute exists :" // trim(aname)
      return
    end if
    if (.not. exists) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Attribute does not exist :" // trim(aname)
      return
    end if

! open the attribute
!
    call h5aopen_by_name_f(gid, '.', aname, aid, ierr)
    if (ierr /= 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot open attribute :" // trim(aname)
      return
    end if

! read attribute value
!
    call h5aread_f(aid, H5T_NATIVE_DOUBLE, avalue, am(:), ierr)
    if (ierr /= 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot read attribute :" // trim(aname)
    end if

! close the attribute
!
    call h5aclose_f(aid, ierr)
    if (ierr /= 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot close attribute :" // trim(aname)
      return
    end if

!-------------------------------------------------------------------------------
!
  end subroutine read_scalar_attribute_double_h5
!
!===============================================================================
!
! subroutine READ_VECTOR_ATTRIBUTE_INTEGER_H5:
! -------------------------------------------
!
!   Subroutine reads a vector of the integer attribute provided by the group
!   identifier to which it is linked and its name.
!
!   Arguments:
!
!     gid    - the group identifier to which the attribute is linked;
!     aname  - the attribute name;
!     avalue - the attribute value;
!
!===============================================================================
!
  subroutine read_vector_attribute_integer_h5(gid, aname, avalue)

! import procedures and variables from other modules
!
    use hdf5           , only : H5T_NATIVE_INTEGER
    use hdf5           , only : hid_t, hsize_t
    use hdf5           , only : h5aexists_by_name_f, h5aget_space_f
    use hdf5           , only : h5aopen_by_name_f, h5aread_f, h5aclose_f
    use hdf5           , only : h5sclose_f, h5sget_simple_extent_dims_f
    use iso_fortran_env, only : error_unit

! local variables are not implicit by default
!
    implicit none

! attribute arguments
!
    integer(hid_t)                , intent(in)    :: gid
    character(len=*)              , intent(in)    :: aname
    integer(kind=4) , dimension(:), intent(inout) :: avalue

! local variables
!
    logical                         :: exists = .false.
    integer(hid_t)                  :: aid, sid
    integer(hsize_t), dimension(1)  :: am, bm
    integer(hsize_t)                :: alen
    integer                         :: ierr

! local parameters
!
    character(len=*), parameter :: loc = 'IO::read_vector_attribute_integer_h5()'
!
!-------------------------------------------------------------------------------
!
! check if the attribute exists in the group provided by gid
!
    call h5aexists_by_name_f(gid, '.', aname, exists, ierr)
    if (ierr /= 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot check if attribute exists :" // trim(aname)
      return
    end if
    if (.not. exists) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Attribute does not exist :" // trim(aname)
      return
    end if

! open the attribute
!
    call h5aopen_by_name_f(gid, '.', aname, aid, ierr)
    if (ierr /= 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot open attribute :" // trim(aname)
      return
    end if

! get the attribute space
!
    call h5aget_space_f(aid, sid, ierr)
    if (ierr == 0) then
      call h5sget_simple_extent_dims_f(sid, am, bm, ierr)
      if (ierr /= 1) then
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot get attribute dimensions :" // trim(aname)
      end if
      call h5sclose_f(sid, ierr)
      if (ierr /= 0) then
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot close the attribute space :" // trim(aname)
      end if
    else
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot get the attribute space :" // trim(aname)
      return
    end if

! check if the output array is large enough
!
    if (am(1) > size(avalue)) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Attribute too large for output argument :" // trim(aname)
      return
    end if

! read attribute value
!
    call h5aread_f(aid, H5T_NATIVE_INTEGER, avalue, am(:), ierr)
    if (ierr /= 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot read attribute :" // trim(aname)
    end if

! close the attribute
!
    call h5aclose_f(aid, ierr)
    if (ierr /= 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot close attribute :" // trim(aname)
      return
    end if

!-------------------------------------------------------------------------------
!
  end subroutine read_vector_attribute_integer_h5
!
!===============================================================================
!
! subroutine READ_VECTOR_ATTRIBUTE_DOUBLE_H5:
! ------------------------------------------
!
!   Subroutine reads a vector of the double precision attribute provided by
!   the group identifier to which it is linked and its name.
!
!   Arguments:
!
!     gid    - the group identifier to which the attribute is linked;
!     aname  - the attribute name;
!     avalue - the attribute value;
!
!===============================================================================
!
  subroutine read_vector_attribute_double_h5(gid, aname, avalue)

! import procedures and variables from other modules
!
    use hdf5           , only : H5T_NATIVE_DOUBLE
    use hdf5           , only : hid_t, hsize_t
    use hdf5           , only : h5aexists_by_name_f, h5aget_space_f
    use hdf5           , only : h5aopen_by_name_f, h5aread_f, h5aclose_f
    use hdf5           , only : h5sclose_f, h5sget_simple_extent_dims_f
    use iso_fortran_env, only : error_unit

! local variables are not implicit by default
!
    implicit none

! attribute arguments
!
    integer(hid_t)                , intent(in)    :: gid
    character(len=*)              , intent(in)    :: aname
    real(kind=8)    , dimension(:), intent(inout) :: avalue

! local variables
!
    logical                         :: exists = .false.
    integer(hid_t)                  :: aid, sid
    integer(hsize_t), dimension(1)  :: am, bm
    integer(hsize_t)                :: alen
    integer                         :: ierr

! local parameters
!
    character(len=*), parameter :: loc = 'IO::read_vector_attribute_double_h5()'
!
!-------------------------------------------------------------------------------
!
! check if the attribute exists in the group provided by gid
!
    call h5aexists_by_name_f(gid, '.', aname, exists, ierr)
    if (ierr /= 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot check if attribute exists :" // trim(aname)
      return
    end if
    if (.not. exists) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Attribute does not exist :" // trim(aname)
      return
    end if

! open the attribute
!
    call h5aopen_by_name_f(gid, '.', aname, aid, ierr)
    if (ierr /= 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot open attribute :" // trim(aname)
      return
    end if

! get the attribute space
!
    call h5aget_space_f(aid, sid, ierr)
    if (ierr == 0) then
      call h5sget_simple_extent_dims_f(sid, am, bm, ierr)
      if (ierr /= 1) then
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot get attribute dimensions :" // trim(aname)
      end if
      call h5sclose_f(sid, ierr)
      if (ierr /= 0) then
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot close the attribute space :" // trim(aname)
      end if
    else
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot get the attribute space :" // trim(aname)
      return
    end if

! check if the output array is large enough
!
    if (am(1) > size(avalue)) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Attribute too large for output argument :" // trim(aname)
      return
    end if

! read attribute value
!
    call h5aread_f(aid, H5T_NATIVE_DOUBLE, avalue, am(:), ierr)
    if (ierr /= 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot read attribute :" // trim(aname)
    end if

! close the attribute
!
    call h5aclose_f(aid, ierr)
    if (ierr /= 0) then
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot close attribute :" // trim(aname)
      return
    end if

!-------------------------------------------------------------------------------
!
  end subroutine read_vector_attribute_double_h5
!
!===============================================================================
!
! WRITE_ARRAY SUBROUTINES
!
!===============================================================================
!
! subroutine WRITE_1D_ARRAY_INTEGER_H5:
! ------------------------------------
!
!   Subroutine stores a one-dimensional integer array in a group specified by
!   identifier.
!
!   Arguments:
!
!     gid   - the HDF5 group identifier
!     name  - the string name describing the array
!     dm    - the array dimensions
!     value - the array values
!
!===============================================================================
!
  subroutine write_1d_array_integer_h5(gid, name, ln, var)

! import procedures and variables from other modules
!
    use hdf5           , only : H5T_NATIVE_INTEGER
    use hdf5           , only : hid_t, hsize_t
    use hdf5           , only : h5screate_simple_f, h5sclose_f
    use hdf5           , only : h5dcreate_f, h5dwrite_f, h5dclose_f
    use hdf5           , only : h5pset_chunk_f
    use iso_fortran_env, only : error_unit

! local variables are not implicit by default
!
    implicit none

! subroutine arguments
!
    integer(hid_t)                , intent(in) :: gid
    character(len=*)              , intent(in) :: name
    integer(hsize_t)              , intent(in) :: ln
    integer(kind=4) , dimension(:), intent(in) :: var

! HDF5 object identifiers
!
    integer(hid_t) :: sid, did

! array dimensions
!
    integer(hsize_t), dimension(1) :: dm

! procedure return value
!
    integer        :: iret

! local parameters
!
    character(len=*), parameter :: loc = 'IO::write_1d_array_integer_h5()'
!
!-------------------------------------------------------------------------------
!
! substitute array dimensions
!
    dm(1) = ln

! create a space for the array
!
    call h5screate_simple_f(1, dm(1:1), sid, iret)

! check if the space has been created successfuly, if not quit
!
    if (iret < 0) then

! print error about the problem with creating the space
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create space for dataset: " // trim(name)

! quit the subroutine
!
      return

    end if

! set the chunk size
!
    call h5pset_chunk_f(pid, 1, dm(1:1), iret)

! check if the chunk size has been set properly
!
    if (iret > 0) then

! print error about the problem with setting the chunk size
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot set the size of the chunk!"

    end if

! create the dataset
!
    call h5dcreate_f(gid, name, H5T_NATIVE_INTEGER, sid, did, iret, pid)

! check if the dataset has been created successfuly
!
    if (iret >= 0) then

! write the dataset data
!
      call h5dwrite_f(did, H5T_NATIVE_INTEGER, var(:), dm(1:1), iret, sid)

! check if the dataset has been written successfuly
!
      if (iret > 0) then

! print error about the problem with writing down the dataset
!
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot write dataset: " // trim(name)

      end if

! close the dataset
!
      call h5dclose_f(did, iret)

! check if the dataset has been closed successfuly
!
      if (iret > 0) then

! print error about the problem with closing the dataset
!
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot close dataset: " // trim(name)

      end if

    else

! print error about the problem with creating the dataset
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create dataset: " // trim(name)

    end if

! release the space
!
    call h5sclose_f(sid, iret)

! check if the space has been released successfuly
!
    if (iret > 0) then

! print error about the problem with closing the space
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot close space for dataset: " // trim(name)

    end if

!-------------------------------------------------------------------------------
!
  end subroutine write_1d_array_integer_h5
!
!===============================================================================
!
! subroutine WRITE_2D_ARRAY_INTEGER_H5:
! ------------------------------------
!
!   Subroutine stores a two-dimensional integer array in a group specified by
!   identifier.
!
!   Arguments:
!
!     gid   - the HDF5 group identifier
!     name  - the string name describing the array
!     dm    - the array dimensions
!     value - the array values
!
!===============================================================================
!
  subroutine write_2d_array_integer_h5(gid, name, dm, var)

! import procedures and variables from other modules
!
    use hdf5           , only : H5T_NATIVE_INTEGER
    use hdf5           , only : hid_t, hsize_t
    use hdf5           , only : h5screate_simple_f, h5sclose_f
    use hdf5           , only : h5dcreate_f, h5dwrite_f, h5dclose_f
    use hdf5           , only : h5pset_chunk_f
    use iso_fortran_env, only : error_unit

! local variables are not implicit by default
!
    implicit none

! subroutine arguments
!
    integer(hid_t)                  , intent(in) :: gid
    character(len=*)                , intent(in) :: name
    integer(hsize_t), dimension(2)  , intent(in) :: dm
    integer(kind=4) , dimension(:,:), intent(in) :: var

! HDF5 object identifiers
!
    integer(hid_t) :: sid, did

! procedure return value
!
    integer        :: iret

! local parameters
!
    character(len=*), parameter :: loc = 'IO::write_2d_array_integer_h5()'
!
!-------------------------------------------------------------------------------
!
! create a space for the array
!
    call h5screate_simple_f(2, dm(1:2), sid, iret)

! check if the space has been created successfuly, if not quit
!
    if (iret < 0) then

! print error about the problem with creating the space
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create space for dataset: " // trim(name)

! quit the subroutine
!
      return

    end if

! set the chunk size
!
    call h5pset_chunk_f(pid, 2, dm(1:2), iret)

! check if the chunk size has been set properly
!
    if (iret > 0) then

! print error about the problem with setting the chunk size
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot set the size of the chunk!"

    end if

! create the dataset
!
    call h5dcreate_f(gid, name, H5T_NATIVE_INTEGER, sid, did, iret, pid)

! check if the dataset has been created successfuly
!
    if (iret >= 0) then

! write the dataset data
!
      call h5dwrite_f(did, H5T_NATIVE_INTEGER, var(:,:), dm(1:2), iret, sid)

! check if the dataset has been written successfuly
!
      if (iret > 0) then

! print error about the problem with writing down the dataset
!
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot write dataset: " // trim(name)

      end if

! close the dataset
!
      call h5dclose_f(did, iret)

! check if the dataset has been closed successfuly
!
      if (iret > 0) then

! print error about the problem with closing the dataset
!
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot close dataset: " // trim(name)

      end if

    else

! print error about the problem with creating the dataset
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create dataset: " // trim(name)

    end if

! release the space
!
    call h5sclose_f(sid, iret)

! check if the space has been released successfuly
!
    if (iret > 0) then

! print error about the problem with closing the space
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot close space for dataset: " // trim(name)

    end if

!-------------------------------------------------------------------------------
!
  end subroutine write_2d_array_integer_h5
!
!===============================================================================
!
! subroutine WRITE_3D_ARRAY_INTEGER_H5:
! ------------------------------------
!
!   Subroutine stores a three-dimensional integer array in a group specified by
!   identifier.
!
!   Arguments:
!
!     gid   - the HDF5 group identifier
!     name  - the string name describing the array
!     dm    - the array dimensions
!     value - the array values
!
!===============================================================================
!
  subroutine write_3d_array_integer_h5(gid, name, dm, var)

! import procedures and variables from other modules
!
    use hdf5           , only : H5T_NATIVE_INTEGER
    use hdf5           , only : hid_t, hsize_t
    use hdf5           , only : h5screate_simple_f, h5sclose_f
    use hdf5           , only : h5dcreate_f, h5dwrite_f, h5dclose_f
    use hdf5           , only : h5pset_chunk_f
    use iso_fortran_env, only : error_unit

! local variables are not implicit by default
!
    implicit none

! subroutine arguments
!
    integer(hid_t)                    , intent(in) :: gid
    character(len=*)                  , intent(in) :: name
    integer(hsize_t), dimension(3)    , intent(in) :: dm
    integer(kind=4) , dimension(:,:,:), intent(in) :: var

! HDF5 object identifiers
!
    integer(hid_t) :: sid, did

! procedure return value
!
    integer        :: iret

! local parameters
!
    character(len=*), parameter :: loc = 'IO::write_3d_array_integer_h5()'
!
!-------------------------------------------------------------------------------
!
! create a space for the array
!
    call h5screate_simple_f(3, dm(1:3), sid, iret)

! check if the space has been created successfuly, if not quit
!
    if (iret < 0) then

! print error about the problem with creating the space
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create space for dataset: " // trim(name)

! quit the subroutine
!
      return

    end if

! set the chunk size
!
    call h5pset_chunk_f(pid, 3, dm(1:3), iret)

! check if the chunk size has been set properly
!
    if (iret > 0) then

! print error about the problem with setting the chunk size
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot set the size of the chunk!"

    end if

! create the dataset
!
    call h5dcreate_f(gid, name, H5T_NATIVE_INTEGER, sid, did, iret, pid)

! check if the dataset has been created successfuly
!
    if (iret >= 0) then

! write the dataset data
!
      call h5dwrite_f(did, H5T_NATIVE_INTEGER, var(:,:,:), dm(1:3), iret, sid)

! check if the dataset has been written successfuly
!
      if (iret > 0) then

! print error about the problem with writing down the dataset
!
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot write dataset: " // trim(name)

      end if

! close the dataset
!
      call h5dclose_f(did, iret)

! check if the dataset has been closed successfuly
!
      if (iret > 0) then

! print error about the problem with closing the dataset
!
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot close dataset: " // trim(name)

      end if

    else

! print error about the problem with creating the dataset
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create dataset: " // trim(name)

    end if

! release the space
!
    call h5sclose_f(sid, iret)

! check if the space has been released successfuly
!
    if (iret > 0) then

! print error about the problem with closing the space
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot close space for dataset: " // trim(name)

    end if

!-------------------------------------------------------------------------------
!
  end subroutine write_3d_array_integer_h5
!
!===============================================================================
!
! subroutine WRITE_4D_ARRAY_INTEGER_H5:
! ------------------------------------
!
!   Subroutine stores a four-dimensional integer array in a group specified by
!   identifier.
!
!   Arguments:
!
!     gid   - the HDF5 group identifier
!     name  - the string name describing the array
!     dm    - the array dimensions
!     value - the array values
!
!===============================================================================
!
  subroutine write_4d_array_integer_h5(gid, name, dm, var)

! import procedures and variables from other modules
!
    use hdf5           , only : H5T_NATIVE_INTEGER
    use hdf5           , only : hid_t, hsize_t
    use hdf5           , only : h5screate_simple_f, h5sclose_f
    use hdf5           , only : h5dcreate_f, h5dwrite_f, h5dclose_f
    use hdf5           , only : h5pset_chunk_f
    use iso_fortran_env, only : error_unit

! local variables are not implicit by default
!
    implicit none

! subroutine arguments
!
    integer(hid_t)                      , intent(in) :: gid
    character(len=*)                    , intent(in) :: name
    integer(hsize_t), dimension(4)      , intent(in) :: dm
    integer(kind=4) , dimension(:,:,:,:), intent(in) :: var

! HDF5 object identifiers
!
    integer(hid_t) :: sid, did

! procedure return value
!
    integer        :: iret

! local parameters
!
    character(len=*), parameter :: loc = 'IO::write_4d_array_integer_h5()'
!
!-------------------------------------------------------------------------------
!
! create a space for the array
!
    call h5screate_simple_f(4, dm(1:4), sid, iret)

! check if the space has been created successfuly, if not quit
!
    if (iret < 0) then

! print error about the problem with creating the space
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create space for dataset: " // trim(name)

! quit the subroutine
!
      return

    end if

! set the chunk size
!
    call h5pset_chunk_f(pid, 4, dm(1:4), iret)

! check if the chunk size has been set properly
!
    if (iret > 0) then

! print error about the problem with setting the chunk size
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot set the size of the chunk!"

    end if

! create the dataset
!
    call h5dcreate_f(gid, name, H5T_NATIVE_INTEGER, sid, did, iret, pid)

! check if the dataset has been created successfuly
!
    if (iret >= 0) then

! write the dataset data
!
      call h5dwrite_f(did, H5T_NATIVE_INTEGER, var(:,:,:,:), dm(1:4), iret, sid)

! check if the dataset has been written successfuly
!
      if (iret > 0) then

! print error about the problem with writing down the dataset
!
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot write dataset: " // trim(name)

      end if

! close the dataset
!
      call h5dclose_f(did, iret)

! check if the dataset has been closed successfuly
!
      if (iret > 0) then

! print error about the problem with closing the dataset
!
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot close dataset: " // trim(name)

      end if

    else

! print error about the problem with creating the dataset
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create dataset: " // trim(name)

    end if

! release the space
!
    call h5sclose_f(sid, iret)

! check if the space has been released successfuly
!
    if (iret > 0) then

! print error about the problem with closing the space
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot close space for dataset: " // trim(name)

    end if

!-------------------------------------------------------------------------------
!
  end subroutine write_4d_array_integer_h5
!
!===============================================================================
!
! subroutine WRITE_5D_ARRAY_INTEGER_H5:
! ------------------------------------
!
!   Subroutine stores a five-dimensional integer array in a group specified by
!   identifier.
!
!   Arguments:
!
!     gid   - the HDF5 group identifier
!     name  - the string name describing the array
!     dm    - the array dimensions
!     value - the array values
!
!===============================================================================
!
  subroutine write_5d_array_integer_h5(gid, name, dm, var)

! import procedures and variables from other modules
!
    use hdf5           , only : H5T_NATIVE_INTEGER
    use hdf5           , only : hid_t, hsize_t
    use hdf5           , only : h5screate_simple_f, h5sclose_f
    use hdf5           , only : h5dcreate_f, h5dwrite_f, h5dclose_f
    use hdf5           , only : h5pset_chunk_f
    use iso_fortran_env, only : error_unit

! local variables are not implicit by default
!
    implicit none

! subroutine arguments
!
    integer(hid_t)                        , intent(in) :: gid
    character(len=*)                      , intent(in) :: name
    integer(hsize_t), dimension(5)        , intent(in) :: dm
    integer(kind=4) , dimension(:,:,:,:,:), intent(in) :: var

! HDF5 object identifiers
!
    integer(hid_t) :: sid, did

! procedure return value
!
    integer        :: iret

! local parameters
!
    character(len=*), parameter :: loc = 'IO::write_5d_array_integer_h5()'
!
!-------------------------------------------------------------------------------
!
! create a space for the array
!
    call h5screate_simple_f(5, dm(1:5), sid, iret)

! check if the space has been created successfuly, if not quit
!
    if (iret < 0) then

! print error about the problem with creating the space
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create space for dataset: " // trim(name)

! quit the subroutine
!
      return

    end if

! set the chunk size
!
    call h5pset_chunk_f(pid, 5, dm(1:5), iret)

! check if the chunk size has been set properly
!
    if (iret > 0) then

! print error about the problem with setting the chunk size
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot set the size of the chunk!"

    end if

! create the dataset
!
    call h5dcreate_f(gid, name, H5T_NATIVE_INTEGER, sid, did, iret, pid)

! check if the dataset has been created successfuly
!
    if (iret >= 0) then

! write the dataset data
!
      call h5dwrite_f(did, H5T_NATIVE_INTEGER, var(:,:,:,:,:), dm(1:5)         &
                                                                  , iret, sid)

! check if the dataset has been written successfuly
!
      if (iret > 0) then

! print error about the problem with writing down the dataset
!
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot write dataset: " // trim(name)

      end if

! close the dataset
!
      call h5dclose_f(did, iret)

! check if the dataset has been closed successfuly
!
      if (iret > 0) then

! print error about the problem with closing the dataset
!
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot close dataset: " // trim(name)

      end if

    else

! print error about the problem with creating the dataset
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create dataset: " // trim(name)

    end if

! release the space
!
    call h5sclose_f(sid, iret)

! check if the space has been released successfuly
!
    if (iret > 0) then

! print error about the problem with closing the space
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot close space for dataset: " // trim(name)

    end if

!-------------------------------------------------------------------------------
!
  end subroutine write_5d_array_integer_h5
!
!===============================================================================
!
! subroutine WRITE_1D_ARRAY_DOUBLE_H5:
! -----------------------------------
!
!   Subroutine stores a one-dimensional double precision array in a group
!   specified by identifier.
!
!   Arguments:
!
!     gid   - the HDF5 group identifier
!     name  - the string name describing the array
!     dm    - the array dimensions
!     value - the array values
!
!===============================================================================
!
  subroutine write_1d_array_double_h5(gid, name, ln, var)

! import procedures and variables from other modules
!
    use hdf5           , only : H5T_NATIVE_DOUBLE
    use hdf5           , only : hid_t, hsize_t
    use hdf5           , only : h5screate_simple_f, h5sclose_f
    use hdf5           , only : h5dcreate_f, h5dwrite_f, h5dclose_f
    use hdf5           , only : h5pset_chunk_f
    use iso_fortran_env, only : error_unit

! local variables are not implicit by default
!
    implicit none

! subroutine arguments
!
    integer(hid_t)                , intent(in) :: gid
    character(len=*)              , intent(in) :: name
    integer(hsize_t)              , intent(in) :: ln
    real(kind=8)    , dimension(:), intent(in) :: var

! HDF5 object identifiers
!
    integer(hid_t) :: sid, did

! array dimensions
!
    integer(hsize_t), dimension(1) :: dm

! procedure return value
!
    integer        :: iret

! local parameters
!
    character(len=*), parameter :: loc = 'IO::write_1d_array_double_h5()'
!
!-------------------------------------------------------------------------------
!
! substitute array dimensions
!
    dm(1) = ln

! create a space for the array
!
    call h5screate_simple_f(1, dm(1:1), sid, iret)

! check if the space has been created successfuly, if not quit
!
    if (iret < 0) then

! print error about the problem with creating the space
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create space for dataset: " // trim(name)

! quit the subroutine
!
      return

    end if

! set the chunk size
!
    call h5pset_chunk_f(pid, 1, dm(1:1), iret)

! check if the chunk size has been set properly
!
    if (iret > 0) then

! print error about the problem with setting the chunk size
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot set the size of the chunk!"

    end if

! create the dataset
!
    call h5dcreate_f(gid, name, H5T_NATIVE_DOUBLE, sid, did, iret, pid)

! check if the dataset has been created successfuly
!
    if (iret >= 0) then

! write the dataset data
!
      call h5dwrite_f(did, H5T_NATIVE_DOUBLE, var(:), dm(1:1), iret, sid)

! check if the dataset has been written successfuly
!
      if (iret > 0) then

! print error about the problem with writing down the dataset
!
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot write dataset: " // trim(name)

      end if

! close the dataset
!
      call h5dclose_f(did, iret)

! check if the dataset has been closed successfuly
!
      if (iret > 0) then

! print error about the problem with closing the dataset
!
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot close dataset: " // trim(name)

      end if

    else

! print error about the problem with creating the dataset
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create dataset: " // trim(name)

    end if

! release the space
!
    call h5sclose_f(sid, iret)

! check if the space has been released successfuly
!
    if (iret > 0) then

! print error about the problem with closing the space
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot close space for dataset: " // trim(name)

    end if

!-------------------------------------------------------------------------------
!
  end subroutine write_1d_array_double_h5
!
!===============================================================================
!
! subroutine WRITE_2D_ARRAY_DOUBLE_H5:
! ------------------------------------
!
!   Subroutine stores a two-dimensional double precision array in a group
!   specified by identifier.
!
!   Arguments:
!
!     gid   - the HDF5 group identifier
!     name  - the string name describing the array
!     dm    - the array dimensions
!     value - the array values
!
!===============================================================================
!
  subroutine write_2d_array_double_h5(gid, name, dm, var)

! import procedures and variables from other modules
!
    use hdf5           , only : H5T_NATIVE_DOUBLE
    use hdf5           , only : hid_t, hsize_t
    use hdf5           , only : h5screate_simple_f, h5sclose_f
    use hdf5           , only : h5dcreate_f, h5dwrite_f, h5dclose_f
    use hdf5           , only : h5pset_chunk_f
    use iso_fortran_env, only : error_unit

! local variables are not implicit by default
!
    implicit none

! subroutine arguments
!
    integer(hid_t)                  , intent(in) :: gid
    character(len=*)                , intent(in) :: name
    integer(hsize_t), dimension(2)  , intent(in) :: dm
    real(kind=8)    , dimension(:,:), intent(in) :: var

! HDF5 object identifiers
!
    integer(hid_t) :: sid, did

! procedure return value
!
    integer        :: iret

! local parameters
!
    character(len=*), parameter :: loc = 'IO::write_2d_array_double_h5()'
!
!-------------------------------------------------------------------------------
!
! create a space for the array
!
    call h5screate_simple_f(2, dm(1:2), sid, iret)

! check if the space has been created successfuly, if not quit
!
    if (iret < 0) then

! print error about the problem with creating the space
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create space for dataset: " // trim(name)

! quit the subroutine
!
      return

    end if

! set the chunk size
!
    call h5pset_chunk_f(pid, 2, dm(1:2), iret)

! check if the chunk size has been set properly
!
    if (iret > 0) then

! print error about the problem with setting the chunk size
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot set the size of the chunk!"

    end if

! create the dataset
!
    call h5dcreate_f(gid, name, H5T_NATIVE_DOUBLE, sid, did, iret, pid)

! check if the dataset has been created successfuly
!
    if (iret >= 0) then

! write the dataset data
!
      call h5dwrite_f(did, H5T_NATIVE_DOUBLE, var(:,:), dm(1:2), iret, sid)

! check if the dataset has been written successfuly
!
      if (iret > 0) then

! print error about the problem with writing down the dataset
!
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot write dataset: " // trim(name)

      end if

! close the dataset
!
      call h5dclose_f(did, iret)

! check if the dataset has been closed successfuly
!
      if (iret > 0) then

! print error about the problem with closing the dataset
!
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot close dataset: " // trim(name)

      end if

    else

! print error about the problem with creating the dataset
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create dataset: " // trim(name)

    end if

! release the space
!
    call h5sclose_f(sid, iret)

! check if the space has been released successfuly
!
    if (iret > 0) then

! print error about the problem with closing the space
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot close space for dataset: " // trim(name)

    end if

!-------------------------------------------------------------------------------
!
  end subroutine write_2d_array_double_h5
!
!===============================================================================
!
! subroutine WRITE_3D_ARRAY_DOUBLE_H5:
! -----------------------------------
!
!   Subroutine stores a three-dimensional double precision array in a group
!   specified by identifier.
!
!   Arguments:
!
!     gid   - the HDF5 group identifier
!     name  - the string name describing the array
!     dm    - the array dimensions
!     value - the array values
!
!===============================================================================
!
  subroutine write_3d_array_double_h5(gid, name, dm, var)

! import procedures and variables from other modules
!
    use hdf5           , only : H5T_NATIVE_DOUBLE
    use hdf5           , only : hid_t, hsize_t
    use hdf5           , only : h5screate_simple_f, h5sclose_f
    use hdf5           , only : h5dcreate_f, h5dwrite_f, h5dclose_f
    use hdf5           , only : h5pset_chunk_f
    use iso_fortran_env, only : error_unit

! local variables are not implicit by default
!
    implicit none

! subroutine arguments
!
    integer(hid_t)                    , intent(in) :: gid
    character(len=*)                  , intent(in) :: name
    integer(hsize_t), dimension(3)    , intent(in) :: dm
    real(kind=8)    , dimension(:,:,:), intent(in) :: var

! HDF5 object identifiers
!
    integer(hid_t) :: sid, did

! procedure return value
!
    integer        :: iret

! local parameters
!
    character(len=*), parameter :: loc = 'IO::write_3d_array_double_h5()'
!
!-------------------------------------------------------------------------------
!
! create a space for the array
!
    call h5screate_simple_f(3, dm(1:3), sid, iret)

! check if the space has been created successfuly, if not quit
!
    if (iret < 0) then

! print error about the problem with creating the space
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create space for dataset: " // trim(name)

! quit the subroutine
!
      return

    end if

! set the chunk size
!
    call h5pset_chunk_f(pid, 3, dm(1:3), iret)

! check if the chunk size has been set properly
!
    if (iret > 0) then

! print error about the problem with setting the chunk size
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot set the size of the chunk!"

    end if

! create the dataset
!
    call h5dcreate_f(gid, name, H5T_NATIVE_DOUBLE, sid, did, iret, pid)

! check if the dataset has been created successfuly
!
    if (iret >= 0) then

! write the dataset data
!
      call h5dwrite_f(did, H5T_NATIVE_DOUBLE, var(:,:,:), dm(1:3), iret, sid)

! check if the dataset has been written successfuly
!
      if (iret > 0) then

! print error about the problem with writing down the dataset
!
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot write dataset: " // trim(name)

      end if

! close the dataset
!
      call h5dclose_f(did, iret)

! check if the dataset has been closed successfuly
!
      if (iret > 0) then

! print error about the problem with closing the dataset
!
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot close dataset: " // trim(name)

      end if

    else

! print error about the problem with creating the dataset
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create dataset: " // trim(name)

    end if

! release the space
!
    call h5sclose_f(sid, iret)

! check if the space has been released successfuly
!
    if (iret > 0) then

! print error about the problem with closing the space
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot close space for dataset: " // trim(name)

    end if

!-------------------------------------------------------------------------------
!
  end subroutine write_3d_array_double_h5
!
!===============================================================================
!
! subroutine WRITE_4D_ARRAY_DOUBLE_H5:
! ------------------------------------
!
!   Subroutine stores a four-dimensional double precision array in a group
!   specified by identifier.
!
!   Arguments:
!
!     gid   - the HDF5 group identifier
!     name  - the string name describing the array
!     dm    - the array dimensions
!     value - the array values
!
!===============================================================================
!
  subroutine write_4d_array_double_h5(gid, name, dm, var)

! import procedures and variables from other modules
!
    use hdf5           , only : H5T_NATIVE_DOUBLE
    use hdf5           , only : hid_t, hsize_t
    use hdf5           , only : h5screate_simple_f, h5sclose_f
    use hdf5           , only : h5dcreate_f, h5dwrite_f, h5dclose_f
    use hdf5           , only : h5pset_chunk_f
    use iso_fortran_env, only : error_unit

! local variables are not implicit by default
!
    implicit none

! subroutine arguments
!
    integer(hid_t)                      , intent(in) :: gid
    character(len=*)                    , intent(in) :: name
    integer(hsize_t), dimension(4)      , intent(in) :: dm
    real(kind=8)    , dimension(:,:,:,:), intent(in) :: var

! HDF5 object identifiers
!
    integer(hid_t) :: sid, did

! procedure return value
!
    integer        :: iret

! local parameters
!
    character(len=*), parameter :: loc = 'IO::write_4d_array_double_h5()'
!
!-------------------------------------------------------------------------------
!
! create a space for the array
!
    call h5screate_simple_f(4, dm(1:4), sid, iret)

! check if the space has been created successfuly, if not quit
!
    if (iret < 0) then

! print error about the problem with creating the space
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create space for dataset: " // trim(name)

! quit the subroutine
!
      return

    end if

! set the chunk size
!
    call h5pset_chunk_f(pid, 4, dm(1:4), iret)

! check if the chunk size has been set properly
!
    if (iret > 0) then

! print error about the problem with setting the chunk size
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot set the size of the chunk!"

    end if

! create the dataset
!
    call h5dcreate_f(gid, name, H5T_NATIVE_DOUBLE, sid, did, iret, pid)

! check if the dataset has been created successfuly
!
    if (iret >= 0) then

! write the dataset data
!
      call h5dwrite_f(did, H5T_NATIVE_DOUBLE, var(:,:,:,:), dm(1:4), iret, sid)

! check if the dataset has been written successfuly
!
      if (iret > 0) then

! print error about the problem with writing down the dataset
!
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot write dataset: " // trim(name)

      end if

! close the dataset
!
      call h5dclose_f(did, iret)

! check if the dataset has been closed successfuly
!
      if (iret > 0) then

! print error about the problem with closing the dataset
!
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot close dataset: " // trim(name)

      end if

    else

! print error about the problem with creating the dataset
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create dataset: " // trim(name)

    end if

! release the space
!
    call h5sclose_f(sid, iret)

! check if the space has been released successfuly
!
    if (iret > 0) then

! print error about the problem with closing the space
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot close space for dataset: " // trim(name)

    end if

!-------------------------------------------------------------------------------
!
  end subroutine write_4d_array_double_h5
!
!===============================================================================
!
! subroutine WRITE_5D_ARRAY_DOUBLE_H5:
! -----------------------------------
!
!   Subroutine stores a five-dimensional double precision array in a group
!   specified by identifier.
!
!   Arguments:
!
!     gid   - the HDF5 group identifier
!     name  - the string name describing the array
!     dm    - the array dimensions
!     value - the array values
!
!===============================================================================
!
  subroutine write_5d_array_double_h5(gid, name, dm, var)

! import procedures and variables from other modules
!
    use hdf5           , only : H5T_NATIVE_DOUBLE
    use hdf5           , only : hid_t, hsize_t
    use hdf5           , only : h5screate_simple_f, h5sclose_f
    use hdf5           , only : h5dcreate_f, h5dwrite_f, h5dclose_f
    use hdf5           , only : h5pset_chunk_f
    use iso_fortran_env, only : error_unit

! local variables are not implicit by default
!
    implicit none

! subroutine arguments
!
    integer(hid_t)                        , intent(in) :: gid
    character(len=*)                      , intent(in) :: name
    integer(hsize_t), dimension(5)        , intent(in) :: dm
    real(kind=8)    , dimension(:,:,:,:,:), intent(in) :: var

! HDF5 object identifiers
!
    integer(hid_t) :: sid, did

! procedure return value
!
    integer        :: iret

! local parameters
!
    character(len=*), parameter :: loc = 'IO::write_5d_array_double_h5()'
!
!-------------------------------------------------------------------------------
!
! create a space for the array
!
    call h5screate_simple_f(5, dm(1:5), sid, iret)

! check if the space has been created successfuly, if not quit
!
    if (iret < 0) then

! print error about the problem with creating the space
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create space for dataset: " // trim(name)

! quit the subroutine
!
      return

    end if

! set the chunk size
!
    call h5pset_chunk_f(pid, 5, dm(1:5), iret)

! check if the chunk size has been set properly
!
    if (iret > 0) then

! print error about the problem with setting the chunk size
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot set the size of the chunk!"

    end if

! create the dataset
!
    call h5dcreate_f(gid, name, H5T_NATIVE_DOUBLE, sid, did, iret, pid)

! check if the dataset has been created successfuly
!
    if (iret >= 0) then

! write the dataset data
!
      call h5dwrite_f(did, H5T_NATIVE_DOUBLE, var(:,:,:,:,:), dm(1:5)          &
                                                                  , iret, sid)

! check if the dataset has been written successfuly
!
      if (iret > 0) then

! print error about the problem with writing down the dataset
!
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot write dataset: " // trim(name)

      end if

! close the dataset
!
      call h5dclose_f(did, iret)

! check if the dataset has been closed successfuly
!
      if (iret > 0) then

! print error about the problem with closing the dataset
!
        write(error_unit,"('[',a,']: ',a)") trim(loc)                          &
                        , "Cannot close dataset: " // trim(name)

      end if

    else

! print error about the problem with creating the dataset
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot create dataset: " // trim(name)

    end if

! release the space
!
    call h5sclose_f(sid, iret)

! check if the space has been released successfuly
!
    if (iret > 0) then

! print error about the problem with closing the space
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot close space for dataset: " // trim(name)

    end if

!-------------------------------------------------------------------------------
!
  end subroutine write_5d_array_double_h5

!===============================================================================
!
! READ_ARRAY SUBROUTINES
!
!===============================================================================
!
! subroutine READ_1D_ARRAY_INTEGER_H5:
! -----------------------------------
!
!   Subroutine restores a one-dimensional integer array from a group specified
!   by identifier.
!
!   Arguments:
!
!     gid   - the HDF5 group identifier
!     name  - the string name describing the array
!     dm    - the array dimensions
!     value - the array values
!
!===============================================================================
!
  subroutine read_1d_array_integer_h5(gid, name, dm, var)

! import procedures and variables from other modules
!
    use hdf5           , only : H5T_NATIVE_INTEGER
    use hdf5           , only : hid_t, hsize_t
    use hdf5           , only : h5dopen_f, h5dread_f, h5dclose_f
    use iso_fortran_env, only : error_unit

! local variables are not implicit by default
!
    implicit none

! subroutine arguments
!
    integer(hid_t)                , intent(in)    :: gid
    character(len=*)              , intent(in)    :: name
    integer(hsize_t), dimension(1), intent(inout) :: dm
    integer(kind=4) , dimension(:), intent(inout) :: var

! local variables
!
    integer(hid_t) :: did
    integer        :: iret

! local parameters
!
    character(len=*), parameter :: loc = 'IO::read_1d_array_integer_h5()'
!
!-------------------------------------------------------------------------------
!
! open the dataset
!
    call h5dopen_f(gid, name, did, iret)

! check if the dataset has been opened successfuly
!
    if (iret < 0) then

! print error about the problem with opening the data space
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot open dataset: " // trim(name)

! quit the subroutine
!
      return

    end if

! read dataset data
!
    call h5dread_f(did, H5T_NATIVE_INTEGER, var(:), dm(1:1), iret)

! check if the dataset has been read successfuly
!
    if (iret > 0) then

! print error about the problem with reading the dataset
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot read dataset: " // trim(name)

    end if

! close the dataset
!
    call h5dclose_f(did, iret)

! check if the dataset has been closed successfuly
!
    if (iret > 0) then

! print error about the problem with closing the dataset
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot close dataset: " // trim(name)

    end if

!-------------------------------------------------------------------------------
!
  end subroutine read_1d_array_integer_h5
!
!===============================================================================
!
! subroutine READ_2D_ARRAY_INTEGER_H5:
! -----------------------------------
!
!   Subroutine restores a two-dimensional integer array from a group specified
!   by identifier.
!
!   Arguments:
!
!     gid   - the HDF5 group identifier
!     name  - the string name describing the array
!     dm    - the array dimensions
!     value - the array values
!
!===============================================================================
!
  subroutine read_2d_array_integer_h5(gid, name, dm, var)

! import procedures and variables from other modules
!
    use hdf5           , only : H5T_NATIVE_INTEGER
    use hdf5           , only : hid_t, hsize_t
    use hdf5           , only : h5dopen_f, h5dread_f, h5dclose_f
    use iso_fortran_env, only : error_unit

! local variables are not implicit by default
!
    implicit none

! subroutine arguments
!
    integer(hid_t)                  , intent(in)    :: gid
    character(len=*)                , intent(in)    :: name
    integer(hsize_t), dimension(2)  , intent(inout) :: dm
    integer(kind=4) , dimension(:,:), intent(inout) :: var

! local variables
!
    integer(hid_t) :: did
    integer        :: iret

! local parameters
!
    character(len=*), parameter :: loc = 'IO::read_2d_array_integer_h5()'
!
!-------------------------------------------------------------------------------
!
! open the dataset
!
    call h5dopen_f(gid, name, did, iret)

! check if the dataset has been opened successfuly
!
    if (iret < 0) then

! print error about the problem with opening the data space
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot open dataset: " // trim(name)

! quit the subroutine
!
      return

    end if

! read dataset data
!
    call h5dread_f(did, H5T_NATIVE_INTEGER, var(:,:), dm(1:2), iret)

! check if the dataset has been read successfuly
!
    if (iret > 0) then

! print error about the problem with reading the dataset
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot read dataset: " // trim(name)

    end if

! close the dataset
!
    call h5dclose_f(did, iret)

! check if the dataset has been closed successfuly
!
    if (iret > 0) then

! print error about the problem with closing the dataset
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot close dataset: " // trim(name)

    end if

!-------------------------------------------------------------------------------
!
  end subroutine read_2d_array_integer_h5
!
!===============================================================================
!
! subroutine READ_3D_ARRAY_INTEGER_H5:
! -----------------------------------
!
!   Subroutine restores a three-dimensional integer array from a group specified
!   by identifier.
!
!   Arguments:
!
!     gid   - the HDF5 group identifier
!     name  - the string name describing the array
!     dm    - the array dimensions
!     value - the array values
!
!===============================================================================
!
  subroutine read_3d_array_integer_h5(gid, name, dm, var)

! import procedures and variables from other modules
!
    use hdf5           , only : H5T_NATIVE_INTEGER
    use hdf5           , only : hid_t, hsize_t
    use hdf5           , only : h5dopen_f, h5dread_f, h5dclose_f
    use iso_fortran_env, only : error_unit

! local variables are not implicit by default
!
    implicit none

! subroutine arguments
!
    integer(hid_t)                    , intent(in)    :: gid
    character(len=*)                  , intent(in)    :: name
    integer(hsize_t), dimension(3)    , intent(inout) :: dm
    integer(kind=4) , dimension(:,:,:), intent(inout) :: var

! local variables
!
    integer(hid_t) :: did
    integer        :: iret

! local parameters
!
    character(len=*), parameter :: loc = 'IO::read_3d_array_integer_h5()'
!
!-------------------------------------------------------------------------------
!
! open the dataset
!
    call h5dopen_f(gid, name, did, iret)

! check if the dataset has been opened successfuly
!
    if (iret < 0) then

! print error about the problem with opening the data space
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot open dataset: " // trim(name)

! quit the subroutine
!
      return

    end if

! read dataset data
!
    call h5dread_f(did, H5T_NATIVE_INTEGER, var(:,:,:), dm(1:3), iret)

! check if the dataset has been read successfuly
!
    if (iret > 0) then

! print error about the problem with reading the dataset
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot read dataset: " // trim(name)

    end if

! close the dataset
!
    call h5dclose_f(did, iret)

! check if the dataset has been closed successfuly
!
    if (iret > 0) then

! print error about the problem with closing the dataset
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot close dataset: " // trim(name)

    end if

!-------------------------------------------------------------------------------
!
  end subroutine read_3d_array_integer_h5
!
!===============================================================================
!
! subroutine READ_4D_ARRAY_INTEGER_H5:
! -----------------------------------
!
!   Subroutine restores a four-dimensional integer array from a group specified
!   by identifier.
!
!   Arguments:
!
!     gid   - the HDF5 group identifier
!     name  - the string name describing the array
!     dm    - the array dimensions
!     value - the array values
!
!===============================================================================
!
  subroutine read_4d_array_integer_h5(gid, name, dm, var)

! import procedures and variables from other modules
!
    use hdf5           , only : H5T_NATIVE_INTEGER
    use hdf5           , only : hid_t, hsize_t
    use hdf5           , only : h5dopen_f, h5dread_f, h5dclose_f
    use iso_fortran_env, only : error_unit

! local variables are not implicit by default
!
    implicit none

! subroutine arguments
!
    integer(hid_t)                      , intent(in)    :: gid
    character(len=*)                    , intent(in)    :: name
    integer(hsize_t), dimension(4)      , intent(inout) :: dm
    integer(kind=4) , dimension(:,:,:,:), intent(inout) :: var

! local variables
!
    integer(hid_t) :: did
    integer        :: iret

! local parameters
!
    character(len=*), parameter :: loc = 'IO::read_4d_array_integer_h5()'
!
!-------------------------------------------------------------------------------
!
! open the dataset
!
    call h5dopen_f(gid, name, did, iret)

! check if the dataset has been opened successfuly
!
    if (iret < 0) then

! print error about the problem with opening the data space
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot open dataset: " // trim(name)

! quit the subroutine
!
      return

    end if

! read dataset data
!
    call h5dread_f(did, H5T_NATIVE_INTEGER, var(:,:,:,:), dm(1:4), iret)

! check if the dataset has been read successfuly
!
    if (iret > 0) then

! print error about the problem with reading the dataset
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot read dataset: " // trim(name)

    end if

! close the dataset
!
    call h5dclose_f(did, iret)

! check if the dataset has been closed successfuly
!
    if (iret > 0) then

! print error about the problem with closing the dataset
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot close dataset: " // trim(name)

    end if

!-------------------------------------------------------------------------------
!
  end subroutine read_4d_array_integer_h5
!
!===============================================================================
!
! subroutine READ_5D_ARRAY_INTEGER_H5:
! -----------------------------------
!
!   Subroutine restores a five-dimensional integer array from a group specified
!   by identifier.
!
!   Arguments:
!
!     gid   - the HDF5 group identifier
!     name  - the string name describing the array
!     dm    - the array dimensions
!     value - the array values
!
!===============================================================================
!
  subroutine read_5d_array_integer_h5(gid, name, dm, var)

! import procedures and variables from other modules
!
    use hdf5           , only : H5T_NATIVE_INTEGER
    use hdf5           , only : hid_t, hsize_t
    use hdf5           , only : h5dopen_f, h5dread_f, h5dclose_f
    use iso_fortran_env, only : error_unit

! local variables are not implicit by default
!
    implicit none

! subroutine arguments
!
    integer(hid_t)                        , intent(in)    :: gid
    character(len=*)                      , intent(in)    :: name
    integer(hsize_t), dimension(5)        , intent(inout) :: dm
    integer(kind=4) , dimension(:,:,:,:,:), intent(inout) :: var

! local variables
!
    integer(hid_t) :: did
    integer        :: iret

! local parameters
!
    character(len=*), parameter :: loc = 'IO::read_5d_array_integer_h5()'
!
!-------------------------------------------------------------------------------
!
! open the dataset
!
    call h5dopen_f(gid, name, did, iret)

! check if the dataset has been opened successfuly
!
    if (iret < 0) then

! print error about the problem with opening the data space
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot open dataset: " // trim(name)

! quit the subroutine
!
      return

    end if

! read dataset data
!
    call h5dread_f(did, H5T_NATIVE_INTEGER, var(:,:,:,:,:), dm(1:5), iret)

! check if the dataset has been read successfuly
!
    if (iret > 0) then

! print error about the problem with reading the dataset
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot read dataset: " // trim(name)

    end if

! close the dataset
!
    call h5dclose_f(did, iret)

! check if the dataset has been closed successfuly
!
    if (iret > 0) then

! print error about the problem with closing the dataset
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot close dataset: " // trim(name)

    end if

!-------------------------------------------------------------------------------
!
  end subroutine read_5d_array_integer_h5
!
!===============================================================================
!
! subroutine READ_1D_ARRAY_DOUBLE_H5:
! ----------------------------------
!
!   Subroutine restores a one-dimensional double precision array from a group
!   specified by identifier.
!
!   Arguments:
!
!     gid   - the HDF5 group identifier
!     name  - the string name describing the array
!     dm    - the array dimensions
!     value - the array values
!
!===============================================================================
!
  subroutine read_1d_array_double_h5(gid, name, dm, var)

! import procedures and variables from other modules
!
    use hdf5           , only : H5T_NATIVE_DOUBLE
    use hdf5           , only : hid_t, hsize_t
    use hdf5           , only : h5dopen_f, h5dread_f, h5dclose_f
    use iso_fortran_env, only : error_unit

! local variables are not implicit by default
!
    implicit none

! subroutine arguments
!
    integer(hid_t)                , intent(in)    :: gid
    character(len=*)              , intent(in)    :: name
    integer(hsize_t), dimension(1), intent(inout) :: dm
    real(kind=8)    , dimension(:), intent(inout) :: var

! local variables
!
    integer(hid_t) :: did
    integer        :: iret

! local parameters
!
    character(len=*), parameter :: loc = 'IO::read_1d_array_double_h5()'
!
!-------------------------------------------------------------------------------
!
! open the dataset
!
    call h5dopen_f(gid, name, did, iret)

! check if the dataset has been opened successfuly
!
    if (iret < 0) then

! print error about the problem with opening the data space
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot open dataset: " // trim(name)

! quit the subroutine
!
      return

    end if

! read dataset data
!
    call h5dread_f(did, H5T_NATIVE_DOUBLE, var(:), dm(1:1), iret)

! check if the dataset has been read successfuly
!
    if (iret > 0) then

! print error about the problem with reading the dataset
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot read dataset: " // trim(name)

    end if

! close the dataset
!
    call h5dclose_f(did, iret)

! check if the dataset has been closed successfuly
!
    if (iret > 0) then

! print error about the problem with closing the dataset
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot close dataset: " // trim(name)

    end if

!-------------------------------------------------------------------------------
!
  end subroutine read_1d_array_double_h5
!
!===============================================================================
!
! subroutine READ_2D_ARRAY_DOUBLE_H5:
! -----------------------------------
!
!   Subroutine restores a two-dimensional double precision array from a group
!   specified by identifier.
!
!   Arguments:
!
!     gid   - the HDF5 group identifier
!     name  - the string name describing the array
!     dm    - the array dimensions
!     value - the array values
!
!===============================================================================
!
  subroutine read_2d_array_double_h5(gid, name, dm, var)

! import procedures and variables from other modules
!
    use hdf5           , only : H5T_NATIVE_DOUBLE
    use hdf5           , only : hid_t, hsize_t
    use hdf5           , only : h5dopen_f, h5dread_f, h5dclose_f
    use iso_fortran_env, only : error_unit

! local variables are not implicit by default
!
    implicit none

! subroutine arguments
!
    integer(hid_t)                  , intent(in)    :: gid
    character(len=*)                , intent(in)    :: name
    integer(hsize_t), dimension(2)  , intent(inout) :: dm
    real(kind=8)    , dimension(:,:), intent(inout) :: var

! local variables
!
    integer(hid_t) :: did
    integer        :: iret

! local parameters
!
    character(len=*), parameter :: loc = 'IO::read_2d_array_double_h5()'
!
!-------------------------------------------------------------------------------
!
! open the dataset
!
    call h5dopen_f(gid, name, did, iret)

! check if the dataset has been opened successfuly
!
    if (iret < 0) then

! print error about the problem with opening the data space
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot open dataset: " // trim(name)

! quit the subroutine
!
      return

    end if

! read dataset data
!
    call h5dread_f(did, H5T_NATIVE_DOUBLE, var(:,:), dm(1:2), iret)

! check if the dataset has been read successfuly
!
    if (iret > 0) then

! print error about the problem with reading the dataset
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot read dataset: " // trim(name)

    end if

! close the dataset
!
    call h5dclose_f(did, iret)

! check if the dataset has been closed successfuly
!
    if (iret > 0) then

! print error about the problem with closing the dataset
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot close dataset: " // trim(name)

    end if

!-------------------------------------------------------------------------------
!
  end subroutine read_2d_array_double_h5
!
!===============================================================================
!
! subroutine READ_3D_ARRAY_DOUBLE_H5:
! -----------------------------------
!
!   Subroutine restores a three-dimensional double precision array from a group
!   specified by identifier.
!
!   Arguments:
!
!     gid   - the HDF5 group identifier
!     name  - the string name describing the array
!     dm    - the array dimensions
!     value - the array values
!
!===============================================================================
!
  subroutine read_3d_array_double_h5(gid, name, dm, var)

! import procedures and variables from other modules
!
    use hdf5           , only : H5T_NATIVE_DOUBLE
    use hdf5           , only : hid_t, hsize_t
    use hdf5           , only : h5dopen_f, h5dread_f, h5dclose_f
    use iso_fortran_env, only : error_unit

! local variables are not implicit by default
!
    implicit none

! subroutine arguments
!
    integer(hid_t)                    , intent(in)    :: gid
    character(len=*)                  , intent(in)    :: name
    integer(hsize_t), dimension(3)    , intent(inout) :: dm
    real(kind=8)    , dimension(:,:,:), intent(inout) :: var

! local variables
!
    integer(hid_t) :: did
    integer        :: iret

! local parameters
!
    character(len=*), parameter :: loc = 'IO::read_3d_array_double_h5()'
!
!-------------------------------------------------------------------------------
!
! open the dataset
!
    call h5dopen_f(gid, name, did, iret)

! check if the dataset has been opened successfuly
!
    if (iret < 0) then

! print error about the problem with opening the data space
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot open dataset: " // trim(name)

! quit the subroutine
!
      return

    end if

! read dataset data
!
    call h5dread_f(did, H5T_NATIVE_DOUBLE, var(:,:,:), dm(1:3), iret)

! check if the dataset has been read successfuly
!
    if (iret > 0) then

! print error about the problem with reading the dataset
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot read dataset: " // trim(name)

    end if

! close the dataset
!
    call h5dclose_f(did, iret)

! check if the dataset has been closed successfuly
!
    if (iret > 0) then

! print error about the problem with closing the dataset
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot close dataset: " // trim(name)

    end if

!-------------------------------------------------------------------------------
!
  end subroutine read_3d_array_double_h5
!
!===============================================================================
!
! subroutine READ_4D_ARRAY_DOUBLE_H5:
! -----------------------------------
!
!   Subroutine restores a four-dimensional double precision array from a group
!   specified by identifier.
!
!   Arguments:
!
!     gid   - the HDF5 group identifier
!     name  - the string name describing the array
!     dm    - the array dimensions
!     value - the array values
!
!===============================================================================
!
  subroutine read_4d_array_double_h5(gid, name, dm, var)

! import procedures and variables from other modules
!
    use hdf5           , only : H5T_NATIVE_DOUBLE
    use hdf5           , only : hid_t, hsize_t
    use hdf5           , only : h5dopen_f, h5dread_f, h5dclose_f
    use iso_fortran_env, only : error_unit

! local variables are not implicit by default
!
    implicit none

! subroutine arguments
!
    integer(hid_t)                      , intent(in)    :: gid
    character(len=*)                    , intent(in)    :: name
    integer(hsize_t), dimension(4)      , intent(inout) :: dm
    real(kind=8)    , dimension(:,:,:,:), intent(inout) :: var

! local variables
!
    integer(hid_t) :: did
    integer        :: iret

! local parameters
!
    character(len=*), parameter :: loc = 'IO::read_4d_array_double_h5()'
!
!-------------------------------------------------------------------------------
!
! open the dataset
!
    call h5dopen_f(gid, name, did, iret)

! check if the dataset has been opened successfuly
!
    if (iret < 0) then

! print error about the problem with opening the data space
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot open dataset: " // trim(name)

! quit the subroutine
!
      return

    end if

! read dataset data
!
    call h5dread_f(did, H5T_NATIVE_DOUBLE, var(:,:,:,:), dm(1:4), iret)

! check if the dataset has been read successfuly
!
    if (iret > 0) then

! print error about the problem with reading the dataset
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot read dataset: " // trim(name)

    end if

! close the dataset
!
    call h5dclose_f(did, iret)

! check if the dataset has been closed successfuly
!
    if (iret > 0) then

! print error about the problem with closing the dataset
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot close dataset: " // trim(name)

    end if

!-------------------------------------------------------------------------------
!
  end subroutine read_4d_array_double_h5
!
!===============================================================================
!
! subroutine READ_5D_ARRAY_DOUBLE_H5:
! -----------------------------------
!
!   Subroutine restores a five-dimensional double precision array from a group
!   specified by identifier.
!
!   Arguments:
!
!     gid   - the HDF5 group identifier
!     name  - the string name describing the array
!     dm    - the array dimensions
!     value - the array values
!
!===============================================================================
!
  subroutine read_5d_array_double_h5(gid, name, dm, var)

! import procedures and variables from other modules
!
    use hdf5           , only : H5T_NATIVE_DOUBLE
    use hdf5           , only : hid_t, hsize_t
    use hdf5           , only : h5dopen_f, h5dread_f, h5dclose_f
    use iso_fortran_env, only : error_unit

! local variables are not implicit by default
!
    implicit none

! subroutine arguments
!
    integer(hid_t)                        , intent(in)    :: gid
    character(len=*)                      , intent(in)    :: name
    integer(hsize_t), dimension(5)        , intent(inout) :: dm
    real(kind=8)    , dimension(:,:,:,:,:), intent(inout) :: var

! local variables
!
    integer(hid_t) :: did
    integer        :: iret

! local parameters
!
    character(len=*), parameter :: loc = 'IO::read_5d_array_double_h5()'
!
!-------------------------------------------------------------------------------
!
! open the dataset
!
    call h5dopen_f(gid, name, did, iret)

! check if the dataset has been opened successfuly
!
    if (iret < 0) then

! print error about the problem with opening the data space
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot open dataset: " // trim(name)

! quit the subroutine
!
      return

    end if

! read dataset data
!
    call h5dread_f(did, H5T_NATIVE_DOUBLE, var(:,:,:,:,:), dm(1:5), iret)

! check if the dataset has been read successfuly
!
    if (iret > 0) then

! print error about the problem with reading the dataset
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot read dataset: " // trim(name)

    end if

! close the dataset
!
    call h5dclose_f(did, iret)

! check if the dataset has been closed successfuly
!
    if (iret > 0) then

! print error about the problem with closing the dataset
!
      write(error_unit,"('[',a,']: ',a)") trim(loc)                            &
                      , "Cannot close dataset: " // trim(name)

    end if

!-------------------------------------------------------------------------------
!
  end subroutine read_5d_array_double_h5
!
!===============================================================================
!
! subroutine WRITE_SNAPSHOT_XDMF:
! ------------------------------
!
!   Subroutine writes an XDMF file per snapshot per MPI process.
!   XDMF file is just a wrapper which helps to load data in a visualization
!   tools like Paraview or Visit.
!
!   Based on the subroutine by Pierre Kestener
!   (see https://bitbucket.org/pkestene/amun-code).
!
!
!===============================================================================
!
  subroutine write_snapshot_xdmf()

! import procedures and variables from other modules
!
    use blocks         , only : block_data, list_data
    use blocks         , only : get_dblocks
    use equations      , only : nv, pvars
    use mpitools       , only : nproc
    use coordinates    , only : ni => ncells, ng => nghosts
    use coordinates    , only : adx, ady, adz
    use evolution      , only : time

! local variables are not implicit by default
!
    implicit none

! local pointers
!
    type(block_data), pointer  :: pdata

! local variables
!
    character(len=64)  :: fname, hname
    character(len=128) :: stmp, ttmp, sdim, bdim, pdim
    integer(kind=4)    :: l, p
    integer(kind=4)    :: ip, jp, kp

! local arrays
!
    integer, dimension(12) :: slab

! local parameters
!
    integer, parameter :: xdmf = 101
!
!-------------------------------------------------------------------------------
!
! prepare the XDMF and HDF5 file names

    write(fname, "(a1,i6.6,'_',i5.5,'.xdmf')") ftype, isnap, nproc
    write(hname, "(a1,i6.6,'_',i5.5,'.h5'  )") ftype, isnap, nproc

! open the XDMF file
!
    open (unit = xdmf, file = fname, status = 'replace')

! write the header
!
    write(xdmf, "(a)") '<?xml version="1.0" encoding="UTF-8"?>'
    write(xdmf, "(a)") '<Xdmf Version="2.2"'                                   &
                               // ' xmlns:xi="http://www.w3.org/2003/XInclude">'
    write(xdmf, "(a)") ' <Domain>'
    write(stmp, "(1i16)") nproc
    write(xdmf, "(a)") '  <Grid Name="region_' // trim(adjustl(stmp))          &
                          // '" GridType="Collection" CollectionType="Spatial">'
    write(stmp, "(1g15.8)") time
    write(xdmf, "(a)") '   <Time TimeType="Single" Value="'                    &
                                                 // trim(adjustl(stmp)) // '"/>'

! do not proceed if there are not data block on this processor
!
    if (get_dblocks() > 0) then

! prepare dimensions
!
      ip = ni + 1
      jp = ni + 1
#if NDIMS == 3
      kp = ni + 1
#endif /* NDIMS == 3 */

#if NDIMS == 3
      write(stmp, "(1i8)") ni
      write(ttmp, "(1i8)") ni
      stmp = trim(adjustl(stmp)) // ' ' // trim(adjustl(ttmp))
      write(ttmp, "(1i8)") ni
#else /* NDIMS == 3 */
      write(stmp, "(1i8)") ni
      write(ttmp, "(1i8)") ni
#endif /* NDIMS == 3 */
      bdim = trim(adjustl(stmp)) // ' ' // trim(adjustl(ttmp))

      write(stmp, "(1i8)") ni
      write(ttmp, "(1i8)") ni
      stmp = trim(adjustl(stmp)) // ' ' // trim(adjustl(ttmp))
      write(ttmp, "(1i8)") ni
      stmp = trim(adjustl(stmp)) // ' ' // trim(adjustl(ttmp))
      write(ttmp, "(1i8)") get_dblocks()
      sdim = trim(adjustl(stmp)) // ' ' // trim(adjustl(ttmp))

! prepare slab indices
!
#if NDIMS == 3
      if (with_ghosts) then
        slab(:) = (/ ng, ng, ng, -1, 1, 1, 1, 1, ni, ni, ni, 1 /)
      else
        slab(:) = (/  0,  0,  0, -1, 1, 1, 1, 1, ni, ni, ni, 1 /)
      end if
#else /* NDIMS == 3 */
      if (with_ghosts) then
        slab(:) = (/  0, ng, ng, -1, 1, 1, 1, 1,  1, ni, ni, 1 /)
      else
        slab(:) = (/  0,  0,  0, -1, 1, 1, 1, 1,  1, ni, ni, 1 /)
      end if
#endif /* NDIMS == 3 */

! iterate over all data blocks
!
      l = 0
      pdata => list_data
      do while(associated(pdata))

! store block geometry information
!
        write(stmp, "(1i16)") pdata%meta%id
        write(xdmf, "(a)") '    <Grid Name="block_'                            &
                              // trim(adjustl(stmp)) // '">'
#if NDIMS == 3
        write(stmp, "(1i8)") kp
        write(ttmp, "(1i8)") jp
        stmp = trim(adjustl(stmp)) // ' ' // trim(adjustl(ttmp))
        write(ttmp, "(1i8)") ip
        stmp = trim(adjustl(stmp)) // ' ' // trim(adjustl(ttmp))
        write(xdmf, "(a)") '     <Topology TopologyType="3DCoRectMesh"'        &
                              // ' Dimensions="' // trim(adjustl(stmp)) // '"/>'
        write(xdmf, "(a)") '     <Geometry GeometryType="ORIGIN_DXDYDZ">'
        write(stmp, "(3es16.8)") pdata%meta%zmin, pdata%meta%ymin              &
                                                , pdata%meta%xmin
        write(xdmf, "(a)") '      <DataItem Name="Origin" NumberType="Float"'  &
                              // ' Precision="4" Dimensions="3" Format="XML">' &
                              // trim(adjustl(stmp)) // '</DataItem>'
        write(stmp, "(3es16.8)") adz(pdata%meta%level), ady(pdata%meta%level)  &
                                                      , adx(pdata%meta%level)
        write(xdmf, "(a)") '      <DataItem Name="Spacing" NumberType="Float"' &
                              // ' Precision="4" Dimensions="3" Format="XML">' &
                              // trim(adjustl(stmp)) // '</DataItem>'
#else /* NDIMS == 3 */
        write(stmp, "(1i8)") jp
        write(ttmp, "(1i8)") ip
        stmp = trim(adjustl(stmp)) // ' ' // trim(adjustl(ttmp))
        write(xdmf, "(a)") '     <Topology TopologyType="2DCoRectMesh"'        &
                              // ' Dimensions="' // trim(adjustl(stmp)) // '"/>'
        write(xdmf, "(a)") '     <Geometry GeometryType="ORIGIN_DXDY">'
        write(stmp, "(2es16.8)") pdata%meta%ymin, pdata%meta%xmin
        write(xdmf, "(a)") '      <DataItem Name="Origin" NumberType="Float"'  &
                              // ' Precision="4" Dimensions="2" Format="XML">' &
                              // trim(adjustl(stmp)) // '</DataItem>'
        write(stmp, "(2es16.8)") ady(pdata%meta%level), adx(pdata%meta%level)
        write(xdmf, "(a)") '      <DataItem Name="Spacing" NumberType="Float"' &
                              // ' Precision="4" Dimensions="2" Format="XML">' &
                              // trim(adjustl(stmp)) // '</DataItem>'
#endif /* NDIMS == 3 */
        write(xdmf, "(a)") '     </Geometry>'

! convert slab dimensions to string
!
        slab(4) = l
        write(pdim, "(1i8)") slab(1)
        do p = 2, 12
          write(ttmp, "(1i8)") slab(p)
          pdim = trim(adjustl(pdim)) // ' ' // trim(adjustl(ttmp))
        end do

! loop over all variables and store their information
!
        do p = 1, nv

          write(xdmf, "(a)") '     <Attribute Name="'                          &
                               // trim(adjustl(pvars(p)))                      &
                               // '" AttributeType="Scalar" Center="Cell">'
          write(xdmf, "(a)") '      <DataItem ItemType="Hyperslab"'            &
                               // ' Dimensions="'  // trim(adjustl(bdim))      &
                               // '" Type="HyperSlab">'
          write(xdmf, "(a)") '       <DataItem Dimensions="3 4" Format="XML">' &
                               // trim(adjustl(pdim)) // '</DataItem>'
          write(ttmp, "(a,':/variables/',a)") trim(hname), trim(pvars(p))
          write(xdmf, "(a)") '       <DataItem NumberType="Float"'             &
                               // ' Precision="8" Dimensions="'                &
                               // trim(adjustl(sdim)) // '" Format="HDF">'     &
                               // trim(adjustl(ttmp)) // '</DataItem>'
          write(xdmf, "(a)") '      </DataItem>'
          write(xdmf, "(a)") '     </Attribute>'

        end do

! close grid structure for the current block
!
        write(xdmf,"(a)") '    </Grid>'

        l = l + 1
        pdata => pdata%next
      end do

    end if ! get_dblocks() > 0

! close the XDMF structures
!
    write(xdmf, "(a)") '  </Grid>'
    write(xdmf, "(a)") ' </Domain>'
    write(xdmf, "(a)") '</Xdmf>'

! close the XDMF file
!
    close(xdmf)

!-------------------------------------------------------------------------------
!
  end subroutine write_snapshot_xdmf
!
!===============================================================================
!
! subroutine WRITE_SNAPSHOT_XDMF_MASTER:
! -------------------------------------
!
! Subroutine writes one XDMF file per snapshot in root MPI process to gather
! all MPI subdomains.
!
!   Based on the subroutine by Pierre Kestener
!   (see https://bitbucket.org/pkestene/amun-code).
!
!
!===============================================================================
!
  subroutine write_snapshot_xdmf_master()

! import procedures and variables from other modules
!
    use mpitools       , only : npmax

! local variables are not implicit by default
!
    implicit none

! local variables
!
    character(len=64)  :: fname, pname
    integer(kind=4)    :: np

! local parameters
!
    integer, parameter :: xdmf = 102
!
!-------------------------------------------------------------------------------
!
! prepare the XDMF and HDF5 file names

    write(fname, "(a1,i6.6,'.xdmf')") ftype, isnap

! open the XDMF file
!
    open (unit = xdmf, file = fname, status = 'replace')

! write the header
!
    write(xdmf, "(a)") '<?xml version="1.0" encoding="UTF-8"?>'
    write(xdmf, "(a)") '<Xdmf Version="2.2"'                                   &
                               // ' xmlns:xi="http://www.w3.org/2003/XInclude">'
    write(xdmf, "(a)") ' <Domain Name="GatherMPISubDomains">'
    write(xdmf, "(a)") '  <Grid Name="FullDomain"'                             &
                           // ' GridType="Collection" CollectionType="Spatial">'

! write references to MPI subdomain files
!
    do np = 0, npmax
      write(pname, "(a1,i6.6,'_',i5.5,'.xdmf')") ftype, isnap, np
      write(xdmf, "(a)") '   <xi:include href="' // trim(pname)                &
                               // '" xpointer="xpointer(//Xdmf/Domain/Grid)"/>'
    end do

! close the XDMF structures
!
    write(xdmf, "(a)") '  </Grid>'
    write(xdmf, "(a)") ' </Domain>'
    write(xdmf, "(a)") '</Xdmf>'

! close the XDMF file
!
    close(xdmf)

!-------------------------------------------------------------------------------
!
  end subroutine write_snapshot_xdmf_master
#endif /* HDF5 */

!===============================================================================
!
end module