Actual source code: hdf5v.c

  1: #include <petsc/private/viewerhdf5impl.h>

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

  6: static PetscErrorCode PetscViewerHDF5GetAbsolutePath_Internal(PetscViewer viewer, const char path[], char *abspath[])
  7: {
  8:   PetscBool      relative = PETSC_FALSE;
  9:   const char     *group;
 10:   char           buf[PETSC_MAX_PATH_LEN] = "";

 14:   PetscViewerHDF5GetGroup(viewer, &group);
 15:   PetscViewerHDF5PathIsRelative(path, PETSC_TRUE, &relative);
 16:   if (relative) {
 17:     PetscStrcpy(buf, group);
 18:     PetscStrcat(buf, "/");
 19:     PetscStrcat(buf, path);
 20:     PetscStrallocpy(buf, abspath);
 21:   } else {
 22:     PetscStrallocpy(path, abspath);
 23:   }
 24:   return(0);
 25: }

 27: static PetscErrorCode PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer, PetscObject obj)
 28: {
 29:   PetscBool      has;

 33:   PetscViewerHDF5HasObject(viewer, obj, &has);
 34:   if (!has) {
 35:     const char *group;
 36:     PetscViewerHDF5GetGroup(viewer, &group);
 37:     SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) \"%s\" not stored in group %s", obj->name, group ? group : "/");
 38:   }
 39:   return(0);
 40: }

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

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

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

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

 75: static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer)
 76: {
 77:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data;
 78:   PetscErrorCode   ierr;

 81:   PetscFree(hdf5->filename);
 82:   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
 83:   return(0);
 84: }

 86: static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer)
 87: {
 88:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
 89:   PetscErrorCode   ierr;

 92:   PetscStackCallHDF5(H5Pclose,(hdf5->dxpl_id));
 93:   PetscViewerFileClose_HDF5(viewer);
 94:   while (hdf5->groups) {
 95:     PetscViewerHDF5GroupList *tmp = hdf5->groups->next;

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

110: static PetscErrorCode  PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type)
111: {
112:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

115:   hdf5->btype = type;
116:   return(0);
117: }

119: static PetscErrorCode  PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type)
120: {
121:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

124:   *type = hdf5->btype;
125:   return(0);
126: }

128: static PetscErrorCode  PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg)
129: {
130:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

133:   hdf5->basedimension2 = flg;
134:   return(0);
135: }

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

141:     Logically Collective on PetscViewer

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

147:   Options Database:
148: .  -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

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

154:   Level: intermediate

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

158: @*/
159: PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg)
160: {

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

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

173:     Logically Collective on PetscViewer

175:   Input Parameter:
176: .  viewer - the PetscViewer, must be of type HDF5

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

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

185:   Level: intermediate

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

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

196:   *flg = hdf5->basedimension2;
197:   return(0);
198: }

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

205:   hdf5->spoutput = flg;
206:   return(0);
207: }

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

213:     Logically Collective on PetscViewer

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

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

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

226:   Level: intermediate

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

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

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

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

246:     Logically Collective on PetscViewer

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

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

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

258:   Level: intermediate

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

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

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

274: static PetscErrorCode  PetscViewerHDF5SetCollective_HDF5(PetscViewer viewer, PetscBool flg)
275: {
277:   /* H5FD_MPIO_COLLECTIVE is wrong in hdf5 1.10.2, and is the same as H5FD_MPIO_INDEPENDENT in earlier versions
278:      - see e.g. https://gitlab.cosma.dur.ac.uk/swift/swiftsim/issues/431 */
279: #if H5_VERSION_GE(1,10,3) && defined(H5_HAVE_PARALLEL)
280:   {
281:     PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
282:     PetscStackCallHDF5(H5Pset_dxpl_mpio,(hdf5->dxpl_id, flg ? H5FD_MPIO_COLLECTIVE : H5FD_MPIO_INDEPENDENT));
283:   }
284: #else
285:   if (flg) {
287:     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");
288:   }
289: #endif
290:   return(0);
291: }

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

296:   Logically Collective; flg must contain common value

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

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

305:   Notes:
306:   Collective mode gives the MPI-IO layer underneath HDF5 a chance to do some additional collective optimizations and hence can perform better.
307:   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.

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

314:   Level: intermediate

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

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

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

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

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

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

350:   Not Collective

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

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

358:   Level: intermediate

360:   Notes:
361:   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.
362:   For more details, see PetscViewerHDF5SetCollective().

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

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


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

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

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

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

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

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

438:   return(0);
439: }

441: /*MC
442:    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().

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

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

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

542:   Not collective

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

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

550:   Level: intermediate

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

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

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

567:   Not collective

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

573:   Level: intermediate

575:   Notes:
576:   This is designed to mnemonically resemble the Unix cd command.
577:   + If name begins with '/', it is interpreted as an absolute path fully replacing current group, otherwise it is taken as relative to the current group.
578:   . NULL, empty string, or any sequence of all slashes (e.g. "///") is interpreted as the root group "/".
579:   - "." means the current group is pushed again.

581:   Example:
582:   Suppose the current group is "/a".
583:   + If name is NULL, empty string, or a sequence of all slashes (e.g. "///"), then the new group will be "/".
584:   . If name is ".", then the new group will be "/a".
585:   . If name is "b", then the new group will be "/a/b".
586:   - If name is "/b", then the new group will be "/b".

588:   Developer Notes:
589:   The root group "/" is internally stored as NULL.

591: .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
592: @*/
593: PetscErrorCode  PetscViewerHDF5PushGroup(PetscViewer viewer, const char name[])
594: {
595:   PetscViewer_HDF5          *hdf5 = (PetscViewer_HDF5*) viewer->data;
596:   PetscViewerHDF5GroupList  *groupNode;
597:   size_t                    i,len;
598:   char                      buf[PETSC_MAX_PATH_LEN];
599:   const char                *gname;
600:   PetscErrorCode            ierr;

605:   PetscStrlen(name, &len);
606:   gname = NULL;
607:   if (len) {
608:     if (len == 1 && name[0] == '.') {
609:       /* use current name */
610:       gname = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : NULL;
611:     } else if (name[0] == '/') {
612:       /* absolute */
613:       for (i=1; i<len; i++) {
614:         if (name[i] != '/') {
615:           gname = name;
616:           break;
617:         }
618:       }
619:     } else {
620:       /* relative */
621:       const char *parent = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : "";
622:       PetscSNPrintf(buf, sizeof(buf), "%s/%s", parent, name);
623:       gname = buf;
624:     }
625:   }
626:   PetscNew(&groupNode);
627:   PetscStrallocpy(gname, (char**) &groupNode->name);
628:   groupNode->next = hdf5->groups;
629:   hdf5->groups    = groupNode;
630:   return(0);
631: }

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

636:   Not collective

638:   Input Parameter:
639: . viewer - the PetscViewer

641:   Level: intermediate

643: .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
644: @*/
645: PetscErrorCode  PetscViewerHDF5PopGroup(PetscViewer viewer)
646: {
647:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
648:   PetscViewerHDF5GroupList *groupNode;
649:   PetscErrorCode   ierr;

653:   if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
654:   groupNode    = hdf5->groups;
655:   hdf5->groups = hdf5->groups->next;
656:   PetscFree(groupNode->name);
657:   PetscFree(groupNode);
658:   return(0);
659: }

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

665:   Not collective

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

670:   Output Parameter:
671: . name - The group name

673:   Level: intermediate

675: .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup()
676: @*/
677: PetscErrorCode  PetscViewerHDF5GetGroup(PetscViewer viewer, const char *name[])
678: {
679:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data;

684:   if (hdf5->groups) *name = hdf5->groups->name;
685:   else *name = NULL;
686:   return(0);
687: }

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

694:   Not collective

696:   Input Parameter:
697: . viewer - the PetscViewer

699:   Output Parameters:
700: + fileId - The HDF5 file ID
701: - groupId - The HDF5 group ID

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

706:   Level: intermediate

708: .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
709: @*/
710: PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId)
711: {
712:   hid_t          file_id;
713:   H5O_type_t     type;
714:   const char     *groupName = NULL, *fileName = NULL;
715:   PetscBool      writable, has;

719:   PetscViewerWritable(viewer, &writable);
720:   PetscViewerHDF5GetFileId(viewer, &file_id);
721:   PetscViewerFileGetName(viewer, &fileName);
722:   PetscViewerHDF5GetGroup(viewer, &groupName);
723:   PetscViewerHDF5Traverse_Internal(viewer, groupName, writable, &has, &type);
724:   if (!has) {
725:     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);
726:     else           SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_LIB, "HDF5 failed to create group %s although file %s is open for writing", groupName, fileName);
727:   }
728:   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);
729:   PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT));
730:   *fileId  = file_id;
731:   return(0);
732: }

734: /*@
735:   PetscViewerHDF5PushTimestepping - Activate timestepping mode for subsequent HDF5 reading and writing.

737:   Not collective

739:   Input Parameter:
740: . viewer - the PetscViewer

742:   Level: intermediate

744:   Notes:
745:   On first PetscViewerHDF5PushTimestepping(), the initial time step is set to 0.
746:   Next timesteps can then be set using PetscViewerHDF5IncrementTimestep() or PetscViewerHDF5SetTimestep().
747:   Current timestep value determines which timestep is read from or written to any dataset on the next HDF5 I/O operation [e.g. VecView()].
748:   Use PetscViewerHDF5PopTimestepping() to deactivate timestepping mode; calling it by the end of the program is NOT mandatory.
749:   Current timestep is remembered between PetscViewerHDF5PopTimestepping() and the next PetscViewerHDF5PushTimestepping().

751:   If a dataset was stored with timestepping, it can be loaded only in the timestepping mode again.
752:   Loading a timestepped dataset with timestepping disabled, or vice-versa results in an error.

754:   Developer notes:
755:   Timestepped HDF5 dataset has an extra dimension and attribute "timestepping" set to true.

757: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PopTimestepping(), PetscViewerHDF5IsTimestepping(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
758: @*/
759: PetscErrorCode  PetscViewerHDF5PushTimestepping(PetscViewer viewer)
760: {
761:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

765:   if (hdf5->timestepping) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping is already pushed");
766:   hdf5->timestepping = PETSC_TRUE;
767:   if (hdf5->timestep < 0) hdf5->timestep = 0;
768:   return(0);
769: }

771: /*@
772:   PetscViewerHDF5PopTimestepping - Deactivate timestepping mode for subsequent HDF5 reading and writing.

774:   Not collective

776:   Input Parameter:
777: . viewer - the PetscViewer

779:   Level: intermediate

781:   Notes:
782:   See PetscViewerHDF5PushTimestepping() for details.

784: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5IsTimestepping(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
785: @*/
786: PetscErrorCode  PetscViewerHDF5PopTimestepping(PetscViewer viewer)
787: {
788:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

792:   if (!hdf5->timestepping) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
793:   hdf5->timestepping = PETSC_FALSE;
794:   return(0);
795: }

797: /*@
798:   PetscViewerHDF5IsTimestepping - Ask the viewer whether it is in timestepping mode currently.

800:   Not collective

802:   Input Parameter:
803: . viewer - the PetscViewer

805:   Output Parameter:
806: . flg - is timestepping active?

808:   Level: intermediate

810:   Notes:
811:   See PetscViewerHDF5PushTimestepping() for details.

813: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5PopTimestepping(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
814: @*/
815: PetscErrorCode  PetscViewerHDF5IsTimestepping(PetscViewer viewer, PetscBool *flg)
816: {
817:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

821:   *flg = hdf5->timestepping;
822:   return(0);
823: }

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

828:   Not collective

830:   Input Parameter:
831: . viewer - the PetscViewer

833:   Level: intermediate

835:   Notes:
836:   This can be called only if the viewer is in timestepping mode. See PetscViewerHDF5PushTimestepping() for details.

838: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep()
839: @*/
840: PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
841: {
842:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

846:   if (!hdf5->timestepping) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
847:   ++hdf5->timestep;
848:   return(0);
849: }

851: /*@
852:   PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time.

854:   Logically collective

856:   Input Parameters:
857: + viewer - the PetscViewer
858: - timestep - The timestep

860:   Level: intermediate

862:   Notes:
863:   This can be called only if the viewer is in timestepping mode. See PetscViewerHDF5PushTimestepping() for details.

865: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
866: @*/
867: PetscErrorCode  PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
868: {
869:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

874:   if (timestep < 0) SETERRQ1(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestep %D is negative", timestep);
875:   if (!hdf5->timestepping) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
876:   hdf5->timestep = timestep;
877:   return(0);
878: }

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

883:   Not collective

885:   Input Parameter:
886: . viewer - the PetscViewer

888:   Output Parameter:
889: . timestep - The timestep

891:   Level: intermediate

893:   Notes:
894:   This can be called only if the viewer is in the timestepping mode. See PetscViewerHDF5PushTimestepping() for details.

896: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep()
897: @*/
898: PetscErrorCode  PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
899: {
900:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

905:   if (!hdf5->timestepping) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
906:   *timestep = hdf5->timestep;
907:   return(0);
908: }

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

913:   Not collective

915:   Input Parameter:
916: . ptype - the PETSc datatype name (for example PETSC_DOUBLE)

918:   Output Parameter:
919: . mtype - the MPI datatype (for example MPI_DOUBLE, ...)

921:   Level: advanced

923: .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
924: @*/
925: PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
926: {
928:   if (ptype == PETSC_INT)
929: #if defined(PETSC_USE_64BIT_INDICES)
930:                                        *htype = H5T_NATIVE_LLONG;
931: #else
932:                                        *htype = H5T_NATIVE_INT;
933: #endif
934:   else if (ptype == PETSC_DOUBLE)      *htype = H5T_NATIVE_DOUBLE;
935:   else if (ptype == PETSC_LONG)        *htype = H5T_NATIVE_LONG;
936:   else if (ptype == PETSC_SHORT)       *htype = H5T_NATIVE_SHORT;
937:   else if (ptype == PETSC_ENUM)        *htype = H5T_NATIVE_INT;
938:   else if (ptype == PETSC_BOOL)        *htype = H5T_NATIVE_INT;
939:   else if (ptype == PETSC_FLOAT)       *htype = H5T_NATIVE_FLOAT;
940:   else if (ptype == PETSC_CHAR)        *htype = H5T_NATIVE_CHAR;
941:   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
942:   else if (ptype == PETSC_STRING)      *htype = H5Tcopy(H5T_C_S1);
943:   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
944:   return(0);
945: }

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

950:   Not collective

952:   Input Parameter:
953: . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...)

955:   Output Parameter:
956: . ptype - the PETSc datatype name (for example PETSC_DOUBLE)

958:   Level: advanced

960: .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
961: @*/
962: PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
963: {
965: #if defined(PETSC_USE_64BIT_INDICES)
966:   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_LONG;
967:   else if (htype == H5T_NATIVE_LLONG)  *ptype = PETSC_INT;
968: #else
969:   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_INT;
970: #endif
971:   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
972:   else if (htype == H5T_NATIVE_LONG)   *ptype = PETSC_LONG;
973:   else if (htype == H5T_NATIVE_SHORT)  *ptype = PETSC_SHORT;
974:   else if (htype == H5T_NATIVE_FLOAT)  *ptype = PETSC_FLOAT;
975:   else if (htype == H5T_NATIVE_CHAR)   *ptype = PETSC_CHAR;
976:   else if (htype == H5T_NATIVE_UCHAR)  *ptype = PETSC_CHAR;
977:   else if (htype == H5T_C_S1)          *ptype = PETSC_STRING;
978:   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
979:   return(0);
980: }

982: /*@C
983:  PetscViewerHDF5WriteAttribute - Write an attribute

985:   Collective

987:   Input Parameters:
988: + viewer - The HDF5 viewer
989: . parent - The parent dataset/group name
990: . name   - The attribute name
991: . datatype - The attribute type
992: - value    - The attribute value

994:   Level: advanced

996:   Notes:
997:   If parent starts with '/', it is taken as an absolute path overriding currently pushed group, else parent is relative to the current pushed group. NULL means the current pushed group.

999: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1000: @*/
1001: PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value)
1002: {
1003:   char           *parentAbsPath;
1004:   hid_t          h5, dataspace, obj, attribute, dtype;
1005:   PetscBool      has;

1014:   PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath);
1015:   PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_TRUE, NULL, NULL);
1016:   PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has);
1017:   PetscDataTypeToHDF5DataType(datatype, &dtype);
1018:   if (datatype == PETSC_STRING) {
1019:     size_t len;
1020:     PetscStrlen((const char *) value, &len);
1021:     PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
1022:   }
1023:   PetscViewerHDF5GetFileId(viewer, &h5);
1024:   PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR));
1025:   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parentAbsPath, H5P_DEFAULT));
1026:   if (has) {
1027:     PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
1028:   } else {
1029:     PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
1030:   }
1031:   PetscStackCallHDF5(H5Awrite,(attribute, dtype, value));
1032:   if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype));
1033:   PetscStackCallHDF5(H5Aclose,(attribute));
1034:   PetscStackCallHDF5(H5Oclose,(obj));
1035:   PetscStackCallHDF5(H5Sclose,(dataspace));
1036:   PetscFree(parentAbsPath);
1037:   return(0);
1038: }

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

1043:   Collective

1045:   Input Parameters:
1046: + viewer   - The HDF5 viewer
1047: . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
1048: . name     - The attribute name
1049: . datatype - The attribute type
1050: - value    - The attribute value

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

1056:   Level: advanced

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

1069:   PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);
1070:   PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);
1071:   return(0);
1072: }

1074: /*@C
1075:  PetscViewerHDF5ReadAttribute - Read an attribute

1077:   Collective

1079:   Input Parameters:
1080: + viewer - The HDF5 viewer
1081: . parent - The parent dataset/group name
1082: . name   - The attribute name
1083: . datatype - The attribute type
1084: - defaultValue - The pointer to the default value

1086:   Output Parameter:
1087: . value    - The pointer to the read HDF5 attribute value

1089:   Notes:
1090:   If defaultValue is NULL and the attribute is not found, an error occurs.
1091:   If defaultValue is not NULL and the attribute is not found, defaultValue is copied to value.
1092:   The pointers defaultValue and value can be the same; for instance
1093: $  flg = PETSC_FALSE;
1094: $  PetscViewerHDF5ReadAttribute(viewer,name,"attr",PETSC_BOOL,&flg,&flg);
1095:   is valid, but make sure the default value is initialized.

1097:   If the datatype is PETSC_STRING, the output string is newly allocated so one must PetscFree() it when no longer needed.

1099:   If parent starts with '/', it is taken as an absolute path overriding currently pushed group, else parent is relative to the current pushed group. NULL means the current pushed group.

1101:   Level: advanced

1103: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1104: @*/
1105: PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *defaultValue, void *value)
1106: {
1107:   char           *parentAbsPath;
1108:   hid_t          h5, obj, attribute, dtype;
1109:   PetscBool      has;

1118:   PetscDataTypeToHDF5DataType(datatype, &dtype);
1119:   PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath);
1120:   PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, &has, NULL);
1121:   if (has) {PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has);}
1122:   if (!has) {
1123:     if (defaultValue) {
1124:       if (defaultValue != value) {
1125:         if (datatype == PETSC_STRING) {
1126:           PetscStrallocpy(*(char**)defaultValue, (char**)value);
1127:         } else {
1128:           size_t len;
1129:           PetscStackCallHDF5Return(len,H5Tget_size,(dtype));
1130:           PetscMemcpy(value, defaultValue, len);
1131:         }
1132:       }
1133:       PetscFree(parentAbsPath);
1134:       return(0);
1135:     } else SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist and default value not provided", parentAbsPath, name);
1136:   }
1137:   PetscViewerHDF5GetFileId(viewer, &h5);
1138:   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parentAbsPath, H5P_DEFAULT));
1139:   PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
1140:   if (datatype == PETSC_STRING) {
1141:     size_t len;
1142:     hid_t  atype;
1143:     PetscStackCallHDF5Return(atype,H5Aget_type,(attribute));
1144:     PetscStackCallHDF5Return(len,H5Tget_size,(atype));
1145:     PetscMalloc((len+1) * sizeof(char), value);
1146:     PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
1147:     PetscStackCallHDF5(H5Aread,(attribute, dtype, *(char**)value));
1148:   } else {
1149:     PetscStackCallHDF5(H5Aread,(attribute, dtype, value));
1150:   }
1151:   PetscStackCallHDF5(H5Aclose,(attribute));
1152:   /* H5Oclose can be used to close groups, datasets, or committed datatypes */
1153:   PetscStackCallHDF5(H5Oclose,(obj));
1154:   PetscFree(parentAbsPath);
1155:   return(0);
1156: }

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

1161:   Collective

1163:   Input Parameters:
1164: + viewer   - The HDF5 viewer
1165: . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
1166: . name     - The attribute name
1167: - datatype - The attribute type

1169:   Output Parameter:
1170: . value    - The attribute value

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

1176:   Level: advanced

1178: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute() PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1179: @*/
1180: PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *defaultValue, void *value)
1181: {

1189:   PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);
1190:   PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, defaultValue, value);
1191:   return(0);
1192: }

1194: PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_)
1195: {
1196:   htri_t exists;
1197:   hid_t group;

1200:   PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT));
1201:   if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT));
1202:   if (!exists && createGroup) {
1203:     PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT));
1204:     PetscStackCallHDF5(H5Gclose,(group));
1205:     exists = PETSC_TRUE;
1206:   }
1207:   *exists_ = (PetscBool) exists;
1208:   return(0);
1209: }

1211: static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype)
1212: {
1213:   const char     rootGroupName[] = "/";
1214:   hid_t          h5;
1215:   PetscBool      exists=PETSC_FALSE;
1216:   PetscInt       i;
1217:   int            n;
1218:   char           **hierarchy;
1219:   char           buf[PETSC_MAX_PATH_LEN]="";

1225:   else name = rootGroupName;
1226:   if (has) {
1228:     *has = PETSC_FALSE;
1229:   }
1230:   if (otype) {
1232:     *otype = H5O_TYPE_UNKNOWN;
1233:   }
1234:   PetscViewerHDF5GetFileId(viewer, &h5);

1236:   /*
1237:      Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing.
1238:      Hence, each of them needs to be tested separately:
1239:      1) whether it's a valid link
1240:      2) whether this link resolves to an object
1241:      See H5Oexists_by_name() documentation.
1242:   */
1243:   PetscStrToArray(name,'/',&n,&hierarchy);
1244:   if (!n) {
1245:     /*  Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */
1246:     if (has)   *has   = PETSC_TRUE;
1247:     if (otype) *otype = H5O_TYPE_GROUP;
1248:     PetscStrToArrayDestroy(n,hierarchy);
1249:     return(0);
1250:   }
1251:   for (i=0; i<n; i++) {
1252:     PetscStrcat(buf,"/");
1253:     PetscStrcat(buf,hierarchy[i]);
1254:     PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);
1255:     if (!exists) break;
1256:   }
1257:   PetscStrToArrayDestroy(n,hierarchy);

1259:   /* If the object exists, get its type */
1260:   if (exists && otype) {
1261:     H5O_info_t info;

1263:     /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */
1264:     PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT));
1265:     *otype = info.type;
1266:   }
1267:   if (has) *has = exists;
1268:   return(0);
1269: }

1271: /*@C
1272:  PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file

1274:   Collective

1276:   Input Parameters:
1277: + viewer - The HDF5 viewer
1278: - path - The group path

1280:   Output Parameter:
1281: . has - Flag for group existence

1283:   Level: advanced

1285:   Notes:
1286:   If path starts with '/', it is taken as an absolute path overriding currently pushed group, else path is relative to the current pushed group.
1287:   So NULL or empty path means the current pushed group.
1288:   If path exists but is not a group, PETSC_FALSE is returned.

1290: .seealso: PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasDataset(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5GetGroup(), PetscViewerHDF5OpenGroup()
1291: @*/
1292: PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, const char path[], PetscBool *has)
1293: {
1294:   H5O_type_t     type;
1295:   char           *abspath;

1302:   PetscViewerHDF5GetAbsolutePath_Internal(viewer, path, &abspath);
1303:   PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type);
1304:   *has = (PetscBool)(type == H5O_TYPE_GROUP);
1305:   PetscFree(abspath);
1306:   return(0);
1307: }

1309: /*@C
1310:  PetscViewerHDF5HasDataset - Check whether a given dataset exists in the HDF5 file

1312:   Collective

1314:   Input Parameters:
1315: + viewer - The HDF5 viewer
1316: - path - The dataset path

1318:   Output Parameter:
1319: . has - Flag whether dataset exists

1321:   Level: advanced

1323:   Notes:
1324:   If path starts with '/', it is taken as an absolute path overriding currently pushed group, else path is relative to the current pushed group.
1325:   If path is NULL or empty, has is set to PETSC_FALSE.
1326:   If path exists but is not a dataset, has is set to PETSC_FALSE as well.

1328: .seealso: PetscViewerHDF5HasObject(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasGroup(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5GetGroup()
1329: @*/
1330: PetscErrorCode PetscViewerHDF5HasDataset(PetscViewer viewer, const char path[], PetscBool *has)
1331: {
1332:   H5O_type_t     type;
1333:   char           *abspath;

1340:   PetscViewerHDF5GetAbsolutePath_Internal(viewer, path, &abspath);
1341:   PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type);
1342:   *has = (PetscBool)(type == H5O_TYPE_DATASET);
1343:   PetscFree(abspath);
1344:   return(0);
1345: }

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

1350:   Collective

1352:   Input Parameters:
1353: + viewer - The HDF5 viewer
1354: - obj    - The named object

1356:   Output Parameter:
1357: . has    - Flag for dataset existence

1359:   Notes:
1360:   If the object is unnamed, an error occurs.
1361:   If the path current_group/object_name exists but is not a dataset, has is set to PETSC_FALSE as well.

1363:   Level: advanced

1365: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasDataset(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5GetGroup()
1366: @*/
1367: PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has)
1368: {
1369:   size_t         len;

1376:   PetscStrlen(obj->name, &len);
1377:   if (!len) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named");
1378:   PetscViewerHDF5HasDataset(viewer, obj->name, has);
1379:   return(0);
1380: }

1382: /*@C
1383:  PetscViewerHDF5HasAttribute - Check whether an attribute exists

1385:   Collective

1387:   Input Parameters:
1388: + viewer - The HDF5 viewer
1389: . parent - The parent dataset/group name
1390: - name   - The attribute name

1392:   Output Parameter:
1393: . has    - Flag for attribute existence

1395:   Level: advanced

1397:   Notes:
1398:   If parent starts with '/', it is taken as an absolute path overriding currently pushed group, else parent is relative to the current pushed group. NULL means the current pushed group.

1400: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1401: @*/
1402: PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
1403: {
1404:   char           *parentAbsPath;

1412:   PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath);
1413:   PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, has, NULL);
1414:   if (*has) {PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, has);}
1415:   PetscFree(parentAbsPath);
1416:   return(0);
1417: }

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

1422:   Collective

1424:   Input Parameters:
1425: + viewer - The HDF5 viewer
1426: . obj    - The object whose name is used to lookup the parent dataset, relative to the current group.
1427: - name   - The attribute name

1429:   Output Parameter:
1430: . has    - Flag for attribute existence

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

1436:   Level: advanced

1438: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1439: @*/
1440: PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has)
1441: {

1449:   PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);
1450:   PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);
1451:   return(0);
1452: }

1454: static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
1455: {
1456:   hid_t          h5;
1457:   htri_t         hhas;

1461:   PetscViewerHDF5GetFileId(viewer, &h5);
1462:   PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT));
1463:   *has = hhas ? PETSC_TRUE : PETSC_FALSE;
1464:   return(0);
1465: }

1467: /*
1468:   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
1469:   is attached to a communicator, in this case the attribute is a PetscViewer.
1470: */
1471: PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;

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

1476:   Collective

1478:   Input Parameter:
1479: . comm - the MPI communicator to share the HDF5 PetscViewer

1481:   Level: intermediate

1483:   Options Database Keys:
1484: . -viewer_hdf5_filename <name> - name of the HDF5 file

1486:   Environmental variables:
1487: . PETSC_VIEWER_HDF5_FILENAME - name of the HDF5 file

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

1494: .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy()
1495: @*/
1496: PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
1497: {
1499:   PetscBool      flg;
1500:   PetscViewer    viewer;
1501:   char           fname[PETSC_MAX_PATH_LEN];
1502:   MPI_Comm       ncomm;

1505:   PetscCommDuplicate(comm,&ncomm,NULL);if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1506:   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
1507:     MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,NULL);
1508:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1509:   }
1510:   MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg);
1511:   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1512:   if (!flg) { /* PetscViewer not yet created */
1513:     PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg);
1514:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1515:     if (!flg) {
1516:       PetscStrcpy(fname,"output.h5");
1517:       if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1518:     }
1519:     PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer);
1520:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1521:     PetscObjectRegisterDestroy((PetscObject)viewer);
1522:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1523:     MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer);
1524:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1525:   }
1526:   PetscCommDestroy(&ncomm);
1527:   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1528:   PetscFunctionReturn(viewer);
1529: }