[petsc-dev] hdf5 viewer question
    Dominic Meiser 
    dmeiser at txcorp.com
       
    Tue Jun 23 09:38:06 CDT 2015
    
    
  
Hi,
I'm running into an issue with the Hdf5 viewer. For parallel vectors of 
length larger than 2^16 I get the error message below. The vector is 
created from a 1D DMDA with 2 dofs. Serial vectors work fine as do 
vectors of length <= 2^16. The error message is reproduced below. The 
issue can be reproduced with the attached example with 'mpirun -np 2 
./testcase -da_grid_x 65537'.
Do I need to configure the viewer in some way before using it?
Thanks,
Dominic
HDF5-DIAG: Error detected in HDF5 (1.8.12) MPI-process 0:
   #000: H5D.c line 170 in H5Dcreate2(): unable to create dataset
     major: Dataset
     minor: Unable to initialize object
   #001: H5Dint.c line 439 in H5D__create_named(): unable to create and 
link to dataset
     major: Dataset
     minor: Unable to initialize object
   #002: H5L.c line 1638 in H5L_link_object(): unable to create new link 
to object
     major: Links
     minor: Unable to initialize object
   #003: H5L.c line 1882 in H5L_create_real(): can't insert link
     major: Symbol table
     minor: Unable to insert object
   #004: H5Gtraverse.c line 861 in H5G_traverse(): internal path 
traversal failed
     major: Symbol table
     minor: Object not found
   #005: H5Gtraverse.c line 641 in H5G_traverse_real(): traversal 
operator failed
     major: Symbol table
     minor: Callback failed
   #006: H5L.c line 1685 in H5L_link_cb(): unable to create object
     major: Object header
     minor: Unable to initialize object
   #007: H5O.c line 3015 in H5O_obj_create(): unable to open object
     major: Object header
HDF5-DIAG: Error detected in HDF5 (1.8.12) MPI-process 1:
   #000: H5D.c line 170 in H5Dcreate2(): unable to create dataset
     major: Dataset
     minor: Unable to initialize object
   #001: H5Dint.c line 439 in H5D__create_named(): unable to create and 
link to dataset
     major: Dataset
     minor: Unable to initialize object
   #002: H5L.c line 1638 in H5L_link_object(): unable to create new link 
to object
     major: Links
     minor: Unable to initialize object
   #003: H5L.c line 1882 in H5L_create_real(): can't insert link
     major: Symbol table
     minor: Unable to insert object
   #004: H5Gtraverse.c line 861 in H5G_traverse(): internal path 
traversal failed
     major: Symbol table
     minor: Object not found
   #005: H5Gtraverse.c line 641 in H5G_traverse_real(): traversal 
operator failed
     major: Symbol table
     minor: Callback failed
   #006: H5L.c line 1685 in H5L_link_cb(): unable to create object
     major: Object header
     minor: Unable to initialize object
   #007: H5O.c line 3015 in H5O_obj_create(): unable to open object
     major: Object header
     minor: Can't open object
   #008: H5Doh.c line 293 in H5O__dset_create(): unable to create dataset
     major: Dataset
     minor: Unable to initialize object
   #009: H5Dint.c line 1044 in H5D__create(): unable to construct layout 
information
     major: Dataset
     minor: Unable to initialize object
   #010: H5Dchunk.c line 539 in H5D__chunk_construct(): chunk size must 
be <= maximum dimension size for fixed-sized dimensions
     major: Dataset
     minor: Unable to initialize object
     minor: Can't open object
   #008: H5Doh.c line 293 in H5O__dset_create(): unable to create dataset
     major: Dataset
     minor: Unable to initialize object
   #009: H5Dint.c line 1044 in H5D__create(): unable to construct layout 
information
     major: Dataset
     minor: Unable to initialize object
   #010: H5Dchunk.c line 539 in H5D__chunk_construct(): chunk size must 
be <= maximum dimension size for fixed-sized dimensions
     major: Dataset
     minor: Unable to initialize object
[0]PETSC ERROR: [1]PETSC ERROR: --------------------- Error Message 
--------------------------------------------------------------
--------------------- Error Message 
--------------------------------------------------------------
[0]PETSC ERROR: [1]PETSC ERROR: Error in external library
Error in external library
[0]PETSC ERROR: Error in HDF5 call H5Dcreate2() Status -1
[0]PETSC ERROR: [1]PETSC ERROR: Error in HDF5 call H5Dcreate2() Status -1
[1]PETSC ERROR: See http://www.mcs.anl.gov/petsc/documentation/faq.html 
for trouble shooting.
See http://www.mcs.anl.gov/petsc/documentation/faq.html for trouble 
shooting.
[0]PETSC ERROR: Petsc Release Version 3.6.0, unknown
[0]PETSC ERROR: [1]PETSC ERROR: Petsc Release Version 3.6.0, unknown
[1]PETSC ERROR: ./testcase on a fftw_test named longspeak.txcorp.com by 
dmeiser Mon Jun 22 16:59:21 2015
./testcase on a fftw_test named longspeak.txcorp.com by dmeiser Mon Jun 
22 16:59:21 2015
[0]PETSC ERROR: Configure options --download-hdf5=1 --download-fftw=1 
--with-mpi=1 PETSC_ARCH=fftw_test
[0]PETSC ERROR: [1]PETSC ERROR: Configure options --download-hdf5=1 
--download-fftw=1 --with-mpi=1 PETSC_ARCH=fftw_test
[1]PETSC ERROR: #1 VecView_MPI_HDF5_DA() line 535 in 
/home/scratch/petsc/src/dm/impls/da/gr2.c
#1 VecView_MPI_HDF5_DA() line 535 in 
/home/scratch/petsc/src/dm/impls/da/gr2.c
[0]PETSC ERROR: #2 VecView_MPI_DA() line 713 in 
/home/scratch/petsc/src/dm/impls/da/gr2.c
[1]PETSC ERROR: #2 VecView_MPI_DA() line 713 in 
/home/scratch/petsc/src/dm/impls/da/gr2.c
[0]PETSC ERROR: #3 VecView() line 618 in 
/home/scratch/petsc/src/vec/vec/interface/vector.c
[1]PETSC ERROR: #3 VecView() line 618 in 
/home/scratch/petsc/src/vec/vec/interface/vector.c
[0]PETSC ERROR: [1]PETSC ERROR: #4 dump() line 75 in testcase.c
#4 dump() line 75 in testcase.c
[0]PETSC ERROR: #5 main() line 31 in testcase.c
[1]PETSC ERROR: #5 main() line 31 in testcase.c
[0]PETSC ERROR: [1]PETSC ERROR: PETSc Option Table entries:
[1]PETSC ERROR: PETSc Option Table entries:
[0]PETSC ERROR: -da_grid_x 65537
[0]PETSC ERROR: -da_grid_x 65537
[1]PETSC ERROR: ----------------End of Error Message -------send entire 
error message to petsc-maint at mcs.anl.gov----------
----------------End of Error Message -------send entire error message to 
petsc-maint at mcs.anl.gov----------
[cli_0]: aborting job:
application called MPI_Abort(MPI_COMM_WORLD, 76) - process 0
[cli_1]: aborting job:
application called MPI_Abort(MPI_COMM_WORLD, 76) - process 1
===================================================================================
=   BAD TERMINATION OF ONE OF YOUR APPLICATION PROCESSES
=   PID 15847 RUNNING AT longspeak.txcorp.com
=   EXIT CODE: 76
=   CLEANING UP REMAINING PROCESSES
=   YOU CAN IGNORE THE BELOW CLEANUP MESSAGES
===================================================================================
-- 
Dominic Meiser
Tech-X Corporation
5621 Arapahoe Avenue
Boulder, CO 80303
USA
Telephone: 303-996-2036
Fax: 303-448-7756
www.txcorp.com
-------------- next part --------------
#include <petscvec.h>
#include <petscdmda.h>
#include <petscdm.h>
#include <petscviewer.h>
#include <petscviewerhdf5.h>
struct Field {
  PetscScalar E;
  PetscScalar H;
};
static const char *fieldNames[] = {"E", "H"};
PetscErrorCode setInitialState(Vec *field);
PetscErrorCode dump(Vec field, const char *fieldname, const char *filename);
int main(int argn, char **argv) {
  PetscErrorCode ierr;
  DM             da;
  Vec            field;
  PetscFunctionBegin;
  ierr = PetscInitialize(&argn, &argv, NULL, NULL);CHKERRQ(ierr);
  ierr = DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, -1000, 2, 1, NULL, &da);CHKERRQ(ierr);
  ierr = DMDASetUniformCoordinates(da, -1.0, 1.0, 0, 0, 0, 0);CHKERRQ(ierr);
  ierr = DMDASetCoordinateName(da, 0, "x");CHKERRQ(ierr);
  ierr = DMDASetFieldNames(da, fieldNames);CHKERRQ(ierr);
  ierr = DMSetFromOptions(da);CHKERRQ(ierr);
  ierr = DMCreateLocalVector(da, &field);CHKERRQ(ierr);
  ierr = PetscObjectSetName((PetscObject)field, "field");CHKERRQ(ierr);
  ierr = setInitialState(&field);CHKERRQ(ierr);
  ierr = dump(field, "field", "field_initial.h5");CHKERRQ(ierr);
  PetscFunctionReturn(0);
}
PetscErrorCode setInitialState(Vec *field) {
  PetscErrorCode ierr;
  struct Field   *f;
  DM             da;
  PetscInt       i, x, m;
  PetscFunctionBegin;
  ierr = VecGetDM(*field, &da);CHKERRQ(ierr);
  ierr = DMDAVecGetArray(da, *field, &f);CHKERRQ(ierr);
  ierr = DMDAGetCorners(da, &x, NULL, NULL, &m, NULL, NULL);CHKERRQ(ierr);
  for (i = x; i < x + m; ++i) {
    f[i].E = 0.0;
    f[i].H = 0.0;
  }
  ierr = DMDAVecRestoreArray(da, *field, &f);CHKERRQ(ierr);
  DMLocalToLocalBegin(da, *field, INSERT_VALUES, *field);CHKERRQ(ierr);
  DMLocalToLocalEnd(da, *field, INSERT_VALUES, *field);CHKERRQ(ierr);
  PetscFunctionReturn(0);
}
PetscErrorCode dump(Vec field, const char *fieldname, const char *filename) {
  PetscErrorCode ierr;
  PetscViewer    v;
  DM             dm;
  Vec            tmp;
  PetscInt       numProcs;
  PetscFunctionBegin;
  ierr = VecGetDM(field, &dm);CHKERRQ(ierr);
  ierr = DMDAGetInfo(dm, 0, 0, 0, 0, &numProcs, 0, 0, 0, 0, 0, 0, 0, 0);
  if (numProcs > 1) {
    ierr = DMGetGlobalVector(dm, &tmp);CHKERRQ(ierr);
    ierr = DMLocalToGlobalBegin(dm, field, INSERT_VALUES, tmp);CHKERRQ(ierr);
    ierr = DMLocalToGlobalEnd(dm, field, INSERT_VALUES, tmp);CHKERRQ(ierr);
  } else {
    tmp = field;
  }
  ierr = PetscObjectSetName((PetscObject)tmp, fieldname);CHKERRQ(ierr);
  ierr = PetscViewerHDF5Open(MPI_COMM_WORLD, filename, FILE_MODE_WRITE, &v);CHKERRQ(ierr);
  ierr = VecView(tmp, v);CHKERRQ(ierr);
  ierr = PetscViewerDestroy(&v);CHKERRQ(ierr);
  PetscFunctionReturn(0);
}
    
    
More information about the petsc-dev
mailing list