Actual source code: hdf5v.c

  1: #include <petsc/private/viewerimpl.h>
  2: #include <petsc/private/viewerhdf5impl.h>
  3: #include <petscviewerhdf5.h>

  5: static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer, const char[], PetscBool, PetscBool*, H5O_type_t*);
  6: static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer, const char[], const char[], PetscBool*);

  8: static PetscErrorCode PetscViewerHDF5GetAbsolutePath_Internal(PetscViewer viewer, const char objname[], char **fullpath)
  9: {
 10:   const char *group;
 11:   char buf[PETSC_MAX_PATH_LEN]="";

 15:   PetscViewerHDF5GetGroup(viewer, &group);
 16:   PetscStrcat(buf, group);
 17:   PetscStrcat(buf, "/");
 18:   PetscStrcat(buf, objname);
 19:   PetscStrallocpy(buf, fullpath);
 20:   return(0);
 21: }

 23: static PetscErrorCode PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer, PetscObject obj)
 24: {
 25:   PetscBool has;
 26:   const char *group;

 30:   if (!obj->name) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named");
 31:   PetscViewerHDF5HasObject(viewer, obj, &has);
 32:   if (!has) {
 33:     PetscViewerHDF5GetGroup(viewer, &group);
 34:     SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) %s not stored in group %s", obj->name, group);
 35:   }
 36:   return(0);
 37: }

 39: static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptionItems *PetscOptionsObject,PetscViewer v)
 40: {
 41:   PetscErrorCode   ierr;
 42:   PetscBool        flg = PETSC_FALSE, set;
 43:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data;

 46:   PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");
 47:   PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);
 48:   PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL);
 49:   PetscOptionsBool("-viewer_hdf5_collective","Enable collective transfer mode","PetscViewerHDF5SetCollective",flg,&flg,&set);
 50:   if (set) {PetscViewerHDF5SetCollective(v,flg);}
 51:   PetscOptionsTail();
 52:   return(0);
 53: }

 55: static PetscErrorCode PetscViewerView_HDF5(PetscViewer v,PetscViewer viewer)
 56: {
 57:   PetscViewer_HDF5  *hdf5 = (PetscViewer_HDF5*)v->data;
 58:   PetscBool         flg;
 59:   PetscErrorCode    ierr;

 62:   if (hdf5->filename) {
 63:     PetscViewerASCIIPrintf(viewer,"Filename: %s\n",hdf5->filename);
 64:   }
 65:   PetscViewerASCIIPrintf(viewer,"Vectors with blocksize 1 saved as 2D datasets: %s\n",PetscBools[hdf5->basedimension2]);
 66:   PetscViewerASCIIPrintf(viewer,"Enforce single precision storage: %s\n",PetscBools[hdf5->spoutput]);
 67:   PetscViewerHDF5GetCollective(v,&flg);
 68:   PetscViewerASCIIPrintf(viewer,"MPI-IO transfer mode: %s\n",flg ? "collective" : "independent");
 69:   return(0);
 70: }

 72: static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer)
 73: {
 74:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data;
 75:   PetscErrorCode   ierr;

 78:   PetscFree(hdf5->filename);
 79:   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
 80:   return(0);
 81: }

 83: static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer)
 84: {
 85:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
 86:   PetscErrorCode   ierr;

 89:   PetscStackCallHDF5(H5Pclose,(hdf5->dxpl_id));
 90:   PetscViewerFileClose_HDF5(viewer);
 91:   while (hdf5->groups) {
 92:     PetscViewerHDF5GroupList *tmp = hdf5->groups->next;

 94:     PetscFree(hdf5->groups->name);
 95:     PetscFree(hdf5->groups);
 96:     hdf5->groups = tmp;
 97:   }
 98:   PetscFree(hdf5);
 99:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);
100:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);
101:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);
102:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetBaseDimension2_C",NULL);
103:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetSPOutput_C",NULL);
104:   return(0);
105: }

107: static PetscErrorCode  PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type)
108: {
109:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

112:   hdf5->btype = type;
113:   return(0);
114: }

116: static PetscErrorCode  PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type)
117: {
118:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

121:   *type = hdf5->btype;
122:   return(0);
123: }

125: static PetscErrorCode  PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg)
126: {
127:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

130:   hdf5->basedimension2 = flg;
131:   return(0);
132: }

134: /*@
135:      PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
136:        dimension of 2.

138:     Logically Collective on PetscViewer

140:   Input Parameters:
141: +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
142: -  flg - if PETSC_TRUE the vector will always have at least a dimension of 2 even if that first dimension is of size 1

144:   Options Database:
145: .  -viewer_hdf5_base_dimension2 - turns on (true) or off (false) using a dimension of 2 in the HDF5 file even if the bs/dof of the vector is 1


148:   Notes:
149:     Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof
150:          of one when the dimension is lower. Others think the option is crazy.

152:   Level: intermediate

154: .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()

156: @*/
157: PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg)
158: {

163:   PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));
164:   return(0);
165: }

167: /*@
168:      PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
169:        dimension of 2.

171:     Logically Collective on PetscViewer

173:   Input Parameter:
174: .  viewer - the PetscViewer, must be of type HDF5

176:   Output Parameter:
177: .  flg - if PETSC_TRUE the vector will always have at least a dimension of 2 even if that first dimension is of size 1

179:   Notes:
180:     Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof
181:          of one when the dimension is lower. Others think the option is crazy.

183:   Level: intermediate

185: .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()

187: @*/
188: PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg)
189: {
190:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

194:   *flg = hdf5->basedimension2;
195:   return(0);
196: }

198: static PetscErrorCode  PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg)
199: {
200:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

203:   hdf5->spoutput = flg;
204:   return(0);
205: }

207: /*@
208:      PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is
209:        compiled with double precision PetscReal.

211:     Logically Collective on PetscViewer

213:   Input Parameters:
214: +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
215: -  flg - if PETSC_TRUE the data will be written to disk with single precision

217:   Options Database:
218: .  -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision


221:   Notes:
222:     Setting this option does not make any difference if PETSc is compiled with single precision
223:          in the first place. It does not affect reading datasets (HDF5 handle this internally).

225:   Level: intermediate

227: .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
228:           PetscReal

230: @*/
231: PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg)
232: {

237:   PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));
238:   return(0);
239: }

241: /*@
242:      PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is
243:        compiled with double precision PetscReal.

245:     Logically Collective on PetscViewer

247:   Input Parameter:
248: .  viewer - the PetscViewer, must be of type HDF5

250:   Output Parameter:
251: .  flg - if PETSC_TRUE the data will be written to disk with single precision

253:   Notes:
254:     Setting this option does not make any difference if PETSc is compiled with single precision
255:          in the first place. It does not affect reading datasets (HDF5 handle this internally).

257:   Level: intermediate

259: .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
260:           PetscReal

262: @*/
263: PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg)
264: {
265:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

269:   *flg = hdf5->spoutput;
270:   return(0);
271: }

273: static PetscErrorCode  PetscViewerHDF5SetCollective_HDF5(PetscViewer viewer, PetscBool flg)
274: {
276:   /* H5FD_MPIO_COLLECTIVE is wrong in hdf5 1.10.2, and is the same as H5FD_MPIO_INDEPENDENT in earlier versions
277:      - see e.g. https://gitlab.cosma.dur.ac.uk/swift/swiftsim/issues/431 */
278: #if H5_VERSION_GE(1,10,3) && defined(H5_HAVE_PARALLEL)
279:   {
280:     PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
281:     PetscStackCallHDF5(H5Pset_dxpl_mpio,(hdf5->dxpl_id, flg ? H5FD_MPIO_COLLECTIVE : H5FD_MPIO_INDEPENDENT));
282:   }
283: #else
284:   if (flg) {
286:     PetscPrintf(PetscObjectComm((PetscObject)viewer), "Warning: PetscViewerHDF5SetCollective(viewer,PETSC_TRUE) is ignored for HDF5 versions prior to 1.10.3 or if built without MPI support\n");
287:   }
288: #endif
289:   return(0);
290: }

292: /*@
293:   PetscViewerHDF5SetCollective - Use collective MPI-IO transfer mode for HDF5 reads and writes.

295:   Logically Collective; flg must contain common value

297:   Input Parameters:
298: + viewer - the PetscViewer; if it is not hdf5 then this command is ignored
299: - flg - PETSC_TRUE for collective mode; PETSC_FALSE for independent mode (default)

301:   Options Database:
302: . -viewer_hdf5_collective - turns on (true) or off (false) collective transfers

304:   Notes:
305:   Collective mode gives the MPI-IO layer underneath HDF5 a chance to do some additional collective optimizations and hence can perform better.
306:   However, this works correctly only since HDF5 1.10.3 and if HDF5 is installed for MPI; hence, we ignore this setting for older versions.

308:   Developer notes:
309:   In the HDF5 layer, PETSC_TRUE / PETSC_FALSE means H5Pset_dxpl_mpio() is called with H5FD_MPIO_COLLECTIVE / H5FD_MPIO_INDEPENDENT, respectively.
310:   This in turn means use of MPI_File_{read,write}_all /  MPI_File_{read,write} in the MPI-IO layer, respectively.
311:   See HDF5 documentation and MPI-IO documentation for details.

313:   Level: intermediate

315: .seealso: PetscViewerHDF5GetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open()

317: @*/
318: PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer viewer,PetscBool flg)
319: {

325:   PetscTryMethod(viewer,"PetscViewerHDF5SetCollective_C",(PetscViewer,PetscBool),(viewer,flg));
326:   return(0);
327: }

329: static PetscErrorCode  PetscViewerHDF5GetCollective_HDF5(PetscViewer viewer, PetscBool *flg)
330: {
331: #if defined(H5_HAVE_PARALLEL)
332:   PetscViewer_HDF5  *hdf5 = (PetscViewer_HDF5*) viewer->data;
333:   H5FD_mpio_xfer_t  mode;
334: #endif

337: #if !defined(H5_HAVE_PARALLEL)
338:   *flg = PETSC_FALSE;
339: #else
340:   PetscStackCallHDF5(H5Pget_dxpl_mpio,(hdf5->dxpl_id, &mode));
341:   *flg = (mode == H5FD_MPIO_COLLECTIVE) ? PETSC_TRUE : PETSC_FALSE;
342: #endif
343:   return(0);
344: }

346: /*@
347:   PetscViewerHDF5GetCollective - Return flag whether collective MPI-IO transfer mode is used for HDF5 reads and writes.

349:   Not Collective

351:   Input Parameters:
352: . viewer - the HDF5 PetscViewer

354:   Output Parameters:
355: . flg - the flag

357:   Level: intermediate

359:   Notes:
360:   This setting works correctly only since HDF5 1.10.3 and if HDF5 was installed for MPI. For older versions, PETSC_FALSE will be always returned.
361:   For more details, see PetscViewerHDF5SetCollective().

363: .seealso: PetscViewerHDF5SetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open()

365: @*/
366: PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer viewer,PetscBool *flg)
367: {


374:   PetscUseMethod(viewer,"PetscViewerHDF5GetCollective_C",(PetscViewer,PetscBool*),(viewer,flg));
375:   return(0);
376: }

378: static PetscErrorCode  PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[])
379: {
380:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
381:   hid_t             plist_id;
382:   PetscErrorCode    ierr;

385:   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
386:   if (hdf5->filename) {PetscFree(hdf5->filename);}
387:   PetscStrallocpy(name, &hdf5->filename);
388:   /* Set up file access property list with parallel I/O access */
389:   PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS));
390: #if defined(H5_HAVE_PARALLEL)
391:   PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), MPI_INFO_NULL));
392: #endif
393:   /* Create or open the file collectively */
394:   switch (hdf5->btype) {
395:   case FILE_MODE_READ:
396:     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id));
397:     break;
398:   case FILE_MODE_APPEND:
399:   case FILE_MODE_UPDATE:
400:   {
401:     PetscBool flg;
402:     PetscTestFile(hdf5->filename, 'r', &flg);
403:     if (flg) PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id));
404:     else     PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_EXCL, H5P_DEFAULT, plist_id));
405:     break;
406:   }
407:   case FILE_MODE_WRITE:
408:     PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id));
409:     break;
410:   case FILE_MODE_UNDEFINED:
411:     SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
412:   default:
413:     SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP, "Unsupported file mode %s",PetscFileModes[hdf5->btype]);
414:   }
415:   if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name);
416:   PetscStackCallHDF5(H5Pclose,(plist_id));
417:   return(0);
418: }

420: static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name)
421: {
422:   PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data;

425:   *name = vhdf5->filename;
426:   return(0);
427: }

429: static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer)
430: {
431:   /*
432:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
433:   PetscErrorCode   ierr;
434:   */

437:   return(0);
438: }

440: /*MC
441:    PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file


444: .seealso:  PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET,
445:            PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING,
446:            PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB,
447:            PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType()

449:   Level: beginner
450: M*/

452: PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v)
453: {
454:   PetscViewer_HDF5 *hdf5;
455:   PetscErrorCode   ierr;

458: #if !defined(H5_HAVE_PARALLEL)
459:   {
460:     PetscMPIInt size;
461:     MPI_Comm_size(PetscObjectComm((PetscObject)v), &size);
462:     if (size > 1) SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_SUP, "Cannot use parallel HDF5 viewer since the given HDF5 does not support parallel I/O (H5_HAVE_PARALLEL is unset)");
463:   }
464: #endif

466:   PetscNewLog(v,&hdf5);

468:   v->data                = (void*) hdf5;
469:   v->ops->destroy        = PetscViewerDestroy_HDF5;
470:   v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5;
471:   v->ops->setup          = PetscViewerSetUp_HDF5;
472:   v->ops->view           = PetscViewerView_HDF5;
473:   v->ops->flush          = NULL;
474:   hdf5->btype            = FILE_MODE_UNDEFINED;
475:   hdf5->filename         = NULL;
476:   hdf5->timestep         = -1;
477:   hdf5->groups           = NULL;

479:   PetscStackCallHDF5Return(hdf5->dxpl_id,H5Pcreate,(H5P_DATASET_XFER));

481:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);
482:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);
483:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);
484:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_HDF5);
485:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);
486:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);
487:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetCollective_C",PetscViewerHDF5SetCollective_HDF5);
488:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetCollective_C",PetscViewerHDF5GetCollective_HDF5);
489:   return(0);
490: }

492: /*@C
493:    PetscViewerHDF5Open - Opens a file for HDF5 input/output.

495:    Collective

497:    Input Parameters:
498: +  comm - MPI communicator
499: .  name - name of file
500: -  type - type of file

502:    Output Parameter:
503: .  hdf5v - PetscViewer for HDF5 input/output to use with the specified file

505:   Options Database:
506: +  -viewer_hdf5_base_dimension2 - turns on (true) or off (false) using a dimension of 2 in the HDF5 file even if the bs/dof of the vector is 1
507: -  -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal

509:    Level: beginner

511:    Notes:
512:    Reading is always available, regardless of the mode. Available modes are
513: +  FILE_MODE_READ - open existing HDF5 file for read only access, fail if file does not exist [H5Fopen() with H5F_ACC_RDONLY]
514: .  FILE_MODE_WRITE - if file exists, fully overwrite it, else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_TRUNC]
515: .  FILE_MODE_APPEND - if file exists, keep existing contents [H5Fopen() with H5F_ACC_RDWR], else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_EXCL]
516: -  FILE_MODE_UPDATE - same as FILE_MODE_APPEND

518:    In case of FILE_MODE_APPEND / FILE_MODE_UPDATE, any stored object (dataset, attribute) can be selectively ovewritten if the same fully qualified name (/group/path/to/object) is specified.

520:    This PetscViewer should be destroyed with PetscViewerDestroy().


523: .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(),
524:           PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(),
525:           MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName()
526: @*/
527: PetscErrorCode  PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
528: {

532:   PetscViewerCreate(comm, hdf5v);
533:   PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);
534:   PetscViewerFileSetMode(*hdf5v, type);
535:   PetscViewerFileSetName(*hdf5v, name);
536:   PetscViewerSetFromOptions(*hdf5v);
537:   return(0);
538: }

540: /*@C
541:   PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls

543:   Not collective

545:   Input Parameter:
546: . viewer - the PetscViewer

548:   Output Parameter:
549: . file_id - The file id

551:   Level: intermediate

553: .seealso: PetscViewerHDF5Open()
554: @*/
555: PetscErrorCode  PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
556: {
557:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

561:   if (file_id) *file_id = hdf5->file_id;
562:   return(0);
563: }

565: /*@C
566:   PetscViewerHDF5PushGroup - Set the current HDF5 group for output

568:   Not collective

570:   Input Parameters:
571: + viewer - the PetscViewer
572: - name - The group name

574:   Level: intermediate

576:   Note: The group name being NULL, empty string, or a sequence of all slashes (e.g. "///") is always internally stored as NULL and interpreted as "/".

578: .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
579: @*/
580: PetscErrorCode  PetscViewerHDF5PushGroup(PetscViewer viewer, const char name[])
581: {
582:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
583:   PetscViewerHDF5GroupList *groupNode;
584:   PetscErrorCode   ierr;

589:   if (name && name[0]) {
590:      size_t i,len;
591:      PetscStrlen(name, &len);
592:      for (i=0; i<len; i++) if (name[i] != '/') break;
593:      if (i == len) name = NULL;
594:   } else name = NULL;
595:   PetscNew(&groupNode);
596:   PetscStrallocpy(name, (char**) &groupNode->name);
597:   groupNode->next = hdf5->groups;
598:   hdf5->groups    = groupNode;
599:   return(0);
600: }

602: /*@
603:   PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value

605:   Not collective

607:   Input Parameter:
608: . viewer - the PetscViewer

610:   Level: intermediate

612: .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
613: @*/
614: PetscErrorCode  PetscViewerHDF5PopGroup(PetscViewer viewer)
615: {
616:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
617:   PetscViewerHDF5GroupList *groupNode;
618:   PetscErrorCode   ierr;

622:   if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
623:   groupNode    = hdf5->groups;
624:   hdf5->groups = hdf5->groups->next;
625:   PetscFree(groupNode->name);
626:   PetscFree(groupNode);
627:   return(0);
628: }

630: /*@C
631:   PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup().
632:   If none has been assigned, returns NULL.

634:   Not collective

636:   Input Parameter:
637: . viewer - the PetscViewer

639:   Output Parameter:
640: . name - The group name

642:   Level: intermediate

644: .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup()
645: @*/
646: PetscErrorCode  PetscViewerHDF5GetGroup(PetscViewer viewer, const char *name[])
647: {
648:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data;

653:   if (hdf5->groups) *name = hdf5->groups->name;
654:   else *name = NULL;
655:   return(0);
656: }

658: /*@
659:   PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(),
660:   and return this group's ID and file ID.
661:   If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID.

663:   Not collective

665:   Input Parameter:
666: . viewer - the PetscViewer

668:   Output Parameter:
669: + fileId - The HDF5 file ID
670: - groupId - The HDF5 group ID

672:   Notes:
673:   If the viewer is writable, the group is created if it doesn't exist yet.

675:   Level: intermediate

677: .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
678: @*/
679: PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId)
680: {
681:   hid_t          file_id;
682:   H5O_type_t     type;
683:   const char     *groupName = NULL, *fileName = NULL;
684:   PetscBool      writable, has;

688:   PetscViewerWritable(viewer, &writable);
689:   PetscViewerHDF5GetFileId(viewer, &file_id);
690:   PetscViewerFileGetName(viewer, &fileName);
691:   PetscViewerHDF5GetGroup(viewer, &groupName);
692:   PetscViewerHDF5Traverse_Internal(viewer, groupName, writable, &has, &type);
693:   if (!has) {
694:     if (!writable) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Group %s does not exist and file %s is not open for writing", groupName, fileName);
695:     else           SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_LIB, "HDF5 failed to create group %s although file %s is open for writing", groupName, fileName);
696:   }
697:   if (type != H5O_TYPE_GROUP) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Path %s in file %s resolves to something which is not a group", groupName, fileName);
698:   PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT));
699:   *fileId  = file_id;
700:   return(0);
701: }

703: /*@
704:   PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time.

706:   Not collective

708:   Input Parameter:
709: . viewer - the PetscViewer

711:   Level: intermediate

713: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep()
714: @*/
715: PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
716: {
717:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

721:   ++hdf5->timestep;
722:   return(0);
723: }

725: /*@
726:   PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep
727:   of -1 disables blocking with timesteps.

729:   Not collective

731:   Input Parameters:
732: + viewer - the PetscViewer
733: - timestep - The timestep number

735:   Level: intermediate

737: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
738: @*/
739: PetscErrorCode  PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
740: {
741:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

745:   hdf5->timestep = timestep;
746:   return(0);
747: }

749: /*@
750:   PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.

752:   Not collective

754:   Input Parameter:
755: . viewer - the PetscViewer

757:   Output Parameter:
758: . timestep - The timestep number

760:   Level: intermediate

762: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep()
763: @*/
764: PetscErrorCode  PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
765: {
766:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

771:   *timestep = hdf5->timestep;
772:   return(0);
773: }

775: /*@C
776:   PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name.

778:   Not collective

780:   Input Parameter:
781: . ptype - the PETSc datatype name (for example PETSC_DOUBLE)

783:   Output Parameter:
784: . mtype - the MPI datatype (for example MPI_DOUBLE, ...)

786:   Level: advanced

788: .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
789: @*/
790: PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
791: {
793:   if (ptype == PETSC_INT)
794: #if defined(PETSC_USE_64BIT_INDICES)
795:                                        *htype = H5T_NATIVE_LLONG;
796: #else
797:                                        *htype = H5T_NATIVE_INT;
798: #endif
799:   else if (ptype == PETSC_DOUBLE)      *htype = H5T_NATIVE_DOUBLE;
800:   else if (ptype == PETSC_LONG)        *htype = H5T_NATIVE_LONG;
801:   else if (ptype == PETSC_SHORT)       *htype = H5T_NATIVE_SHORT;
802:   else if (ptype == PETSC_ENUM)        *htype = H5T_NATIVE_INT;
803:   else if (ptype == PETSC_BOOL)        *htype = H5T_NATIVE_INT;
804:   else if (ptype == PETSC_FLOAT)       *htype = H5T_NATIVE_FLOAT;
805:   else if (ptype == PETSC_CHAR)        *htype = H5T_NATIVE_CHAR;
806:   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
807:   else if (ptype == PETSC_STRING)      *htype = H5Tcopy(H5T_C_S1);
808:   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
809:   return(0);
810: }

812: /*@C
813:   PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name

815:   Not collective

817:   Input Parameter:
818: . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...)

820:   Output Parameter:
821: . ptype - the PETSc datatype name (for example PETSC_DOUBLE)

823:   Level: advanced

825: .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
826: @*/
827: PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
828: {
830: #if defined(PETSC_USE_64BIT_INDICES)
831:   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_LONG;
832:   else if (htype == H5T_NATIVE_LLONG)  *ptype = PETSC_INT;
833: #else
834:   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_INT;
835: #endif
836:   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
837:   else if (htype == H5T_NATIVE_LONG)   *ptype = PETSC_LONG;
838:   else if (htype == H5T_NATIVE_SHORT)  *ptype = PETSC_SHORT;
839:   else if (htype == H5T_NATIVE_FLOAT)  *ptype = PETSC_FLOAT;
840:   else if (htype == H5T_NATIVE_CHAR)   *ptype = PETSC_CHAR;
841:   else if (htype == H5T_NATIVE_UCHAR)  *ptype = PETSC_CHAR;
842:   else if (htype == H5T_C_S1)          *ptype = PETSC_STRING;
843:   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
844:   return(0);
845: }

847: /*@C
848:  PetscViewerHDF5WriteAttribute - Write an attribute

850:   Input Parameters:
851: + viewer - The HDF5 viewer
852: . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
853: . name   - The attribute name
854: . datatype - The attribute type
855: - value    - The attribute value

857:   Level: advanced

859: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
860: @*/
861: PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, const void *value)
862: {
863:   char           *parent;
864:   hid_t          h5, dataspace, obj, attribute, dtype;
865:   PetscBool      has;

873:   PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);
874:   PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_TRUE, NULL, NULL);
875:   PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);
876:   PetscDataTypeToHDF5DataType(datatype, &dtype);
877:   if (datatype == PETSC_STRING) {
878:     size_t len;
879:     PetscStrlen((const char *) value, &len);
880:     PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
881:   }
882:   PetscViewerHDF5GetFileId(viewer, &h5);
883:   PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR));
884:   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
885:   if (has) {
886:     PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
887:   } else {
888:     PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
889:   }
890:   PetscStackCallHDF5(H5Awrite,(attribute, dtype, value));
891:   if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype));
892:   PetscStackCallHDF5(H5Aclose,(attribute));
893:   PetscStackCallHDF5(H5Oclose,(obj));
894:   PetscStackCallHDF5(H5Sclose,(dataspace));
895:   PetscFree(parent);
896:   return(0);
897: }

899: /*@C
900:  PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name

902:   Input Parameters:
903: + viewer   - The HDF5 viewer
904: . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
905: . name     - The attribute name
906: . datatype - The attribute type
907: - value    - The attribute value

909:   Notes:
910:   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
911:   You might want to check first if it does using PetscViewerHDF5HasObject().

913:   Level: advanced

915: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
916: @*/
917: PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value)
918: {

926:   PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);
927:   PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);
928:   return(0);
929: }

931: /*@C
932:  PetscViewerHDF5ReadAttribute - Read an attribute

934:   Input Parameters:
935: + viewer - The HDF5 viewer
936: . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
937: . name   - The attribute name
938: - datatype - The attribute type

940:   Output Parameter:
941: . value    - The attribute value

943:   Notes: If the datatype is PETSC_STRING one must PetscFree() the obtained value when it is no longer needed.

945:   Level: advanced

947: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
948: @*/
949: PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, void *value)
950: {
951:   char           *parent;
952:   hid_t          h5, obj, attribute, atype, dtype;
953:   PetscBool      has;

961:   PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);
962:   PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, &has, NULL);
963:   if (has) {PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);}
964:   if (!has) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist", parent, name);
965:   PetscDataTypeToHDF5DataType(datatype, &dtype);
966:   PetscViewerHDF5GetFileId(viewer, &h5);
967:   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
968:   PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
969:   if (datatype == PETSC_STRING) {
970:     size_t len;
971:     PetscStackCallHDF5Return(atype,H5Aget_type,(attribute));
972:     PetscStackCall("H5Tget_size",len = H5Tget_size(atype));
973:     PetscMalloc((len+1) * sizeof(char), value);
974:     PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
975:     PetscStackCallHDF5(H5Aread,(attribute, dtype, *(char**)value));
976:   } else {
977:     PetscStackCallHDF5(H5Aread,(attribute, dtype, value));
978:   }
979:   PetscStackCallHDF5(H5Aclose,(attribute));
980:   /* H5Oclose can be used to close groups, datasets, or committed datatypes */
981:   PetscStackCallHDF5(H5Oclose,(obj));
982:   PetscFree(parent);
983:   return(0);
984: }

986: /*@C
987:  PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name

989:   Input Parameters:
990: + viewer   - The HDF5 viewer
991: . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
992: . name     - The attribute name
993: - datatype - The attribute type

995:   Output Parameter:
996: . value    - The attribute value

998:   Notes:
999:   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1000:   You might want to check first if it does using PetscViewerHDF5HasObject().

1002:   Level: advanced

1004: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute() PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1005: @*/
1006: PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *value)
1007: {

1015:   PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);
1016:   PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, value);
1017:   return(0);
1018: }

1020: PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_)
1021: {
1022:   htri_t exists;
1023:   hid_t group;

1026:   PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT));
1027:   if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT));
1028:   if (!exists && createGroup) {
1029:     PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT));
1030:     PetscStackCallHDF5(H5Gclose,(group));
1031:     exists = PETSC_TRUE;
1032:   }
1033:   *exists_ = (PetscBool) exists;
1034:   return(0);
1035: }

1037: static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype)
1038: {
1039:   const char     rootGroupName[] = "/";
1040:   hid_t          h5;
1041:   PetscBool      exists=PETSC_FALSE;
1042:   PetscInt       i;
1043:   int            n;
1044:   char           **hierarchy;
1045:   char           buf[PETSC_MAX_PATH_LEN]="";

1051:   else name = rootGroupName;
1052:   if (has) {
1054:     *has = PETSC_FALSE;
1055:   }
1056:   if (otype) {
1058:     *otype = H5O_TYPE_UNKNOWN;
1059:   }
1060:   PetscViewerHDF5GetFileId(viewer, &h5);

1062:   /*
1063:      Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing.
1064:      Hence, each of them needs to be tested separately:
1065:      1) whether it's a valid link
1066:      2) whether this link resolves to an object
1067:      See H5Oexists_by_name() documentation.
1068:   */
1069:   PetscStrToArray(name,'/',&n,&hierarchy);
1070:   if (!n) {
1071:     /*  Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */
1072:     if (has)   *has   = PETSC_TRUE;
1073:     if (otype) *otype = H5O_TYPE_GROUP;
1074:     PetscStrToArrayDestroy(n,hierarchy);
1075:     return(0);
1076:   }
1077:   for (i=0; i<n; i++) {
1078:     PetscStrcat(buf,"/");
1079:     PetscStrcat(buf,hierarchy[i]);
1080:     PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);
1081:     if (!exists) break;
1082:   }
1083:   PetscStrToArrayDestroy(n,hierarchy);

1085:   /* If the object exists, get its type */
1086:   if (exists && otype) {
1087:     H5O_info_t info;

1089:     /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */
1090:     PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT));
1091:     *otype = info.type;
1092:   }
1093:   if (has) *has = exists;
1094:   return(0);
1095: }

1097: /*@
1098:  PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file

1100:   Input Parameters:
1101: . viewer - The HDF5 viewer

1103:   Output Parameter:
1104: . has    - Flag for group existence

1106:   Notes:
1107:   If the path exists but is not a group, this returns PETSC_FALSE as well.

1109:   Level: advanced

1111: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5OpenGroup()
1112: @*/
1113: PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, PetscBool *has)
1114: {
1115:   H5O_type_t type;
1116:   const char *name;

1122:   PetscViewerHDF5GetGroup(viewer, &name);
1123:   PetscViewerHDF5Traverse_Internal(viewer, name, PETSC_FALSE, has, &type);
1124:   *has = (type == H5O_TYPE_GROUP) ? PETSC_TRUE : PETSC_FALSE;
1125:   return(0);
1126: }

1128: /*@
1129:  PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group

1131:   Input Parameters:
1132: + viewer - The HDF5 viewer
1133: - obj    - The named object

1135:   Output Parameter:
1136: . has    - Flag for dataset existence; PETSC_FALSE for unnamed object

1138:   Notes:
1139:   If the path exists but is not a dataset, this returns PETSC_FALSE as well.

1141:   Level: advanced

1143: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1144: @*/
1145: PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has)
1146: {
1147:   H5O_type_t type;
1148:   char *path;

1155:   *has = PETSC_FALSE;
1156:   if (!obj->name) return(0);
1157:   PetscViewerHDF5GetAbsolutePath_Internal(viewer, obj->name, &path);
1158:   PetscViewerHDF5Traverse_Internal(viewer, path, PETSC_FALSE, has, &type);
1159:   *has = (type == H5O_TYPE_DATASET) ? PETSC_TRUE : PETSC_FALSE;
1160:   PetscFree(path);
1161:   return(0);
1162: }

1164: /*@C
1165:  PetscViewerHDF5HasAttribute - Check whether an attribute exists

1167:   Input Parameters:
1168: + viewer - The HDF5 viewer
1169: . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
1170: - name   - The attribute name

1172:   Output Parameter:
1173: . has    - Flag for attribute existence

1175:   Level: advanced

1177: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1178: @*/
1179: PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscBool *has)
1180: {
1181:   char           *parent;

1189:   PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);
1190:   PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, has, NULL);
1191:   if (*has) {PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, has);}
1192:   PetscFree(parent);
1193:   return(0);
1194: }

1196: /*@C
1197:  PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name

1199:   Input Parameters:
1200: + viewer - The HDF5 viewer
1201: . obj    - The object whose name is used to lookup the parent dataset, relative to the current group.
1202: - name   - The attribute name

1204:   Output Parameter:
1205: . has    - Flag for attribute existence

1207:   Notes:
1208:   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1209:   You might want to check first if it does using PetscViewerHDF5HasObject().

1211:   Level: advanced

1213: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1214: @*/
1215: PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has)
1216: {

1224:   PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);
1225:   PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);
1226:   return(0);
1227: }

1229: static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
1230: {
1231:   hid_t          h5;
1232:   htri_t         hhas;

1236:   PetscViewerHDF5GetFileId(viewer, &h5);
1237:   PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT));
1238:   *has = hhas ? PETSC_TRUE : PETSC_FALSE;
1239:   return(0);
1240: }

1242: /*
1243:   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
1244:   is attached to a communicator, in this case the attribute is a PetscViewer.
1245: */
1246: PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;

1248: /*@C
1249:   PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator.

1251:   Collective

1253:   Input Parameter:
1254: . comm - the MPI communicator to share the HDF5 PetscViewer

1256:   Level: intermediate

1258:   Options Database Keys:
1259: . -viewer_hdf5_filename <name> - name of the HDF5 file

1261:   Environmental variables:
1262: . PETSC_VIEWER_HDF5_FILENAME - name of the HDF5 file

1264:   Notes:
1265:   Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return
1266:   an error code.  The HDF5 PetscViewer is usually used in the form
1267: $       XXXView(XXX object, PETSC_VIEWER_HDF5_(comm));

1269: .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy()
1270: @*/
1271: PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
1272: {
1274:   PetscBool      flg;
1275:   PetscViewer    viewer;
1276:   char           fname[PETSC_MAX_PATH_LEN];
1277:   MPI_Comm       ncomm;

1280:   PetscCommDuplicate(comm,&ncomm,NULL);if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1281:   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
1282:     MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,NULL);
1283:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1284:   }
1285:   MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg);
1286:   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1287:   if (!flg) { /* PetscViewer not yet created */
1288:     PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg);
1289:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1290:     if (!flg) {
1291:       PetscStrcpy(fname,"output.h5");
1292:       if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1293:     }
1294:     PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer);
1295:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1296:     PetscObjectRegisterDestroy((PetscObject)viewer);
1297:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1298:     MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer);
1299:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1300:   }
1301:   PetscCommDestroy(&ncomm);
1302:   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1303:   PetscFunctionReturn(viewer);
1304: }