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] = "";

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

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

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

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

 43:   PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");
 44:   PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);
 45:   PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL);
 46:   PetscOptionsBool("-viewer_hdf5_collective","Enable collective transfer mode","PetscViewerHDF5SetCollective",flg,&flg,&set);
 47:   if (set) PetscViewerHDF5SetCollective(v,flg);
 48:   flg  = PETSC_FALSE;
 49:   PetscOptionsBool("-viewer_hdf5_default_timestepping","Set default timestepping state","PetscViewerHDF5SetDefaultTimestepping",flg,&flg,&set);
 50:   if (set) PetscViewerHDF5SetDefaultTimestepping(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;

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

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

 75:   PetscFree(hdf5->filename);
 76:   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
 77:   return 0;
 78: }

 80: static PetscErrorCode PetscViewerFlush_HDF5(PetscViewer viewer)
 81: {
 82:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data;

 84:   if (hdf5->file_id) PetscStackCallHDF5(H5Fflush,(hdf5->file_id, H5F_SCOPE_LOCAL));
 85:   return 0;
 86: }

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

 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:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5GetDefaultTimestepping_C",NULL);
108:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetDefaultTimestepping_C",NULL);
109:   return 0;
110: }

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

116:   hdf5->btype = type;
117:   return 0;
118: }

120: static PetscErrorCode  PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type)
121: {
122:   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;

132:   hdf5->basedimension2 = flg;
133:   return 0;
134: }

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

140:     Logically Collective on PetscViewer

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

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

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

153:   Level: intermediate

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

157: @*/
158: PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg)
159: {
161:   PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));
162:   return 0;
163: }

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

169:     Logically Collective on PetscViewer

171:   Input Parameter:
172: .  viewer - the PetscViewer, must be of type HDF5

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

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

181:   Level: intermediate

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

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

191:   *flg = hdf5->basedimension2;
192:   return 0;
193: }

195: static PetscErrorCode  PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg)
196: {
197:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

199:   hdf5->spoutput = flg;
200:   return 0;
201: }

203: /*@
204:      PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is
205:        compiled with double precision PetscReal.

207:     Logically Collective on PetscViewer

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

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

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

220:   Level: intermediate

222: .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
223:           PetscReal

225: @*/
226: PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg)
227: {
229:   PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));
230:   return 0;
231: }

233: /*@
234:      PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is
235:        compiled with double precision PetscReal.

237:     Logically Collective on PetscViewer

239:   Input Parameter:
240: .  viewer - the PetscViewer, must be of type HDF5

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

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

249:   Level: intermediate

251: .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
252:           PetscReal

254: @*/
255: PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg)
256: {
257:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

260:   *flg = hdf5->spoutput;
261:   return 0;
262: }

264: static PetscErrorCode  PetscViewerHDF5SetCollective_HDF5(PetscViewer viewer, PetscBool flg)
265: {
266:   /* H5FD_MPIO_COLLECTIVE is wrong in hdf5 1.10.2, and is the same as H5FD_MPIO_INDEPENDENT in earlier versions
267:      - see e.g. https://gitlab.cosma.dur.ac.uk/swift/swiftsim/issues/431 */
268: #if H5_VERSION_GE(1,10,3) && defined(H5_HAVE_PARALLEL)
269:   {
270:     PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
271:     PetscStackCallHDF5(H5Pset_dxpl_mpio,(hdf5->dxpl_id, flg ? H5FD_MPIO_COLLECTIVE : H5FD_MPIO_INDEPENDENT));
272:   }
273: #else
274:   if (flg) 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");
275: #endif
276:   return 0;
277: }

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

282:   Logically Collective; flg must contain common value

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

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

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

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

300:   Level: intermediate

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

304: @*/
305: PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer viewer,PetscBool flg)
306: {
309:   PetscTryMethod(viewer,"PetscViewerHDF5SetCollective_C",(PetscViewer,PetscBool),(viewer,flg));
310:   return 0;
311: }

313: static PetscErrorCode  PetscViewerHDF5GetCollective_HDF5(PetscViewer viewer, PetscBool *flg)
314: {
315: #if defined(H5_HAVE_PARALLEL)
316:   PetscViewer_HDF5  *hdf5 = (PetscViewer_HDF5*) viewer->data;
317:   H5FD_mpio_xfer_t  mode;
318: #endif

320: #if !defined(H5_HAVE_PARALLEL)
321:   *flg = PETSC_FALSE;
322: #else
323:   PetscStackCallHDF5(H5Pget_dxpl_mpio,(hdf5->dxpl_id, &mode));
324:   *flg = (mode == H5FD_MPIO_COLLECTIVE) ? PETSC_TRUE : PETSC_FALSE;
325: #endif
326:   return 0;
327: }

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

332:   Not Collective

334:   Input Parameters:
335: . viewer - the HDF5 PetscViewer

337:   Output Parameters:
338: . flg - the flag

340:   Level: intermediate

342:   Notes:
343:   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.
344:   For more details, see PetscViewerHDF5SetCollective().

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

348: @*/
349: PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer viewer,PetscBool *flg)
350: {

354:   PetscUseMethod(viewer,"PetscViewerHDF5GetCollective_C",(PetscViewer,PetscBool*),(viewer,flg));
355:   return 0;
356: }

358: static PetscErrorCode  PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[])
359: {
360:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
361:   hid_t             plist_id;

363:   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
364:   if (hdf5->filename) PetscFree(hdf5->filename);
365:   PetscStrallocpy(name, &hdf5->filename);
366:   /* Set up file access property list with parallel I/O access */
367:   PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS));
368: #if defined(H5_HAVE_PARALLEL)
369:   PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), MPI_INFO_NULL));
370: #endif
371:   /* Create or open the file collectively */
372:   switch (hdf5->btype) {
373:   case FILE_MODE_READ:
374:     if (PetscDefined(USE_DEBUG)) {
375:       PetscMPIInt rank;
376:       PetscBool   flg;

378:       MPI_Comm_rank(PetscObjectComm((PetscObject)viewer),&rank);
379:       if (rank == 0) {
380:         PetscTestFile(hdf5->filename, 'r', &flg);
382:       }
383:       MPI_Barrier(PetscObjectComm((PetscObject)viewer));
384:     }
385:     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id));
386:     break;
387:   case FILE_MODE_APPEND:
388:   case FILE_MODE_UPDATE:
389:   {
390:     PetscBool flg;
391:     PetscTestFile(hdf5->filename, 'r', &flg);
392:     if (flg) PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id));
393:     else     PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_EXCL, H5P_DEFAULT, plist_id));
394:     break;
395:   }
396:   case FILE_MODE_WRITE:
397:     PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id));
398:     break;
399:   case FILE_MODE_UNDEFINED:
400:     SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
401:   default:
402:     SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP, "Unsupported file mode %s",PetscFileModes[hdf5->btype]);
403:   }
405:   PetscStackCallHDF5(H5Pclose,(plist_id));
406:   return 0;
407: }

409: static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name)
410: {
411:   PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data;

413:   *name = vhdf5->filename;
414:   return 0;
415: }

417: static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer)
418: {
419:   /*
420:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
421:   */

423:   return 0;
424: }

426: static PetscErrorCode PetscViewerHDF5SetDefaultTimestepping_HDF5(PetscViewer viewer, PetscBool flg)
427: {
428:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

430:   hdf5->defTimestepping = flg;
431:   return 0;
432: }

434: static PetscErrorCode PetscViewerHDF5GetDefaultTimestepping_HDF5(PetscViewer viewer, PetscBool *flg)
435: {
436:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

438:   *flg = hdf5->defTimestepping;
439:   return 0;
440: }

442: /*@
443:   PetscViewerHDF5SetDefaultTimestepping - Set the flag for default timestepping

445:   Logically Collective on PetscViewer

447:   Input Parameters:
448: + viewer - the PetscViewer; if it is not hdf5 then this command is ignored
449: - flg    - if PETSC_TRUE we will assume that timestepping is on

451:   Options Database:
452: . -viewer_hdf5_default_timestepping - turns on (true) or off (false) default timestepping

454:   Notes:
455:   If the timestepping attribute is not found for an object, then the default timestepping is used

457:   Level: intermediate

459: .seealso: PetscViewerHDF5GetDefaultTimestepping(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5GetTimestep()
460: @*/
461: PetscErrorCode PetscViewerHDF5SetDefaultTimestepping(PetscViewer viewer, PetscBool flg)
462: {
464:   PetscTryMethod(viewer, "PetscViewerHDF5SetDefaultTimestepping_C", (PetscViewer,PetscBool), (viewer,flg));
465:   return 0;
466: }

468: /*@
469:   PetscViewerHDF5GetDefaultTimestepping - Get the flag for default timestepping

471:   Not collective

473:   Input Parameter:
474: . viewer - the PetscViewer

476:   Output Parameter:
477: . flg    - if PETSC_TRUE we will assume that timestepping is on

479:   Notes:
480:   If the timestepping attribute is not found for an object, then the default timestepping is used

482:   Level: intermediate

484: .seealso: PetscViewerHDF5SetDefaultTimestepping(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5GetTimestep()
485: @*/
486: PetscErrorCode PetscViewerHDF5GetDefaultTimestepping(PetscViewer viewer, PetscBool *flg)
487: {
489:   PetscUseMethod(viewer, "PetscViewerHDF5GetDefaultTimestepping_C", (PetscViewer,PetscBool*), (viewer,flg));
490:   return 0;
491: }

493: /*MC
494:    PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file

496: .seealso:  PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET,
497:            PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING,
498:            PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB,
499:            PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType()

501:   Level: beginner
502: M*/

504: PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v)
505: {
506:   PetscViewer_HDF5 *hdf5;

508: #if !defined(H5_HAVE_PARALLEL)
509:   {
510:     PetscMPIInt size;
511:     MPI_Comm_size(PetscObjectComm((PetscObject)v), &size);
513:   }
514: #endif

516:   PetscNewLog(v,&hdf5);

518:   v->data                = (void*) hdf5;
519:   v->ops->destroy        = PetscViewerDestroy_HDF5;
520:   v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5;
521:   v->ops->setup          = PetscViewerSetUp_HDF5;
522:   v->ops->view           = PetscViewerView_HDF5;
523:   v->ops->flush          = PetscViewerFlush_HDF5;
524:   hdf5->btype            = FILE_MODE_UNDEFINED;
525:   hdf5->filename         = NULL;
526:   hdf5->timestep         = -1;
527:   hdf5->groups           = NULL;

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

531:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);
532:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);
533:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);
534:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_HDF5);
535:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);
536:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);
537:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetCollective_C",PetscViewerHDF5SetCollective_HDF5);
538:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetCollective_C",PetscViewerHDF5GetCollective_HDF5);
539:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetDefaultTimestepping_C",PetscViewerHDF5GetDefaultTimestepping_HDF5);
540:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetDefaultTimestepping_C",PetscViewerHDF5SetDefaultTimestepping_HDF5);
541:   return 0;
542: }

544: /*@C
545:    PetscViewerHDF5Open - Opens a file for HDF5 input/output.

547:    Collective

549:    Input Parameters:
550: +  comm - MPI communicator
551: .  name - name of file
552: -  type - type of file

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

557:   Options Database:
558: +  -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
559: -  -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal

561:    Level: beginner

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

570:    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.

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

574: .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(),
575:           PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(),
576:           MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName()
577: @*/
578: PetscErrorCode  PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
579: {
580:   PetscViewerCreate(comm, hdf5v);
581:   PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);
582:   PetscViewerFileSetMode(*hdf5v, type);
583:   PetscViewerFileSetName(*hdf5v, name);
584:   PetscViewerSetFromOptions(*hdf5v);
585:   return 0;
586: }

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

591:   Not collective

593:   Input Parameter:
594: . viewer - the PetscViewer

596:   Output Parameter:
597: . file_id - The file id

599:   Level: intermediate

601: .seealso: PetscViewerHDF5Open()
602: @*/
603: PetscErrorCode  PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
604: {
605:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

608:   if (file_id) *file_id = hdf5->file_id;
609:   return 0;
610: }

612: /*@C
613:   PetscViewerHDF5PushGroup - Set the current HDF5 group for output

615:   Not collective

617:   Input Parameters:
618: + viewer - the PetscViewer
619: - name - The group name

621:   Level: intermediate

623:   Notes:
624:   This is designed to mnemonically resemble the Unix cd command.
625:   + 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.
626:   . NULL, empty string, or any sequence of all slashes (e.g. "///") is interpreted as the root group "/".
627:   - "." means the current group is pushed again.

629:   Example:
630:   Suppose the current group is "/a".
631:   + If name is NULL, empty string, or a sequence of all slashes (e.g. "///"), then the new group will be "/".
632:   . If name is ".", then the new group will be "/a".
633:   . If name is "b", then the new group will be "/a/b".
634:   - If name is "/b", then the new group will be "/b".

636:   Developer Notes:
637:   The root group "/" is internally stored as NULL.

639: .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
640: @*/
641: PetscErrorCode  PetscViewerHDF5PushGroup(PetscViewer viewer, const char name[])
642: {
643:   PetscViewer_HDF5          *hdf5 = (PetscViewer_HDF5*) viewer->data;
644:   PetscViewerHDF5GroupList  *groupNode;
645:   size_t                    i,len;
646:   char                      buf[PETSC_MAX_PATH_LEN];
647:   const char                *gname;

651:   PetscStrlen(name, &len);
652:   gname = NULL;
653:   if (len) {
654:     if (len == 1 && name[0] == '.') {
655:       /* use current name */
656:       gname = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : NULL;
657:     } else if (name[0] == '/') {
658:       /* absolute */
659:       for (i=1; i<len; i++) {
660:         if (name[i] != '/') {
661:           gname = name;
662:           break;
663:         }
664:       }
665:     } else {
666:       /* relative */
667:       const char *parent = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : "";
668:       PetscSNPrintf(buf, sizeof(buf), "%s/%s", parent, name);
669:       gname = buf;
670:     }
671:   }
672:   PetscNew(&groupNode);
673:   PetscStrallocpy(gname, (char**) &groupNode->name);
674:   groupNode->next = hdf5->groups;
675:   hdf5->groups    = groupNode;
676:   return 0;
677: }

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

682:   Not collective

684:   Input Parameter:
685: . viewer - the PetscViewer

687:   Level: intermediate

689: .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
690: @*/
691: PetscErrorCode  PetscViewerHDF5PopGroup(PetscViewer viewer)
692: {
693:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
694:   PetscViewerHDF5GroupList *groupNode;

698:   groupNode    = hdf5->groups;
699:   hdf5->groups = hdf5->groups->next;
700:   PetscFree(groupNode->name);
701:   PetscFree(groupNode);
702:   return 0;
703: }

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

709:   Not collective

711:   Input Parameter:
712: . viewer - the PetscViewer

714:   Output Parameter:
715: . name - The group name

717:   Level: intermediate

719: .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup()
720: @*/
721: PetscErrorCode  PetscViewerHDF5GetGroup(PetscViewer viewer, const char *name[])
722: {
723:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data;

727:   if (hdf5->groups) *name = hdf5->groups->name;
728:   else *name = NULL;
729:   return 0;
730: }

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

737:   Not collective

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

742:   Output Parameters:
743: + fileId - The HDF5 file ID
744: - groupId - The HDF5 group ID

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

749:   Level: intermediate

751: .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
752: @*/
753: PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId)
754: {
755:   hid_t          file_id;
756:   H5O_type_t     type;
757:   const char     *groupName = NULL, *fileName = NULL;
758:   PetscBool      writable, has;

760:   PetscViewerWritable(viewer, &writable);
761:   PetscViewerHDF5GetFileId(viewer, &file_id);
762:   PetscViewerFileGetName(viewer, &fileName);
763:   PetscViewerHDF5GetGroup(viewer, &groupName);
764:   PetscViewerHDF5Traverse_Internal(viewer, groupName, writable, &has, &type);
765:   if (!has) {
767:     else           SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_LIB, "HDF5 failed to create group %s although file %s is open for writing", groupName, fileName);
768:   }
770:   PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT));
771:   *fileId  = file_id;
772:   return 0;
773: }

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

778:   Not collective

780:   Input Parameter:
781: . viewer - the PetscViewer

783:   Level: intermediate

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

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

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

798: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PopTimestepping(), PetscViewerHDF5IsTimestepping(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
799: @*/
800: PetscErrorCode  PetscViewerHDF5PushTimestepping(PetscViewer viewer)
801: {
802:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

806:   hdf5->timestepping = PETSC_TRUE;
807:   if (hdf5->timestep < 0) hdf5->timestep = 0;
808:   return 0;
809: }

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

814:   Not collective

816:   Input Parameter:
817: . viewer - the PetscViewer

819:   Level: intermediate

821:   Notes:
822:   See PetscViewerHDF5PushTimestepping() for details.

824: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5IsTimestepping(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
825: @*/
826: PetscErrorCode  PetscViewerHDF5PopTimestepping(PetscViewer viewer)
827: {
828:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

832:   hdf5->timestepping = PETSC_FALSE;
833:   return 0;
834: }

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

839:   Not collective

841:   Input Parameter:
842: . viewer - the PetscViewer

844:   Output Parameter:
845: . flg - is timestepping active?

847:   Level: intermediate

849:   Notes:
850:   See PetscViewerHDF5PushTimestepping() for details.

852: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5PopTimestepping(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
853: @*/
854: PetscErrorCode  PetscViewerHDF5IsTimestepping(PetscViewer viewer, PetscBool *flg)
855: {
856:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

859:   *flg = hdf5->timestepping;
860:   return 0;
861: }

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

866:   Not collective

868:   Input Parameter:
869: . viewer - the PetscViewer

871:   Level: intermediate

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

876: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep()
877: @*/
878: PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
879: {
880:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

884:   ++hdf5->timestep;
885:   return 0;
886: }

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

891:   Logically collective

893:   Input Parameters:
894: + viewer - the PetscViewer
895: - timestep - The timestep

897:   Level: intermediate

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

902: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
903: @*/
904: PetscErrorCode  PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
905: {
906:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

912:   hdf5->timestep = timestep;
913:   return 0;
914: }

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

919:   Not collective

921:   Input Parameter:
922: . viewer - the PetscViewer

924:   Output Parameter:
925: . timestep - The timestep

927:   Level: intermediate

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

932: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep()
933: @*/
934: PetscErrorCode  PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
935: {
936:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

941:   *timestep = hdf5->timestep;
942:   return 0;
943: }

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

948:   Not collective

950:   Input Parameter:
951: . ptype - the PETSc datatype name (for example PETSC_DOUBLE)

953:   Output Parameter:
954: . mtype - the MPI datatype (for example MPI_DOUBLE, ...)

956:   Level: advanced

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

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

984:   Not collective

986:   Input Parameter:
987: . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...)

989:   Output Parameter:
990: . ptype - the PETSc datatype name (for example PETSC_DOUBLE)

992:   Level: advanced

994: .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
995: @*/
996: PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
997: {
998: #if defined(PETSC_USE_64BIT_INDICES)
999:   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_LONG;
1000:   else if (htype == H5T_NATIVE_LLONG)  *ptype = PETSC_INT;
1001: #else
1002:   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_INT;
1003: #endif
1004:   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
1005:   else if (htype == H5T_NATIVE_LONG)   *ptype = PETSC_LONG;
1006:   else if (htype == H5T_NATIVE_SHORT)  *ptype = PETSC_SHORT;
1007:   else if (htype == H5T_NATIVE_FLOAT)  *ptype = PETSC_FLOAT;
1008:   else if (htype == H5T_NATIVE_CHAR)   *ptype = PETSC_CHAR;
1009:   else if (htype == H5T_NATIVE_UCHAR)  *ptype = PETSC_CHAR;
1010:   else if (htype == H5T_C_S1)          *ptype = PETSC_STRING;
1011:   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
1012:   return 0;
1013: }

1015: /*@C
1016:  PetscViewerHDF5WriteAttribute - Write an attribute

1018:   Collective

1020:   Input Parameters:
1021: + viewer - The HDF5 viewer
1022: . parent - The parent dataset/group name
1023: . name   - The attribute name
1024: . datatype - The attribute type
1025: - value    - The attribute value

1027:   Level: advanced

1029:   Notes:
1030:   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.

1032: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1033: @*/
1034: PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value)
1035: {
1036:   char           *parentAbsPath;
1037:   hid_t          h5, dataspace, obj, attribute, dtype;
1038:   PetscBool      has;

1045:   PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath);
1046:   PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_TRUE, NULL, NULL);
1047:   PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has);
1048:   PetscDataTypeToHDF5DataType(datatype, &dtype);
1049:   if (datatype == PETSC_STRING) {
1050:     size_t len;
1051:     PetscStrlen((const char *) value, &len);
1052:     PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
1053:   }
1054:   PetscViewerHDF5GetFileId(viewer, &h5);
1055:   PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR));
1056:   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parentAbsPath, H5P_DEFAULT));
1057:   if (has) {
1058:     PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
1059:   } else {
1060:     PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
1061:   }
1062:   PetscStackCallHDF5(H5Awrite,(attribute, dtype, value));
1063:   if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype));
1064:   PetscStackCallHDF5(H5Aclose,(attribute));
1065:   PetscStackCallHDF5(H5Oclose,(obj));
1066:   PetscStackCallHDF5(H5Sclose,(dataspace));
1067:   PetscFree(parentAbsPath);
1068:   return 0;
1069: }

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

1074:   Collective

1076:   Input Parameters:
1077: + viewer   - The HDF5 viewer
1078: . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
1079: . name     - The attribute name
1080: . datatype - The attribute type
1081: - value    - The attribute value

1083:   Notes:
1084:   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).
1085:   You might want to check first if it does using PetscViewerHDF5HasObject().

1087:   Level: advanced

1089: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1090: @*/
1091: PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value)
1092: {
1097:   PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);
1098:   PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);
1099:   return 0;
1100: }

1102: /*@C
1103:  PetscViewerHDF5ReadAttribute - Read an attribute

1105:   Collective

1107:   Input Parameters:
1108: + viewer - The HDF5 viewer
1109: . parent - The parent dataset/group name
1110: . name   - The attribute name
1111: . datatype - The attribute type
1112: - defaultValue - The pointer to the default value

1114:   Output Parameter:
1115: . value    - The pointer to the read HDF5 attribute value

1117:   Notes:
1118:   If defaultValue is NULL and the attribute is not found, an error occurs.
1119:   If defaultValue is not NULL and the attribute is not found, defaultValue is copied to value.
1120:   The pointers defaultValue and value can be the same; for instance
1121: $  flg = PETSC_FALSE;
1122: $  PetscViewerHDF5ReadAttribute(viewer,name,"attr",PETSC_BOOL,&flg,&flg);
1123:   is valid, but make sure the default value is initialized.

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

1127:   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.

1129:   Level: advanced

1131: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1132: @*/
1133: PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *defaultValue, void *value)
1134: {
1135:   char           *parentAbsPath;
1136:   hid_t          h5, obj, attribute, dtype;
1137:   PetscBool      has;

1144:   PetscDataTypeToHDF5DataType(datatype, &dtype);
1145:   PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath);
1146:   PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, &has, NULL);
1147:   if (has) PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has);
1148:   if (!has) {
1149:     if (defaultValue) {
1150:       if (defaultValue != value) {
1151:         if (datatype == PETSC_STRING) {
1152:           PetscStrallocpy(*(char**)defaultValue, (char**)value);
1153:         } else {
1154:           size_t len;
1155:           PetscStackCallHDF5ReturnNoCheck(len,H5Tget_size,(dtype));
1156:           PetscMemcpy(value, defaultValue, len);
1157:         }
1158:       }
1159:       PetscFree(parentAbsPath);
1160:       return 0;
1161:     } else SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist and default value not provided", parentAbsPath, name);
1162:   }
1163:   PetscViewerHDF5GetFileId(viewer, &h5);
1164:   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parentAbsPath, H5P_DEFAULT));
1165:   PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
1166:   if (datatype == PETSC_STRING) {
1167:     size_t len;
1168:     hid_t  atype;
1169:     PetscStackCallHDF5Return(atype,H5Aget_type,(attribute));
1170:     PetscStackCallHDF5ReturnNoCheck(len,H5Tget_size,(atype));
1171:     PetscMalloc((len+1) * sizeof(char), value);
1172:     PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
1173:     PetscStackCallHDF5(H5Aread,(attribute, dtype, *(char**)value));
1174:   } else {
1175:     PetscStackCallHDF5(H5Aread,(attribute, dtype, value));
1176:   }
1177:   PetscStackCallHDF5(H5Aclose,(attribute));
1178:   /* H5Oclose can be used to close groups, datasets, or committed datatypes */
1179:   PetscStackCallHDF5(H5Oclose,(obj));
1180:   PetscFree(parentAbsPath);
1181:   return 0;
1182: }

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

1187:   Collective

1189:   Input Parameters:
1190: + viewer   - The HDF5 viewer
1191: . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
1192: . name     - The attribute name
1193: - datatype - The attribute type

1195:   Output Parameter:
1196: . value    - The attribute value

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

1202:   Level: advanced

1204: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute() PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1205: @*/
1206: PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *defaultValue, void *value)
1207: {
1212:   PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);
1213:   PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, defaultValue, value);
1214:   return 0;
1215: }

1217: static inline PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_)
1218: {
1219:   htri_t exists;
1220:   hid_t group;

1222:   PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT));
1223:   if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT));
1224:   if (!exists && createGroup) {
1225:     PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT));
1226:     PetscStackCallHDF5(H5Gclose,(group));
1227:     exists = PETSC_TRUE;
1228:   }
1229:   *exists_ = (PetscBool) exists;
1230:   return 0;
1231: }

1233: static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype)
1234: {
1235:   const char     rootGroupName[] = "/";
1236:   hid_t          h5;
1237:   PetscBool      exists=PETSC_FALSE;
1238:   PetscInt       i;
1239:   int            n;
1240:   char           **hierarchy;
1241:   char           buf[PETSC_MAX_PATH_LEN]="";

1245:   else name = rootGroupName;
1246:   if (has) {
1248:     *has = PETSC_FALSE;
1249:   }
1250:   if (otype) {
1252:     *otype = H5O_TYPE_UNKNOWN;
1253:   }
1254:   PetscViewerHDF5GetFileId(viewer, &h5);

1256:   /*
1257:      Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing.
1258:      Hence, each of them needs to be tested separately:
1259:      1) whether it's a valid link
1260:      2) whether this link resolves to an object
1261:      See H5Oexists_by_name() documentation.
1262:   */
1263:   PetscStrToArray(name,'/',&n,&hierarchy);
1264:   if (!n) {
1265:     /*  Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */
1266:     if (has)   *has   = PETSC_TRUE;
1267:     if (otype) *otype = H5O_TYPE_GROUP;
1268:     PetscStrToArrayDestroy(n,hierarchy);
1269:     return 0;
1270:   }
1271:   for (i=0; i<n; i++) {
1272:     PetscStrcat(buf,"/");
1273:     PetscStrcat(buf,hierarchy[i]);
1274:     PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);
1275:     if (!exists) break;
1276:   }
1277:   PetscStrToArrayDestroy(n,hierarchy);

1279:   /* If the object exists, get its type */
1280:   if (exists && otype) {
1281:     H5O_info_t info;

1283:     /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */
1284:     PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT));
1285:     *otype = info.type;
1286:   }
1287:   if (has) *has = exists;
1288:   return 0;
1289: }

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

1294:   Collective

1296:   Input Parameters:
1297: + viewer - The HDF5 viewer
1298: - path - The group path

1300:   Output Parameter:
1301: . has - Flag for group existence

1303:   Level: advanced

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

1310: .seealso: PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasDataset(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5GetGroup(), PetscViewerHDF5OpenGroup()
1311: @*/
1312: PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, const char path[], PetscBool *has)
1313: {
1314:   H5O_type_t     type;
1315:   char           *abspath;

1320:   PetscViewerHDF5GetAbsolutePath_Internal(viewer, path, &abspath);
1321:   PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type);
1322:   *has = (PetscBool)(type == H5O_TYPE_GROUP);
1323:   PetscFree(abspath);
1324:   return 0;
1325: }

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

1330:   Collective

1332:   Input Parameters:
1333: + viewer - The HDF5 viewer
1334: - path - The dataset path

1336:   Output Parameter:
1337: . has - Flag whether dataset exists

1339:   Level: advanced

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

1346: .seealso: PetscViewerHDF5HasObject(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasGroup(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5GetGroup()
1347: @*/
1348: PetscErrorCode PetscViewerHDF5HasDataset(PetscViewer viewer, const char path[], PetscBool *has)
1349: {
1350:   H5O_type_t     type;
1351:   char           *abspath;

1356:   PetscViewerHDF5GetAbsolutePath_Internal(viewer, path, &abspath);
1357:   PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type);
1358:   *has = (PetscBool)(type == H5O_TYPE_DATASET);
1359:   PetscFree(abspath);
1360:   return 0;
1361: }

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

1366:   Collective

1368:   Input Parameters:
1369: + viewer - The HDF5 viewer
1370: - obj    - The named object

1372:   Output Parameter:
1373: . has    - Flag for dataset existence

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

1379:   Level: advanced

1381: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasDataset(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5GetGroup()
1382: @*/
1383: PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has)
1384: {
1385:   size_t         len;

1390:   PetscStrlen(obj->name, &len);
1392:   PetscViewerHDF5HasDataset(viewer, obj->name, has);
1393:   return 0;
1394: }

1396: /*@C
1397:  PetscViewerHDF5HasAttribute - Check whether an attribute exists

1399:   Collective

1401:   Input Parameters:
1402: + viewer - The HDF5 viewer
1403: . parent - The parent dataset/group name
1404: - name   - The attribute name

1406:   Output Parameter:
1407: . has    - Flag for attribute existence

1409:   Level: advanced

1411:   Notes:
1412:   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.

1414: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1415: @*/
1416: PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
1417: {
1418:   char           *parentAbsPath;

1424:   PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath);
1425:   PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, has, NULL);
1426:   if (*has) PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, has);
1427:   PetscFree(parentAbsPath);
1428:   return 0;
1429: }

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

1434:   Collective

1436:   Input Parameters:
1437: + viewer - The HDF5 viewer
1438: . obj    - The object whose name is used to lookup the parent dataset, relative to the current group.
1439: - name   - The attribute name

1441:   Output Parameter:
1442: . has    - Flag for attribute existence

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

1448:   Level: advanced

1450: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1451: @*/
1452: PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has)
1453: {
1458:   PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);
1459:   PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);
1460:   return 0;
1461: }

1463: static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
1464: {
1465:   hid_t          h5;
1466:   htri_t         hhas;

1468:   PetscViewerHDF5GetFileId(viewer, &h5);
1469:   PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT));
1470:   *has = hhas ? PETSC_TRUE : PETSC_FALSE;
1471:   return 0;
1472: }

1474: /*
1475:   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
1476:   is attached to a communicator, in this case the attribute is a PetscViewer.
1477: */
1478: PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;

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

1483:   Collective

1485:   Input Parameter:
1486: . comm - the MPI communicator to share the HDF5 PetscViewer

1488:   Level: intermediate

1490:   Options Database Keys:
1491: . -viewer_hdf5_filename <name> - name of the HDF5 file

1493:   Environmental variables:
1494: . PETSC_VIEWER_HDF5_FILENAME - name of the HDF5 file

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

1501: .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy()
1502: @*/
1503: PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
1504: {
1506:   PetscBool      flg;
1507:   PetscViewer    viewer;
1508:   char           fname[PETSC_MAX_PATH_LEN];
1509:   MPI_Comm       ncomm;

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