hi

I'm having trouble reading in a vector on a DMPlex from an HDF5 file. The attached minimal program (readdmvec2.F90) illustrates the problem.

What it does is:
- reads in a DMPlex from file (ghotcol.exo) and sets it up (distributes, constructs ghost cells and default section)
- uses DMGetGlobalVector() to get the vector
- sets its object name to 'cell_geometry'
- opens the HDF5 file (hotcol_ss.h5) in a viewer and tries to read the vector using VecLoad()

The VecLoad() dies, saying "object 'cell_geometry' doesn't exist" in the HDF5 file. (Full error message is below.)

But it does exist- it shows up as a dataset in the '/' group in h5dump. I can also read it in from a Python script and h5py, e.g.:

>>> import h5py
>>> f = h5py.File('hotcol_ss.h5')
>>> f['cell_geometry'][:]

Any idea why VecLoad() can't find it? I get the same error running on any number of processors.

The program that writes 'cell_geometry' to the HDF5 file also subsequently writes time-dependent simulation results to the same file (that's what the other datasets are in there) using DMSetOutputSequenceNumber() and VecView(). I am wondering if the problem is related to this, because if I just write a simple vector out to the file by itself it seems to read back in ok.

Cheers, Adrian

--
HDF5-DIAG: Error detected in HDF5 (1.8.8) MPI-process 0:
  #000: ../../../src/H5D.c line 334 in H5Dopen2(): not found
    major: Dataset
    minor: Object not found
  #001: ../../../src/H5Gloc.c line 430 in H5G_loc_find(): can't find object
    major: Symbol table
    minor: Object not found
#002: ../../../src/H5Gtraverse.c line 861 in H5G_traverse(): internal path traversal failed
    major: Symbol table
    minor: Object not found
#003: ../../../src/H5Gtraverse.c line 641 in H5G_traverse_real(): traversal operator failed
    major: Symbol table
    minor: Callback failed
#004: ../../../src/H5Gloc.c line 385 in H5G_loc_find_cb(): object 'cell_geometry' doesn't exist
    major: Symbol table
    minor: Object not found
[0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------
[0]PETSC ERROR: Error in external library
[0]PETSC ERROR: Error in HDF5 call H5Dopen2() Status -1
[0]PETSC ERROR: See http://www.mcs.anl.gov/petsc/documentation/faq.html for trouble shooting. [0]PETSC ERROR: Petsc Development GIT revision: v3.6.3-3170-g725a008 GIT Date: 2016-01-19 11:34:55 -0600
[0]PETSC ERROR: readdmvec2

`6` on a linux-gnu-c-opt named des108 by acro018 Fri Feb 12 17:36:00 2016 [0]PETSC ERROR: Configure options --download-netcdf --download-exodusii --with-hdf5-dir=/usr --download-triangle --download-ptscotch --
download-chaco
[0]PETSC ERROR: #1 VecLoad_HDF5() line 269 in /home/acro018/software/PETSc/code/src/vec/vec/utils/vecio.c [0]PETSC ERROR: #2 VecLoad_Default() line 411 in /home/acro018/software/PETSc/code/src/vec/vec/utils/vecio.c [0]PETSC ERROR: #3 VecLoad_Plex_Local() line 214 in /home/acro018/software/PETSc/code/src/dm/impls/plex/plex.c [0]PETSC ERROR: #4 VecLoad_Plex_HDF5() line 211 in /home/acro018/software/PETSc/code/src/dm/impls/plex/plexhdf5.c [0]PETSC ERROR: #5 VecLoad_Plex() line 238 in /home/acro018/software/PETSc/code/src/dm/impls/plex/plex.c [0]PETSC ERROR: #6 VecLoad() line 975 in /home/acro018/software/PETSc/code/src/vec/vec/interface/vector.c

--
Dr Adrian Croucher
Senior Research Fellow
Department of Engineering Science
University of Auckland, New Zealand
email: [email protected]
tel: +64 (0)9 923 84611

program readdmvec

  implicit none

#include <petsc/finclude/petsc.h90>

  PetscInt, parameter :: overlap = 1
  MPI_Comm :: comm
  DM :: dm, dist_dm, ghost_dm
  PetscInt :: dim
  PetscErrorCode :: ierr
  Vec :: cell_geom
  PetscViewer :: viewer
  character(len = 8) :: open_boundary_label_name = "open bdy"
  PetscInt :: num_fields, num_bc, i
  PetscInt, target :: num_components(2), field_dim(2)
  PetscInt, allocatable, target :: num_dof(:)
  PetscInt, pointer :: pnum_dof(:), pnum_components(:)
  PetscInt, target :: bc_field(1)
  PetscInt, pointer :: pbc_field(:)
  IS, target :: bc_comps(1), bc_points(1)
  IS, pointer :: pbc_comps(:), pbc_points(:)
  PetscSection :: section

  call PetscInitialize(PETSC_NULL_CHARACTER, ierr); CHKERRQ(ierr)
  comm = PETSC_COMM_WORLD

  call DMPlexCreateFromFile(comm, "ghotcol.exo", PETSC_TRUE, dm, ierr); CHKERRQ(ierr)
  call DMPlexSetAdjacencyUseCone(dm, PETSC_TRUE, ierr); CHKERRQ(ierr)
  call DMPlexSetAdjacencyUseClosure(dm, PETSC_FALSE, ierr); CHKERRQ(ierr)

  call DMCreateLabel(dm, open_boundary_label_name, ierr); CHKERRQ(ierr)
  call DMSetLabelValue(dm, open_boundary_label_name, 55, 1, ierr); CHKERRQ(ierr)

  call DMPlexDistribute(dm, overlap, PETSC_NULL_OBJECT, dist_dm, ierr)
  CHKERRQ(ierr)
  if (dist_dm /= 0) then
     call DMDestroy(dm, ierr); CHKERRQ(ierr)
     dm = dist_dm
  end if
  call DMPlexConstructGhostCells(dm, open_boundary_label_name, &
         PETSC_NULL_INTEGER, ghost_dm, ierr); CHKERRQ(ierr)
  if (ghost_dm /= 0) then
     call DMDestroy(dm, ierr); CHKERRQ(ierr);
     dm = ghost_dm
  end if

  ! Set up section:
  field_dim = [3, 3]
  num_fields = size(field_dim)
  num_components = [3, 1]; pnum_components => num_components
  call DMGetDimension(dm, dim, ierr); CHKERRQ(ierr)
  allocate(num_dof(num_fields*(dim+1)))
  num_dof = 0
  do i = 1, num_fields
     num_dof((i-1) * (dim+1) + field_dim(i) + 1) = num_components(i)
  end do
  pnum_dof => num_dof
  num_bc = 0
  bc_field(1) = 0; pbc_field => bc_field
  pbc_comps => bc_comps
  pbc_points => bc_points
  call DMPlexCreateSection(dm, dim, num_fields, pnum_components, &
       pnum_dof, num_bc, pbc_field, pbc_comps, pbc_points, &
       PETSC_NULL_OBJECT, section, ierr); CHKERRQ(ierr)
  call DMSetDefaultSection(dm, section, ierr); CHKERRQ(ierr)

  call DMView(dm, PETSC_VIEWER_STDOUT_WORLD, ierr); CHKERRQ(ierr)

  ! Read vector from HDF5:
  call DMGetGlobalVector(dm, cell_geom, ierr); CHKERRQ(ierr)
  call PetscObjectSetName(cell_geom, "cell_geometry", ierr); CHKERRQ(ierr)
  call PetscViewerHDF5Open(comm, "hotcol_ss.h5", FILE_MODE_READ, &
       viewer, ierr); CHKERRQ(ierr)
  call PetscViewerHDF5PushGroup(viewer, "/", ierr); CHKERRQ(ierr)
  call VecLoad(cell_geom, viewer, ierr); CHKERRQ(ierr)
  call PetscViewerHDF5PopGroup(viewer, ierr); CHKERRQ(ierr)
  call PetscViewerDestroy(viewer, ierr); CHKERRQ(ierr)
  call VecView(cell_geom, PETSC_VIEWER_STDOUT_WORLD, ierr); CHKERRQ(ierr)

  call DMRestoreGlobalVector(dm, cell_geom, ierr); CHKERRQ(ierr)
  call DMDestroy(dm, ierr); CHKERRQ(ierr)
  call PetscFinalize(ierr); CHKERRQ(ierr)

end program readdmvec

Attachment: ghotcol.exo
Description: Binary data

Attachment: hotcol_ss.h5
Description: HDF file

Reply via email to