Actual source code: hdf5v.c

petsc-master 2019-06-15
Report Typos and Errors
  1: #include <petsc/private/viewerimpl.h>
  2: #include <petsc/private/viewerhdf5impl.h>
  3: #include <petscviewerhdf5.h>    /*I   "petscviewerhdf5.h"   I*/

  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)
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\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; 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:   PetscViewer_HDF5  *hdf5 = (PetscViewer_HDF5*) viewer->data;
332:   H5FD_mpio_xfer_t  mode;

335:   PetscStackCallHDF5(H5Pget_dxpl_mpio,(hdf5->dxpl_id, &mode));
336:   *flg = (mode == H5FD_MPIO_COLLECTIVE) ? PETSC_TRUE : PETSC_FALSE;
337:   return(0);
338: }

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

343:   Not Collective

345:   Input Parameters:
346: . viewer - the HDF5 PetscViewer

348:   Output Parameters:
349: . flg - the flag

351:   Level: intermediate

353:   Notes:
354:   This setting works correctly only since HDF5 1.10.3. For older versions, PETSC_FALSE will be always returned.
355:   For more details, see PetscViewerHDF5SetCollective().

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

359: @*/
360: PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer viewer,PetscBool *flg)
361: {


368:   PetscUseMethod(viewer,"PetscViewerHDF5GetCollective_C",(PetscViewer,PetscBool*),(viewer,flg));
369:   return(0);
370: }

372: static PetscErrorCode  PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[])
373: {
374:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
375:   hid_t             plist_id;
376:   PetscErrorCode    ierr;

379:   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
380:   if (hdf5->filename) {PetscFree(hdf5->filename);}
381:   PetscStrallocpy(name, &hdf5->filename);
382:   /* Set up file access property list with parallel I/O access */
383:   PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS));
384:   PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), MPI_INFO_NULL));
385:   /* Create or open the file collectively */
386:   switch (hdf5->btype) {
387:   case FILE_MODE_READ:
388:     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id));
389:     break;
390:   case FILE_MODE_APPEND:
391:     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id));
392:     break;
393:   case FILE_MODE_WRITE:
394:     PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id));
395:     break;
396:   default:
397:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
398:   }
399:   if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name);
400:   PetscStackCallHDF5(H5Pclose,(plist_id));
401:   return(0);
402: }

404: static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name)
405: {
406:   PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data;

409:   *name = vhdf5->filename;
410:   return(0);
411: }

413: static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer)
414: {
415:   /*
416:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
417:   PetscErrorCode   ierr;
418:   */

421:   return(0);
422: }

424: /*MC
425:    PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file


428: .seealso:  PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET,
429:            PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING,
430:            PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB,
431:            PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType()

433:   Level: beginner
434: M*/

436: PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v)
437: {
438:   PetscViewer_HDF5 *hdf5;
439:   PetscErrorCode   ierr;

442:   PetscNewLog(v,&hdf5);

444:   v->data                = (void*) hdf5;
445:   v->ops->destroy        = PetscViewerDestroy_HDF5;
446:   v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5;
447:   v->ops->setup          = PetscViewerSetUp_HDF5;
448:   v->ops->view           = PetscViewerView_HDF5;
449:   v->ops->flush          = 0;
450:   hdf5->btype            = (PetscFileMode) -1;
451:   hdf5->filename         = 0;
452:   hdf5->timestep         = -1;
453:   hdf5->groups           = NULL;

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

457:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);
458:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);
459:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);
460:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_HDF5);
461:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);
462:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);
463:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetCollective_C",PetscViewerHDF5SetCollective_HDF5);
464:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetCollective_C",PetscViewerHDF5GetCollective_HDF5);
465:   return(0);
466: }

468: /*@C
469:    PetscViewerHDF5Open - Opens a file for HDF5 input/output.

471:    Collective

473:    Input Parameters:
474: +  comm - MPI communicator
475: .  name - name of file
476: -  type - type of file
477: $    FILE_MODE_WRITE - create new file for binary output
478: $    FILE_MODE_READ - open existing file for binary input
479: $    FILE_MODE_APPEND - open existing file for binary output

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

484:   Options Database:
485: .  -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
486: .  -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal

488:    Level: beginner

490:    Note:
491:    This PetscViewer should be destroyed with PetscViewerDestroy().


494: .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(),
495:           PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(),
496:           MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName()
497: @*/
498: PetscErrorCode  PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
499: {

503:   PetscViewerCreate(comm, hdf5v);
504:   PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);
505:   PetscViewerFileSetMode(*hdf5v, type);
506:   PetscViewerFileSetName(*hdf5v, name);
507:   PetscViewerSetFromOptions(*hdf5v);
508:   return(0);
509: }

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

514:   Not collective

516:   Input Parameter:
517: . viewer - the PetscViewer

519:   Output Parameter:
520: . file_id - The file id

522:   Level: intermediate

524: .seealso: PetscViewerHDF5Open()
525: @*/
526: PetscErrorCode  PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
527: {
528:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

532:   if (file_id) *file_id = hdf5->file_id;
533:   return(0);
534: }

536: /*@C
537:   PetscViewerHDF5PushGroup - Set the current HDF5 group for output

539:   Not collective

541:   Input Parameters:
542: + viewer - the PetscViewer
543: - name - The group name

545:   Level: intermediate

547: .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
548: @*/
549: PetscErrorCode  PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name)
550: {
551:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
552:   PetscViewerHDF5GroupList *groupNode;
553:   PetscErrorCode   ierr;

558:   PetscNew(&groupNode);
559:   PetscStrallocpy(name, (char**) &groupNode->name);

561:   groupNode->next = hdf5->groups;
562:   hdf5->groups    = groupNode;
563:   return(0);
564: }

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

569:   Not collective

571:   Input Parameter:
572: . viewer - the PetscViewer

574:   Level: intermediate

576: .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
577: @*/
578: PetscErrorCode  PetscViewerHDF5PopGroup(PetscViewer viewer)
579: {
580:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
581:   PetscViewerHDF5GroupList *groupNode;
582:   PetscErrorCode   ierr;

586:   if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
587:   groupNode    = hdf5->groups;
588:   hdf5->groups = hdf5->groups->next;
589:   PetscFree(groupNode->name);
590:   PetscFree(groupNode);
591:   return(0);
592: }

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

598:   Not collective

600:   Input Parameter:
601: . viewer - the PetscViewer

603:   Output Parameter:
604: . name - The group name

606:   Level: intermediate

608: .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup()
609: @*/
610: PetscErrorCode  PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name)
611: {
612:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data;

617:   if (hdf5->groups) *name = hdf5->groups->name;
618:   else *name = NULL;
619:   return(0);
620: }

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

627:   Not collective

629:   Input Parameter:
630: . viewer - the PetscViewer

632:   Output Parameter:
633: + fileId - The HDF5 file ID
634: - groupId - The HDF5 group ID

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

639:   Level: intermediate

641: .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
642: @*/
643: PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId)
644: {
645:   hid_t          file_id;
646:   H5O_type_t     type;
647:   const char     *groupName = NULL;
648:   PetscBool      create;

652:   PetscViewerWritable(viewer, &create);
653:   PetscViewerHDF5GetFileId(viewer, &file_id);
654:   PetscViewerHDF5GetGroup(viewer, &groupName);
655:   PetscViewerHDF5Traverse_Internal(viewer, groupName, create, NULL, &type);
656:   if (type != H5O_TYPE_GROUP) SETERRQ1(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Path %s resolves to something which is not a group", groupName);
657:   PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT));
658:   *fileId  = file_id;
659:   return(0);
660: }

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

665:   Not collective

667:   Input Parameter:
668: . viewer - the PetscViewer

670:   Level: intermediate

672: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep()
673: @*/
674: PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
675: {
676:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

680:   ++hdf5->timestep;
681:   return(0);
682: }

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

688:   Not collective

690:   Input Parameters:
691: + viewer - the PetscViewer
692: - timestep - The timestep number

694:   Level: intermediate

696: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
697: @*/
698: PetscErrorCode  PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
699: {
700:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

704:   hdf5->timestep = timestep;
705:   return(0);
706: }

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

711:   Not collective

713:   Input Parameter:
714: . viewer - the PetscViewer

716:   Output Parameter:
717: . timestep - The timestep number

719:   Level: intermediate

721: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep()
722: @*/
723: PetscErrorCode  PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
724: {
725:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

730:   *timestep = hdf5->timestep;
731:   return(0);
732: }

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

737:   Not collective

739:   Input Parameter:
740: . ptype - the PETSc datatype name (for example PETSC_DOUBLE)

742:   Output Parameter:
743: . mtype - the MPI datatype (for example MPI_DOUBLE, ...)

745:   Level: advanced

747: .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
748: @*/
749: PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
750: {
752:   if (ptype == PETSC_INT)
753: #if defined(PETSC_USE_64BIT_INDICES)
754:                                        *htype = H5T_NATIVE_LLONG;
755: #else
756:                                        *htype = H5T_NATIVE_INT;
757: #endif
758:   else if (ptype == PETSC_DOUBLE)      *htype = H5T_NATIVE_DOUBLE;
759:   else if (ptype == PETSC_LONG)        *htype = H5T_NATIVE_LONG;
760:   else if (ptype == PETSC_SHORT)       *htype = H5T_NATIVE_SHORT;
761:   else if (ptype == PETSC_ENUM)        *htype = H5T_NATIVE_INT;
762:   else if (ptype == PETSC_BOOL)        *htype = H5T_NATIVE_INT;
763:   else if (ptype == PETSC_FLOAT)       *htype = H5T_NATIVE_FLOAT;
764:   else if (ptype == PETSC_CHAR)        *htype = H5T_NATIVE_CHAR;
765:   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
766:   else if (ptype == PETSC_STRING)      *htype = H5Tcopy(H5T_C_S1);
767:   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
768:   return(0);
769: }

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

774:   Not collective

776:   Input Parameter:
777: . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...)

779:   Output Parameter:
780: . ptype - the PETSc datatype name (for example PETSC_DOUBLE)

782:   Level: advanced

784: .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
785: @*/
786: PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
787: {
789: #if defined(PETSC_USE_64BIT_INDICES)
790:   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_LONG;
791:   else if (htype == H5T_NATIVE_LLONG)  *ptype = PETSC_INT;
792: #else
793:   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_INT;
794: #endif
795:   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
796:   else if (htype == H5T_NATIVE_LONG)   *ptype = PETSC_LONG;
797:   else if (htype == H5T_NATIVE_SHORT)  *ptype = PETSC_SHORT;
798:   else if (htype == H5T_NATIVE_FLOAT)  *ptype = PETSC_FLOAT;
799:   else if (htype == H5T_NATIVE_CHAR)   *ptype = PETSC_CHAR;
800:   else if (htype == H5T_NATIVE_UCHAR)  *ptype = PETSC_CHAR;
801:   else if (htype == H5T_C_S1)          *ptype = PETSC_STRING;
802:   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
803:   return(0);
804: }

806: /*@C
807:  PetscViewerHDF5WriteAttribute - Write an attribute

809:   Input Parameters:
810: + viewer - The HDF5 viewer
811: . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
812: . name   - The attribute name
813: . datatype - The attribute type
814: - value    - The attribute value

816:   Level: advanced

818: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
819: @*/
820: PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, const void *value)
821: {
822:   char           *parent;
823:   hid_t          h5, dataspace, obj, attribute, dtype;
824:   PetscBool      has;

832:   PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);
833:   PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_TRUE, NULL, NULL);
834:   PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);
835:   PetscDataTypeToHDF5DataType(datatype, &dtype);
836:   if (datatype == PETSC_STRING) {
837:     size_t len;
838:     PetscStrlen((const char *) value, &len);
839:     PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
840:   }
841:   PetscViewerHDF5GetFileId(viewer, &h5);
842:   PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR));
843:   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
844:   if (has) {
845:     PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
846:   } else {
847:     PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
848:   }
849:   PetscStackCallHDF5(H5Awrite,(attribute, dtype, value));
850:   if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype));
851:   PetscStackCallHDF5(H5Aclose,(attribute));
852:   PetscStackCallHDF5(H5Oclose,(obj));
853:   PetscStackCallHDF5(H5Sclose,(dataspace));
854:   PetscFree(parent);
855:   return(0);
856: }

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

861:   Input Parameters:
862: + viewer   - The HDF5 viewer
863: . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
864: . name     - The attribute name
865: . datatype - The attribute type
866: - value    - The attribute value

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

872:   Level: advanced

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

885:   PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);
886:   PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);
887:   return(0);
888: }

890: /*@C
891:  PetscViewerHDF5ReadAttribute - Read an attribute

893:   Input Parameters:
894: + viewer - The HDF5 viewer
895: . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
896: . name   - The attribute name
897: - datatype - The attribute type

899:   Output Parameter:
900: . value    - The attribute value

902:   Level: advanced

904: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
905: @*/
906: PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, void *value)
907: {
908:   char           *parent;
909:   hid_t          h5, obj, attribute, atype, dtype;
910:   PetscBool      has;

918:   PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);
919:   PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, &has, NULL);
920:   if (has) {PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);}
921:   if (!has) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist", parent, name);
922:   PetscDataTypeToHDF5DataType(datatype, &dtype);
923:   PetscViewerHDF5GetFileId(viewer, &h5);
924:   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
925:   PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
926:   if (datatype == PETSC_STRING) {
927:     size_t len;
928:     PetscStackCallHDF5Return(atype,H5Aget_type,(attribute));
929:     PetscStackCall("H5Tget_size",len = H5Tget_size(atype));
930:     PetscStackCallHDF5(H5Tclose,(atype));
931:     PetscMalloc((len+1) * sizeof(char *), &value);
932:   }
933:   PetscStackCallHDF5(H5Aread,(attribute, dtype, value));
934:   PetscStackCallHDF5(H5Aclose,(attribute));
935:   /* H5Oclose can be used to close groups, datasets, or committed datatypes */
936:   PetscStackCallHDF5(H5Oclose,(obj));
937:   PetscFree(parent);
938:   return(0);
939: }

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

944:   Input Parameters:
945: + viewer   - The HDF5 viewer
946: . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
947: . name     - The attribute name
948: - datatype - The attribute type

950:   Output Parameter:
951: . value    - The attribute value

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

957:   Level: advanced

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

970:   PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);
971:   PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, value);
972:   return(0);
973: }

975: PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_)
976: {
977:   htri_t exists;
978:   hid_t group;

981:   PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT));
982:   if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT));
983:   if (!exists && createGroup) {
984:     PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT));
985:     PetscStackCallHDF5(H5Gclose,(group));
986:     exists = PETSC_TRUE;
987:   }
988:   *exists_ = (PetscBool) exists;
989:   return(0);
990: }

992: static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype)
993: {
994:   const char     rootGroupName[] = "/";
995:   hid_t          h5;
996:   PetscBool      exists=PETSC_FALSE;
997:   PetscInt       i;
998:   int            n;
999:   char           **hierarchy;
1000:   char           buf[PETSC_MAX_PATH_LEN]="";

1006:   else name = rootGroupName;
1007:   if (has) {
1009:     *has = PETSC_FALSE;
1010:   }
1011:   if (otype) {
1013:     *otype = H5O_TYPE_UNKNOWN;
1014:   }
1015:   PetscViewerHDF5GetFileId(viewer, &h5);

1017:   /*
1018:      Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing.
1019:      Hence, each of them needs to be tested separately:
1020:      1) whether it's a valid link
1021:      2) whether this link resolves to an object
1022:      See H5Oexists_by_name() documentation.
1023:   */
1024:   PetscStrToArray(name,'/',&n,&hierarchy);
1025:   if (!n) {
1026:     /*  Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */
1027:     if (has)   *has   = PETSC_TRUE;
1028:     if (otype) *otype = H5O_TYPE_GROUP;
1029:     PetscStrToArrayDestroy(n,hierarchy);
1030:     return(0);
1031:   }
1032:   for (i=0; i<n; i++) {
1033:     PetscStrcat(buf,"/");
1034:     PetscStrcat(buf,hierarchy[i]);
1035:     PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);
1036:     if (!exists) break;
1037:   }
1038:   PetscStrToArrayDestroy(n,hierarchy);

1040:   /* If the object exists, get its type */
1041:   if (exists && otype) {
1042:     H5O_info_t info;

1044:     /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */
1045:     PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT));
1046:     *otype = info.type;
1047:   }
1048:   if (has) *has = exists;
1049:   return(0);
1050: }

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

1055:   Input Parameters:
1056: . viewer - The HDF5 viewer

1058:   Output Parameter:
1059: . has    - Flag for group existence

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

1064:   Level: advanced

1066: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5OpenGroup()
1067: @*/
1068: PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, PetscBool *has)
1069: {
1070:   H5O_type_t type;
1071:   const char *name;

1077:   PetscViewerHDF5GetGroup(viewer, &name);
1078:   PetscViewerHDF5Traverse_Internal(viewer, name, PETSC_FALSE, has, &type);
1079:   *has = (type == H5O_TYPE_GROUP) ? PETSC_TRUE : PETSC_FALSE;
1080:   return(0);
1081: }

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

1086:   Input Parameters:
1087: + viewer - The HDF5 viewer
1088: - obj    - The named object

1090:   Output Parameter:
1091: . has    - Flag for dataset existence; PETSC_FALSE for unnamed object

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

1096:   Level: advanced

1098: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1099: @*/
1100: PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has)
1101: {
1102:   H5O_type_t type;
1103:   char *path;

1110:   *has = PETSC_FALSE;
1111:   if (!obj->name) return(0);
1112:   PetscViewerHDF5GetAbsolutePath_Internal(viewer, obj->name, &path);
1113:   PetscViewerHDF5Traverse_Internal(viewer, path, PETSC_FALSE, has, &type);
1114:   *has = (type == H5O_TYPE_DATASET) ? PETSC_TRUE : PETSC_FALSE;
1115:   PetscFree(path);
1116:   return(0);
1117: }

1119: /*@C
1120:  PetscViewerHDF5HasAttribute - Check whether an attribute exists

1122:   Input Parameters:
1123: + viewer - The HDF5 viewer
1124: . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
1125: - name   - The attribute name

1127:   Output Parameter:
1128: . has    - Flag for attribute existence

1130:   Level: advanced

1132: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1133: @*/
1134: PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscBool *has)
1135: {
1136:   char           *parent;

1144:   PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);
1145:   PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, has, NULL);
1146:   if (*has) {PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, has);}
1147:   PetscFree(parent);
1148:   return(0);
1149: }

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

1154:   Input Parameters:
1155: + viewer - The HDF5 viewer
1156: . obj    - The object whose name is used to lookup the parent dataset, relative to the current group.
1157: - name   - The attribute name

1159:   Output Parameter:
1160: . has    - Flag for attribute existence

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

1166:   Level: advanced

1168: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1169: @*/
1170: PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has)
1171: {

1179:   PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);
1180:   PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);
1181:   return(0);
1182: }

1184: static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
1185: {
1186:   hid_t          h5;
1187:   htri_t         hhas;

1191:   PetscViewerHDF5GetFileId(viewer, &h5);
1192:   PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT));
1193:   *has = hhas ? PETSC_TRUE : PETSC_FALSE;
1194:   return(0);
1195: }

1197: /*
1198:   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
1199:   is attached to a communicator, in this case the attribute is a PetscViewer.
1200: */
1201: PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;

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

1206:   Collective

1208:   Input Parameter:
1209: . comm - the MPI communicator to share the HDF5 PetscViewer

1211:   Level: intermediate

1213:   Options Database Keys:
1214: . -viewer_hdf5_filename <name>

1216:   Environmental variables:
1217: . PETSC_VIEWER_HDF5_FILENAME

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

1224: .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy()
1225: @*/
1226: PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
1227: {
1229:   PetscBool      flg;
1230:   PetscViewer    viewer;
1231:   char           fname[PETSC_MAX_PATH_LEN];
1232:   MPI_Comm       ncomm;

1235:   PetscCommDuplicate(comm,&ncomm,NULL);if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
1236:   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
1237:     MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0);
1238:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
1239:   }
1240:   MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg);
1241:   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
1242:   if (!flg) { /* PetscViewer not yet created */
1243:     PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg);
1244:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
1245:     if (!flg) {
1246:       PetscStrcpy(fname,"output.h5");
1247:       if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
1248:     }
1249:     PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer);
1250:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
1251:     PetscObjectRegisterDestroy((PetscObject)viewer);
1252:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
1253:     MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer);
1254:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
1255:   }
1256:   PetscCommDestroy(&ncomm);
1257:   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
1258:   PetscFunctionReturn(viewer);
1259: }