[petsc-users] Error Writing Large HDF5 File

David Scott d.scott at epcc.ed.ac.uk
Fri Nov 26 09:10:23 CST 2021


Hello,

I am trying to write out HDF5 files. The program I have works all right
up to a point but then fails when I double the size of the file that I
am trying to write out. I have attached a file containing the error
message and another file containing the short program that I used to
generate it.

The first of these files reports that Petsc 3.14.2 was used. I know that
this is old but I have obtained the same result with Petsc 3.16.1. Also,
I have observed the same behaviour on two different machines. I tried
using 64 bit indices but it made no difference.

The program runs successfully with the following command line options

-global_dim_x 1344 -global_dim_y 336 -global_dim_z 336

but fails with these

-global_dim_x 2688 -global_dim_y 336 -global_dim_z 336

All the best,

David

The University of Edinburgh is a charitable body, registered in Scotland, with registration number SC005336. Is e buidheann carthannais a th’ ann an Oilthigh Dhùn Èideann, clàraichte an Alba, àireamh clàraidh SC005336.
-------------- next part --------------
 -global_dim_x                 =        2688
 -global_dim_y                 =         336
 -global_dim_z                 =         336
 Writing 1
 Writing 2
 Writing 3
HDF5-DIAG: Error detected in HDF5 (1.12.0) MPI-process 0:
  #000: ../../src/H5D.c line 151 in H5Dcreate2(): unable to create dataset
    major: Dataset
    minor: Unable to initialize object
  #001: ../../src/H5VLcallback.c line 1869 in H5VL_dataset_create(): dataset create failed
    major: Virtual Object Layer
    minor: Unable to create file
  #002: ../../src/H5VLcallback.c line 1835 in H5VL__dataset_create(): dataset create failed
    major: Virtual Object Layer
    minor: Unable to create file
  #003: ../../src/H5VLnative_dataset.c line 75 in H5VL__native_dataset_create(): unable to create dataset
    major: Dataset
    minor: Unable to initialize object
  #004: ../../src/H5Dint.c line 411 in H5D__create_named(): unable to create and link to dataset
    major: Dataset
    minor: Unable to initialize object
  #005: ../../src/H5L.c line 1804 in H5L_link_object(): unable to create new link to object
    major: Links
    minor: Unable to initialize object
  #006: ../../src/H5L.c line 2045 in H5L__create_real(): can't insert link
    major: Links
    minor: Unable to insert object
  #007: ../../src/H5Gtraverse.c line 855 in H5G_traverse(): internal path traversal failed
    major: Symbol table
    minor: Object not found
  #008: ../../src/H5Gtraverse.c line 630 in H5G__traverse_real(): traversal operator failed
    major: Symbol table
    minor: Callback failed
  #009: ../../src/H5L.c line 1851 in H5L__link_cb(): unable to create object
    major: Links
    minor: Unable to initialize object
  #010: ../../src/H5Oint.c line 2522 in H5O_obj_create(): unable to open object
    major: Object header
    minor: Can't open object
  #011: ../../src/H5Doh.c line 301 in H5O__dset_create(): unable to create dataset
    major: Dataset
    minor: Unable to initialize object
  #012: ../../src/H5Dint.c line 1374 in H5D__create(): can't update the metadata cache
    major: Dataset
    minor: Unable to initialize object
  #013: ../../src/H5Dint.c line 1071 in H5D__update_oh_info(): unable to update layout/pline/efl header message
    major: Dataset
    minor: Unable to initialize object
  #014: ../../src/H5Dlayout.c line 508 in H5D__layout_oh_create(): unable to initialize storage
    major: Dataset
    minor: Unable to initialize object
  #015: ../../src/H5Dint.c line 2387 in H5D__alloc_storage(): unable to initialize dataset with fill value
    major: Dataset
    minor: Unable to initialize object
  #016: ../../src/H5Dint.c line 2474 in H5D__init_storage(): unable to allocate all chunks of dataset
    major: Dataset
    minor: Unable to initialize object
  #017: ../../src/H5Dchunk.c line 4737 in H5D__chunk_allocate(): unable to write raw data to file
    major: Low-level I/O
    minor: Write failed
  #018: ../../src/H5Dchunk.c line 5081 in H5D__chunk_collective_fill(): MPI_Type_free failed
    major: Internal error (too specific to document in detail)
    minor: Some MPI function failed
  #019: ../../src/H5Dchunk.c line 5081 in H5D__chunk_collective_fill(): Invalid datatype, error stack:
PMPI_Type_free(149): MPI_Type_free(datatype_p=0x7fff23a2b938) failed
PMPI_Type_free(83).: Invalid datatype
    major: Internal error (too specific to document in detail)
    minor: MPI Error String
  #020: ../../src/H5Dchunk.c line 5079 in H5D__chunk_collective_fill(): MPI_Type_free failed
    major: Internal error (too specific to document in detail)
    minor: Some MPI function failed
  #021: ../../src/H5Dchunk.c line 5079 in H5D__chunk_collective_fill(): Invalid datatype, error stack:
PMPI_Type_free(149): MPI_Type_free(datatype_p=0x7fff23a2b940) failed
PMPI_Type_free(83).: Invalid datatype
    major: Internal error (too specific to document in detail)
    minor: MPI Error String
  #022: ../../src/H5Dchunk.c line 5038 in H5D__chunk_collective_fill(): MPI_Type_create_hindexed failed
    major: Internal error (too specific to document in detail)
    minor: Some MPI function failed
  #023: ../../src/H5Dchunk.c line 5038 in H5D__chunk_collective_fill(): Invalid argument, error stack:
PMPI_Type_create_hindexed(143): MPI_Type_create_hindexed(count=1, array_of_blocklengths=0x1eba120, array_of_displacements=0x1eba140, MPI_BYTE, newtype=0x7fff23a2b940) failed
PMPI_Type_create_hindexed(96).: Invalid value for blocklength, must be non-negative but is -1860026368
    major: Internal error (too specific to document in detail)
    minor: MPI Error String
[0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------
[0]PETSC ERROR: Error in external library
[0]PETSC ERROR: Error in HDF5 call H5Dcreate2() Status -1
[0]PETSC ERROR: See https://www.mcs.anl.gov/petsc/documentation/faq.html for trouble shooting.
[0]PETSC ERROR: Petsc Release Version 3.14.2, Dec 03, 2020 
[0]PETSC ERROR: hdf5.x on a  named nid005829 by dmse643 Fri Nov 26 12:54:34 2021
[0]PETSC ERROR: Configure options --known-has-attribute-aligned=1 --known-mpi-int64_t=0 --known-bits-per-byte=8 --known-64-bit-blas-indices=0 --known-sdot-returns-double=0 --known-snrm2-returns-double=0 --known-level1-dcache-assoc=4 --known-level1-dcache-linesize=64 --known-level1-dcache-size=16384 --known-memcmp-ok=1 --known-mpi-c-double-complex=1 --known-mpi-long-double=1 --known-mpi-shared-libraries=0 --known-sizeof-MPI_Comm=4 --known-sizeof-MPI_Fint=4 --known-sizeof-char=1 --known-sizeof-double=8 --known-sizeof-float=4 --known-sizeof-int=4 --known-sizeof-long-long=8 --known-sizeof-long=8 --known-sizeof-short=2 --known-sizeof-size_t=8 --known-sizeof-void-p=8 --with-ar=ar --with-batch=0 --with-cc=cc --with-clib-autodetect=0 --with-cxx=CC --with-cxxlib-autodetect=0 --with-debugging=0 --with-dependencies=0 --with-fc=ftn --with-fortran-datatypes=1 --with-fortran-interfaces=1 --with-fortran-bindings=1 --with-fortranlib-autodetect=0 --with-ranlib=ranlib --with-scalar-type=real --with-shared-ld=ar --with-etags=0 --with-x=0 --with-ssl=0 --with-shared-libraries=0 --with-mpi-lib="[]" --with-mpi-include="[]" --with-mpiexec=srun --with-blaslapack=1 --with-superlu=1 --with-superlu_dist=1 --with-parmetis=1 --with-metis=1 --with-hypre=1 --with-scalapack=1 --with-ptscotch=1 --with-ptscotch-include= --with-ptscotch-lib= --with-mumps=1 --with-mumps-include= --with-mumps-lib=-lmpifort --with-hdf5=1 --CFLAGS="-O3 -ffast-math  -fPIC " --CPPFLAGS="-I/work/y07/shared/libs/core/petsc/3.14.2/GNU/9.3/include " --CXXFLAGS="-O3 -ffast-math  -fPIC " --with-cxx-dialect=C++11 --FFLAGS="-O3 -ffast-math  -fPIC " --LDFLAGS="-L/work/y07/shared/libs/core/petsc/3.14.2/GNU/9.3/lib " --LIBS="-lgfortran -lgcc  -lstdc++" --PETSC_ARCH=x86-rome --prefix=/work/y07/shared/libs/core/petsc/3.14.2/GNU/9.3
[0]PETSC ERROR: #1 VecView_MPI_HDF5_DA() line 518 in /mnt/lustre/a2fs-work1/work/y07/y07/cse/pe-scripts/petsc-3.14.2/src/dm/impls/da/gr2.c
[0]PETSC ERROR: #2 VecView_MPI_DA() line 698 in /mnt/lustre/a2fs-work1/work/y07/y07/cse/pe-scripts/petsc-3.14.2/src/dm/impls/da/gr2.c
[0]PETSC ERROR: #3 VecView() line 608 in /mnt/lustre/a2fs-work1/work/y07/y07/cse/pe-scripts/petsc-3.14.2/src/vec/vec/interface/vector.c
 Writing 4
HDF5-DIAG: Error detected in HDF5 (1.12.0) MPI-process 0:
  #000: ../../src/H5F.c line 886 in H5Fclose(): decrementing file ID failed
    major: Object atom
    minor: Unable to close file
  #001: ../../src/H5I.c line 1422 in H5I_dec_app_ref(): can't decrement ID ref count
    major: Object atom
    minor: Unable to decrement reference count
  #002: ../../src/H5F.c line 243 in H5F__close_cb(): unable to close file
    major: File accessibility
    minor: Unable to close file
  #003: ../../src/H5VLcallback.c line 3977 in H5VL_file_close(): file close failed
    major: Virtual Object Layer
    minor: Unable to close file
  #004: ../../src/H5VLcallback.c line 3945 in H5VL__file_close(): file close failed
    major: Virtual Object Layer
    minor: Unable to close file
  #005: ../../src/H5VLnative_file.c line 884 in H5VL__native_file_close(): can't close file
    major: File accessibility
    minor: Unable to decrement reference count
  #006: ../../src/H5Fint.c line 2049 in H5F__close(): can't close file, there are objects still open
    major: File accessibility
    minor: Unable to close file
[0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------
[0]PETSC ERROR: Error in external library
[0]PETSC ERROR: Error in HDF5 call H5Fclose() Status -1
[0]PETSC ERROR: See https://www.mcs.anl.gov/petsc/documentation/faq.html for trouble shooting.
[0]PETSC ERROR: Petsc Release Version 3.14.2, Dec 03, 2020 
[0]PETSC ERROR: hdf5.x on a  named nid005829 by dmse643 Fri Nov 26 12:54:34 2021
[0]PETSC ERROR: Configure options --known-has-attribute-aligned=1 --known-mpi-int64_t=0 --known-bits-per-byte=8 --known-64-bit-blas-indices=0 --known-sdot-returns-double=0 --known-snrm2-returns-double=0 --known-level1-dcache-assoc=4 --known-level1-dcache-linesize=64 --known-level1-dcache-size=16384 --known-memcmp-ok=1 --known-mpi-c-double-complex=1 --known-mpi-long-double=1 --known-mpi-shared-libraries=0 --known-sizeof-MPI_Comm=4 --known-sizeof-MPI_Fint=4 --known-sizeof-char=1 --known-sizeof-double=8 --known-sizeof-float=4 --known-sizeof-int=4 --known-sizeof-long-long=8 --known-sizeof-long=8 --known-sizeof-short=2 --known-sizeof-size_t=8 --known-sizeof-void-p=8 --with-ar=ar --with-batch=0 --with-cc=cc --with-clib-autodetect=0 --with-cxx=CC --with-cxxlib-autodetect=0 --with-debugging=0 --with-dependencies=0 --with-fc=ftn --with-fortran-datatypes=1 --with-fortran-interfaces=1 --with-fortran-bindings=1 --with-fortranlib-autodetect=0 --with-ranlib=ranlib --with-scalar-type=real --with-shared-ld=ar --with-etags=0 --with-x=0 --with-ssl=0 --with-shared-libraries=0 --with-mpi-lib="[]" --with-mpi-include="[]" --with-mpiexec=srun --with-blaslapack=1 --with-superlu=1 --with-superlu_dist=1 --with-parmetis=1 --with-metis=1 --with-hypre=1 --with-scalapack=1 --with-ptscotch=1 --with-ptscotch-include= --with-ptscotch-lib= --with-mumps=1 --with-mumps-include= --with-mumps-lib=-lmpifort --with-hdf5=1 --CFLAGS="-O3 -ffast-math  -fPIC " --CPPFLAGS="-I/work/y07/shared/libs/core/petsc/3.14.2/GNU/9.3/include " --CXXFLAGS="-O3 -ffast-math  -fPIC " --with-cxx-dialect=C++11 --FFLAGS="-O3 -ffast-math  -fPIC " --LDFLAGS="-L/work/y07/shared/libs/core/petsc/3.14.2/GNU/9.3/lib " --LIBS="-lgfortran -lgcc  -lstdc++" --PETSC_ARCH=x86-rome --prefix=/work/y07/shared/libs/core/petsc/3.14.2/GNU/9.3
[0]PETSC ERROR: #4 PetscViewerFileClose_HDF5() line 79 in /mnt/lustre/a2fs-work1/work/y07/y07/cse/pe-scripts/petsc-3.14.2/src/sys/classes/viewer/impls/hdf5/hdf5v.c
[0]PETSC ERROR: #5 PetscViewerDestroy_HDF5() line 90 in /mnt/lustre/a2fs-work1/work/y07/y07/cse/pe-scripts/petsc-3.14.2/src/sys/classes/viewer/impls/hdf5/hdf5v.c
 Writing 5
[0]PETSC ERROR: #6 PetscViewerDestroy() line 118 in /mnt/lustre/a2fs-work1/work/y07/y07/cse/pe-scripts/petsc-3.14.2/src/sys/classes/viewer/interface/view.c
 Writing 6
 Writing 7
 Writing 8
-------------- next part --------------
program hdf5_test_program

  use petsc

  implicit none

#include "petsc/finclude/petsc.h"

  integer :: ierr
  logical :: found
  DM :: da_w
  Vec :: w_vec
  PetscViewer :: viewer
  integer :: num_procs
  PetscInt :: global_dim_x, global_dim_y, global_dim_z
  PetscInt, parameter :: dof = 1
  PetscInt, parameter :: stencil_width = 1
  PetscInt, dimension(0:0) :: x_lengths
  PetscInt, dimension(0:0) :: y_lengths
  PetscInt, dimension(0:0) :: z_lengths_w
  character(len = 30) :: option_name
  character(len = 60) :: msg

  call PetscInitialize(PETSC_NULL_CHARACTER, ierr)

  call MPI_Comm_size(PETSC_COMM_WORLD, num_procs, ierr)
  if (num_procs .ne. 1) then
    write(*, *) 'The number of processes must be 1'
    call MPI_Abort(PETSC_COMM_WORLD, -1, ierr)
  end if

  ! Read command-line options and options file.
  option_name = '-global_dim_x'
  call PetscOptionsGetInt(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, option_name, global_dim_x, found, ierr)
  if (found) then
    write(msg, *) option_name, '=', global_dim_x
    call PetscPrintf(PETSC_COMM_WORLD, trim(msg) // NEW_LINE('a'), ierr)
  else
    write(msg, *) 'Error:', option_name, 'not found.'
    call PetscPrintf(PETSC_COMM_WORLD, trim(msg) // NEW_LINE('a'), ierr)
  end if

  option_name = '-global_dim_y'
  call PetscOptionsGetInt(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, option_name, global_dim_y, found, ierr)
  if (found) then
    write(msg, *) option_name, '=', global_dim_y
    call PetscPrintf(PETSC_COMM_WORLD, trim(msg) // NEW_LINE('a'), ierr)
  else
    write(msg, *) 'Error:', option_name, 'not found.'
    call PetscPrintf(PETSC_COMM_WORLD, trim(msg) // NEW_LINE('a'), ierr)
  end if

  option_name = '-global_dim_z'
  call PetscOptionsGetInt(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, option_name, global_dim_z, found, ierr)
  if (found) then
    write(msg, *) option_name, '=', global_dim_z
    call PetscPrintf(PETSC_COMM_WORLD, trim(msg) // NEW_LINE('a'), ierr)
    if(global_dim_z .le. 1) then
      write(msg, *) 'global_dim_z must be greater than 1.'
      call PetscPrintf(PETSC_COMM_WORLD, trim(msg) // NEW_LINE('a'), ierr)
    end if
  else
    write(msg, *) 'Error:', option_name, 'not found.'
    call PetscPrintf(PETSC_COMM_WORLD, trim(msg) // NEW_LINE('a'), ierr)
  end if

  x_lengths = global_dim_x
  y_lengths = global_dim_y
  z_lengths_w = global_dim_z+1

  call DMDACreate3d(PETSC_COMM_WORLD,                         &
       DM_BOUNDARY_PERIODIC,                                  &
       DM_BOUNDARY_PERIODIC,                                  &
       DM_BOUNDARY_NONE,                                      &
       DMDA_STENCIL_BOX,                                      &
       global_dim_x, global_dim_y, global_dim_z+1,            &
       1, 1, 1,                                               &
       dof, stencil_width,                                    &
       x_lengths,                                             &
       y_lengths,                                             &
       z_lengths_w,                                           &
       da_w, ierr);CHKERRA(ierr)
  call DMSetUp(da_w, ierr);CHKERRA(ierr)

  call DMCreateGlobalVector(da_w, w_vec, ierr)
  call VecZeroEntries(w_vec, ierr)

  call PetscViewerHDF5Open(PETSC_COMM_WORLD, './w.h5', FILE_MODE_WRITE, viewer, ierr)
  call PetscViewerPushFormat(viewer, PETSC_VIEWER_HDF5_XDMF, ierr)

  ! Write the Vec without one extra dimension for BS
  write(*, *) 'Writing 1'
  call PetscViewerHDF5SetBaseDimension2(viewer, PETSC_FALSE, ierr)
  write(*, *) 'Writing 2'
  call PetscObjectSetName(w_vec, 'w', ierr)
  write(*, *) 'Writing 3'
  call VecView(w_vec, viewer, ierr)

  write(*, *) 'Writing 4'
  call PetscViewerDestroy(viewer, ierr)

  write(*, *) 'Writing 5'
  call VecDestroy(w_vec, ierr)

  write(*, *) 'Writing 6'
  call DMDestroy(da_w, ierr);CHKERRA(ierr)
  
  write(*, *) 'Writing 7'
  call PetscFinalize(ierr)

  write(*, *) 'Writing 8'

end program hdf5_test_program


More information about the petsc-users mailing list