47  use mpi_f08, 
only : mpi_info_null, mpi_allreduce, mpi_in_place, &
 
   70    class(*), 
target, 
intent(in) :: data
 
   71    real(kind=
rp), 
intent(in), 
optional :: t
 
   72    type(
mesh_t), 
pointer :: msh
 
   76    real(kind=
rp), 
pointer :: dtlag(:)
 
   77    real(kind=
rp), 
pointer :: tlag(:)
 
   78    integer :: ierr, info, drank, i, j
 
   79    integer(hid_t) :: plist_id, file_id, dset_id, grp_id, attr_id
 
   80    integer(hid_t) :: filespace, memspace
 
   81    integer(hid_t) :: h5t_neko_real
 
   82    integer(hsize_t), 
dimension(1) :: ddim, dcount, doffset
 
   84    character(len=5) :: id_str
 
   85    character(len=1024) :: fname
 
   87    call hdf5_file_determine_data(
data, msh, dof, fp, fsp, dtlag, tlag)
 
   89    if (.not. this%overwrite) 
call this%increment_counter()
 
   90    fname = trim(this%get_fname())
 
   94    call hdf5_file_determine_real(h5t_neko_real)
 
   96    call h5pcreate_f(h5p_file_access_f, plist_id, ierr)
 
   97    info = mpi_info_null%mpi_val
 
   98    call h5pset_fapl_mpio_f(plist_id, 
neko_comm%mpi_val, info, ierr)
 
  100    call h5fcreate_f(fname, h5f_acc_trunc_f, &
 
  101         file_id, ierr, access_prp = plist_id)
 
  103    call h5pcreate_f(h5p_dataset_xfer_f, plist_id, ierr)
 
  104    call h5pset_dxpl_mpio_f(plist_id, h5fd_mpio_collective_f, ierr)
 
  106    call h5screate_f(h5s_scalar_f, filespace, ierr)
 
  110       call h5acreate_f(file_id, 
"Time", h5t_neko_real, filespace, attr_id, &
 
  111            ierr, h5p_default_f, h5p_default_f)
 
  112       call h5awrite_f(attr_id, h5t_neko_real, t, ddim, ierr)
 
  113       call h5aclose_f(attr_id, ierr)
 
  116    if (
associated(dof)) 
then 
  117       call h5acreate_f(file_id, 
"Lx", h5t_native_integer, filespace, attr_id, &
 
  118            ierr, h5p_default_f, h5p_default_f)
 
  119       call h5awrite_f(attr_id, h5t_native_integer, dof%Xh%lx, ddim, ierr)
 
  120       call h5aclose_f(attr_id, ierr)
 
  123    if (
associated(msh)) 
then 
  124       call h5gcreate_f(file_id, 
"Mesh", grp_id, ierr, lcpl_id=h5p_default_f, &
 
  125            gcpl_id=h5p_default_f, gapl_id=h5p_default_f)
 
  127       call h5acreate_f(grp_id, 
"Elements", h5t_native_integer, filespace, attr_id, &
 
  128            ierr, h5p_default_f, h5p_default_f)
 
  129       call h5awrite_f(attr_id, h5t_native_integer, msh%glb_nelv, ddim, ierr)
 
  130       call h5aclose_f(attr_id, ierr)
 
  132       call h5acreate_f(grp_id, 
"Dimension", h5t_native_integer, filespace, attr_id, &
 
  133            ierr, h5p_default_f, h5p_default_f)
 
  134       call h5awrite_f(attr_id, h5t_native_integer, msh%gdim, ddim, ierr)
 
  135       call h5aclose_f(attr_id, ierr)
 
  137       call h5gclose_f(grp_id, ierr)
 
  141    call h5sclose_f(filespace, ierr)
 
  146    if (
associated(tlag) .and. 
associated(dtlag)) 
then 
  147       call h5gcreate_f(file_id, 
"Restart", grp_id, ierr, lcpl_id=h5p_default_f, &
 
  148            gcpl_id=h5p_default_f, gapl_id=h5p_default_f)
 
  159       call h5screate_simple_f(drank, ddim, filespace, ierr)
 
  161       call h5dcreate_f(grp_id,
'tlag', h5t_neko_real, &
 
  162            filespace, dset_id, ierr)
 
  163       call h5dget_space_f(dset_id, filespace, ierr)
 
  164       call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
 
  165            doffset, dcount, ierr)
 
  166       call h5dwrite_f(dset_id, h5t_neko_real, tlag, &
 
  167            ddim, ierr, xfer_prp = plist_id)
 
  168       call h5dclose_f(dset_id, ierr)
 
  170       call h5dcreate_f(grp_id,
'dtlag', h5t_neko_real, &
 
  171            filespace, dset_id, ierr)
 
  172       call h5dget_space_f(dset_id, filespace, ierr)
 
  173       call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
 
  174            doffset, dcount, ierr)
 
  175       call h5dwrite_f(dset_id, h5t_neko_real, dtlag, &
 
  176            ddim, ierr, xfer_prp = plist_id)
 
  177       call h5dclose_f(dset_id, ierr)
 
  179       call h5sclose_f(filespace, ierr)
 
  180       call h5gclose_f(grp_id, ierr)
 
  188    if (
allocated(fp) .or. 
allocated(fsp)) 
then 
  189       call h5gcreate_f(file_id, 
"Fields", grp_id, ierr, lcpl_id=h5p_default_f, &
 
  190            gcpl_id=h5p_default_f, gapl_id=h5p_default_f)
 
  192       dcount(1) = int(dof%size(), 8)
 
  193       doffset(1) = int(msh%offset_el, 8) * int((dof%Xh%lx**3),8)
 
  194       ddim = int(dof%size(), 8)
 
  196       call mpi_allreduce(mpi_in_place, ddim(1), 1, &
 
  199       call h5screate_simple_f(drank, ddim, filespace, ierr)
 
  200       call h5screate_simple_f(drank, dcount, memspace, ierr)
 
  203       if (
allocated(fp)) 
then 
  205             call h5dcreate_f(grp_id, fp(i)%ptr%name, h5t_neko_real, &
 
  206                  filespace, dset_id, ierr)
 
  207             call h5dget_space_f(dset_id, filespace, ierr)
 
  208             call h5sselect_hyperslab_f(filespace, h5s_select_set_f, &
 
  209                  doffset, dcount, ierr)
 
  210             call h5dwrite_f(dset_id, h5t_neko_real, &
 
  211                  fp(i)%ptr%x(1,1,1,1), &
 
  212                  ddim, ierr, file_space_id = filespace, &
 
  213                  mem_space_id = memspace, xfer_prp = plist_id)
 
  214             call h5dclose_f(dset_id, ierr)
 
  219       if (
allocated(fsp)) 
then 
  221             do j = 1, fsp(i)%ptr%size()
 
  222                call h5dcreate_f(grp_id, fsp(i)%ptr%lf(j)%name, &
 
  223                     h5t_neko_real, filespace, dset_id, ierr)
 
  224                call h5dget_space_f(dset_id, filespace, ierr)
 
  225                call h5sselect_hyperslab_f(filespace, h5s_select_set_f, &
 
  226                     doffset, dcount, ierr)
 
  227                call h5dwrite_f(dset_id, h5t_neko_real, &
 
  228                     fsp(i)%ptr%lf(j)%x(1,1,1,1), &
 
  229                     ddim, ierr, file_space_id = filespace, &
 
  230                     mem_space_id = memspace, xfer_prp = plist_id)
 
  231                call h5dclose_f(dset_id, ierr)
 
  237       call h5gclose_f(grp_id, ierr)
 
  238       call h5sclose_f(filespace, ierr)
 
  239       call h5sclose_f(memspace, ierr)
 
  242    call h5pclose_f(plist_id, ierr)
 
  243    call h5fclose_f(file_id, ierr)
 
  252    class(*), 
target, 
intent(inout) :: data
 
  253    integer(hid_t) :: plist_id, file_id, dset_id, grp_id, attr_id
 
  254    integer(hid_t) :: filespace, memspace
 
  255    integer(hid_t) :: h5t_neko_real
 
  256    integer(hsize_t), 
dimension(1) :: ddim, dcount, doffset
 
  257    integer :: i,j, ierr, info, glb_nelv, gdim, lx, drank
 
  258    type(
mesh_t), 
pointer :: msh
 
  262    real(kind=
rp), 
pointer :: dtlag(:)
 
  263    real(kind=
rp), 
pointer :: tlag(:)
 
  265    character(len=1024) :: fname
 
  267    fname = trim(this%get_fname())
 
  269    call hdf5_file_determine_data(
data, msh, dof, fp, fsp, dtlag, tlag)
 
  273    call hdf5_file_determine_real(h5t_neko_real)
 
  275    call h5pcreate_f(h5p_file_access_f, plist_id, ierr)
 
  276    info = mpi_info_null%mpi_val
 
  277    call h5pset_fapl_mpio_f(plist_id, 
neko_comm%mpi_val, info, ierr)
 
  279    call h5fopen_f(fname, h5f_acc_rdonly_f, &
 
  280         file_id, ierr, access_prp = plist_id)
 
  282    call h5pcreate_f(h5p_dataset_xfer_f, plist_id, ierr)
 
  283    call h5pset_dxpl_mpio_f(plist_id, h5fd_mpio_collective_f, ierr)
 
  286    call h5aopen_name_f(file_id, 
'Time', attr_id, ierr)
 
  287    call h5aread_f(attr_id, h5t_neko_real, t, ddim, ierr)
 
  288    call h5aclose_f(attr_id, ierr)
 
  295    call h5aopen_name_f(file_id, 
'Lx', attr_id, ierr)
 
  296    call h5aread_f(attr_id, h5t_native_integer, lx, ddim, ierr)
 
  297    call h5aclose_f(attr_id, ierr)
 
  299    call h5gopen_f(file_id, 
'Mesh', grp_id, ierr, gapl_id=h5p_default_f)
 
  301    call h5aopen_name_f(grp_id, 
'Elements', attr_id, ierr)
 
  302    call h5aread_f(attr_id, h5t_native_integer, glb_nelv, ddim, ierr)
 
  303    call h5aclose_f(attr_id, ierr)
 
  305    call h5aopen_name_f(grp_id, 
'Dimension', attr_id, ierr)
 
  306    call h5aread_f(attr_id, h5t_native_integer, gdim, ddim, ierr)
 
  307    call h5aclose_f(attr_id, ierr)
 
  308    call h5gclose_f(grp_id, ierr)
 
  311    if (
associated(tlag) .and. 
associated(dtlag)) 
then 
  321       call h5gopen_f(file_id, 
'Restart', grp_id, ierr, gapl_id=h5p_default_f)
 
  322       call h5dopen_f(grp_id, 
'tlag', dset_id, ierr)
 
  323       call h5dget_space_f(dset_id, filespace, ierr)
 
  324       call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
 
  325            doffset, dcount, ierr)
 
  326       call h5dread_f(dset_id, h5t_neko_real, tlag, ddim, ierr, xfer_prp=plist_id)
 
  327       call h5dclose_f(dset_id, ierr)
 
  328       call h5sclose_f(filespace, ierr)
 
  330       call h5dopen_f(grp_id, 
'dtlag', dset_id, ierr)
 
  331       call h5dget_space_f(dset_id, filespace, ierr)
 
  332       call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
 
  333            doffset, dcount, ierr)
 
  334       call h5dread_f(dset_id, h5t_neko_real, dtlag, ddim, ierr, xfer_prp=plist_id)
 
  335       call h5dclose_f(dset_id, ierr)
 
  336       call h5sclose_f(filespace, ierr)
 
  338       call h5gclose_f(grp_id, ierr)
 
  341    if (
allocated(fp) .or. 
allocated(fsp)) 
then 
  342       call h5gopen_f(file_id, 
'Fields', grp_id, ierr, gapl_id=h5p_default_f)
 
  344       dcount(1) = int(dof%size(), 8)
 
  345       doffset(1) = int(msh%offset_el, 8) * int((dof%Xh%lx**3),8)
 
  346       ddim = int(dof%size(), 8)
 
  349       dcount(1) = int(dof%size(), 8)
 
  350       doffset(1) = int(msh%offset_el, 8) * int((dof%Xh%lx**3),8)
 
  351       ddim = int(dof%size(), 8)
 
  353       call mpi_allreduce(mpi_in_place, ddim(1), 1, &
 
  356       call h5screate_simple_f(drank, dcount, memspace, ierr)
 
  358       if (
allocated(fp)) 
then 
  360             call h5dopen_f(grp_id, fp(i)%ptr%name, dset_id, ierr)
 
  361             call h5dget_space_f(dset_id, filespace, ierr)
 
  362             call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
 
  363                  doffset, dcount, ierr)
 
  364             call h5dread_f(dset_id, h5t_neko_real, &
 
  365                  fp(i)%ptr%x(1,1,1,1), &
 
  366                  ddim, ierr, file_space_id = filespace, &
 
  367                  mem_space_id = memspace, xfer_prp=plist_id)
 
  368             call h5dclose_f(dset_id, ierr)
 
  369             call h5sclose_f(filespace, ierr)
 
  373       if (
allocated(fsp)) 
then 
  375             do j = 1, fsp(i)%ptr%size()
 
  376                call h5dopen_f(grp_id, fsp(i)%ptr%lf(j)%name, dset_id, ierr)
 
  377                call h5dget_space_f(dset_id, filespace, ierr)
 
  378                call h5sselect_hyperslab_f (filespace, h5s_select_set_f, &
 
  379                     doffset, dcount, ierr)
 
  380                call h5dread_f(dset_id, h5t_neko_real, &
 
  381                     fsp(i)%ptr%lf(j)%x(1,1,1,1), &
 
  382                     ddim, ierr, file_space_id = filespace, &
 
  383                     mem_space_id = memspace, xfer_prp=plist_id)
 
  384                call h5dclose_f(dset_id, ierr)
 
  385                call h5sclose_f(filespace, ierr)
 
  389       call h5sclose_f(memspace, ierr)
 
  390       call h5gclose_f(grp_id, ierr)
 
  393    call h5pclose_f(plist_id, ierr)
 
  394    call h5fclose_f(file_id, ierr)
 
  401  subroutine hdf5_file_determine_data(data, msh, dof, fp, fsp, dtlag, tlag)
 
  402    class(*), 
target, 
intent(in) :: data
 
  403    type(
mesh_t), 
pointer, 
intent(inout) :: msh
 
  404    type(
dofmap_t), 
pointer, 
intent(inout) :: dof
 
  405    type(
field_ptr_t), 
allocatable, 
intent(inout) :: fp(:)
 
  407    real(kind=
rp), 
pointer, 
intent(inout) :: dtlag(:)
 
  408    real(kind=
rp), 
pointer, 
intent(inout) :: tlag(:)
 
  409    integer :: i, j, fp_size, fp_cur, fsp_size, fsp_cur, scalar_count, ab_count
 
  410    character(len=32) :: scalar_name
 
  417       allocate(fp(fp_size))
 
  425       if (data%size() .gt. 0) 
then 
  426          allocate(fp(data%size()))
 
  431          do i = 1, data%size()
 
  432             fp(i)%ptr => data%items(i)%ptr
 
  443       if ( .not. 
associated(data%u) .or. &
 
  444            .not. 
associated(data%v) .or. &
 
  445            .not. 
associated(data%w) .or. &
 
  446            .not. 
associated(data%p) ) 
then 
  452       if (
allocated(data%scalar_lags%items) .and. data%scalar_lags%size() > 0) 
then 
  453          scalar_count = data%scalar_lags%size()
 
  454       else if (
associated(data%s)) 
then 
  460       if (scalar_count .gt. 1) 
then 
  461          fp_size = fp_size + scalar_count
 
  464          fp_size = fp_size + (scalar_count * 2)
 
  465       else if (
associated(data%s)) 
then 
  467          fp_size = fp_size + 1
 
  468          if (
associated(data%abs1)) 
then 
  469             fp_size = fp_size + 2
 
  473       if (
associated(data%abx1)) 
then 
  474          fp_size = fp_size + 6
 
  477       allocate(fp(fp_size))
 
  480       if (
associated(data%ulag)) 
then 
  481          fsp_size = fsp_size + 3
 
  484       if (scalar_count .gt. 1) 
then 
  485          if (
allocated(data%scalar_lags%items)) 
then 
  486             fsp_size = fsp_size + data%scalar_lags%size()
 
  488       else if (
associated(data%slag)) 
then 
  489          fsp_size = fsp_size + 1
 
  492       if (fsp_size .gt. 0) 
then 
  493          allocate(fsp(fsp_size))
 
  507       if (scalar_count .gt. 1) 
then 
  510            do i = 1, scalar_count
 
  511               slag => data%scalar_lags%get(i)
 
  512               fp(fp_cur)%ptr => slag%f
 
  517          do i = 1, scalar_count
 
  518             fp(fp_cur)%ptr => data%scalar_abx1(i)%ptr
 
  520             fp(fp_cur)%ptr => data%scalar_abx2(i)%ptr
 
  523       else if (
associated(data%s)) 
then 
  525          fp(fp_cur)%ptr => data%s
 
  528          if (
associated(data%abs1)) 
then 
  529             fp(fp_cur)%ptr => data%abs1
 
  530             fp(fp_cur+1)%ptr => data%abs2
 
  535       if (
associated(data%abx1)) 
then 
  536          fp(fp_cur)%ptr => data%abx1
 
  537          fp(fp_cur+1)%ptr => data%abx2
 
  538          fp(fp_cur+2)%ptr => data%aby1
 
  539          fp(fp_cur+3)%ptr => data%aby2
 
  540          fp(fp_cur+4)%ptr => data%abz1
 
  541          fp(fp_cur+5)%ptr => data%abz2
 
  545       if (
associated(data%ulag)) 
then 
  546          fsp(fsp_cur)%ptr => data%ulag
 
  547          fsp(fsp_cur+1)%ptr => data%vlag
 
  548          fsp(fsp_cur+2)%ptr => data%wlag
 
  549          fsp_cur = fsp_cur + 3
 
  553       if (scalar_count .gt. 1) 
then 
  554          if (
allocated(data%scalar_lags%items)) 
then 
  555             do j = 1, data%scalar_lags%size()
 
  556                fsp(fsp_cur)%ptr => data%scalar_lags%get(j)
 
  557                fsp_cur = fsp_cur + 1
 
  560       else if (
associated(data%slag)) 
then 
  561          fsp(fsp_cur)%ptr => data%slag
 
  562          fsp_cur = fsp_cur + 1
 
  565       if (
associated(data%tlag)) 
then 
  574  end subroutine hdf5_file_determine_data
 
  579  subroutine hdf5_file_determine_real(H5T_NEKO_REAL)
 
  580    integer(hid_t), 
intent(inout) :: h5t_neko_real
 
  583       h5t_neko_real = h5t_native_double
 
  585       h5t_neko_real = h5t_native_real
 
  587       call neko_error(
"Unsupported real type")
 
  589  end subroutine hdf5_file_determine_real
 
  594    logical, 
intent(in) :: overwrite
 
  595    this%overwrite = overwrite
 
  603    class(*), 
target, 
intent(in) :: data
 
  604    real(kind=rp), 
intent(in), 
optional :: t
 
  605    call neko_error(
'Neko needs to be built with HDF5 support')
 
 
  611    class(*), 
target, 
intent(inout) :: data
 
  612    call neko_error(
'Neko needs to be built with HDF5 support')
 
 
  618    logical, 
intent(in) :: overwrite
 
  619    call neko_error(
'Neko needs to be built with HDF5 support')
 
 
integer, public pe_rank
MPI rank.
 
type(mpi_comm), public neko_comm
MPI communicator.
 
Defines a mapping of the degrees of freedom.
 
Defines a registry for storing solution fields.
 
type(field_registry_t), target, public neko_field_registry
Global field registry.
 
Contains the field_serties_t type.
 
subroutine hdf5_file_read(this, data)
Read data in HDF5 format.
 
subroutine hdf5_file_write(this, data, t)
Write data in HDF5 format.
 
subroutine hdf5_file_set_overwrite(this, overwrite)
Set the overwrite flag for HDF5 files.
 
type(log_t), public neko_log
Global log stream.
 
integer, parameter, public dp
 
integer, parameter, public sp
 
integer, parameter, public rp
Global precision used in computations.
 
pure integer function, public filename_suffix_pos(fname)
Find position (in the string) of a filename's suffix.
 
field_ptr_t, To easily obtain a pointer to a field
 
field_list_t, To be able to group fields together
 
A wrapper for a pointer to a field_series_t.
 
Stores a series (sequence) of fields, logically connected to a base field, and arranged according to ...
 
Interface for HDF5 files.