Actual source code: hdf5v.c

petsc-3.7.5 2017-01-01
Report Typos and Errors
  1: #include <petsc/private/viewerimpl.h>    /*I   "petscsys.h"   I*/
  2: #include <petscviewerhdf5.h>    /*I   "petscviewerhdf5.h"   I*/

  4: typedef struct GroupList {
  5:   const char       *name;
  6:   struct GroupList *next;
  7: } GroupList;

  9: typedef struct {
 10:   char          *filename;
 11:   PetscFileMode btype;
 12:   hid_t         file_id;
 13:   PetscInt      timestep;
 14:   GroupList     *groups;
 15:   PetscBool     basedimension2;  /* save vectors and DMDA vectors with a dimension of at least 2 even if the bs/dof is 1 */
 16:   PetscBool     spoutput;  /* write data in single precision even if PETSc is compiled with double precision PetscReal */
 17: } PetscViewer_HDF5;

 21: static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptionItems *PetscOptionsObject,PetscViewer v)
 22: {
 23:   PetscErrorCode   ierr;
 24:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data;

 27:   PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");
 28:   PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);
 29:   PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL);
 30:   PetscOptionsTail();
 31:   return(0);
 32: }

 36: static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer)
 37: {
 38:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data;
 39:   PetscErrorCode   ierr;

 42:   PetscFree(hdf5->filename);
 43:   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
 44:   return(0);
 45: }

 49: PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer)
 50: {
 51:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
 52:   PetscErrorCode   ierr;

 55:   PetscViewerFileClose_HDF5(viewer);
 56:   while (hdf5->groups) {
 57:     GroupList *tmp = hdf5->groups->next;

 59:     PetscFree(hdf5->groups->name);
 60:     PetscFree(hdf5->groups);
 61:     hdf5->groups = tmp;
 62:   }
 63:   PetscFree(hdf5);
 64:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);
 65:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);
 66:   return(0);
 67: }

 71: PetscErrorCode  PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type)
 72: {
 73:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

 77:   hdf5->btype = type;
 78:   return(0);
 79: }

 83: PetscErrorCode  PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg)
 84: {
 85:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

 88:   hdf5->basedimension2 = flg;
 89:   return(0);
 90: }

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

 98:     Logically Collective on PetscViewer

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

104:   Options Database:
105: .  -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


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

111:   Level: intermediate

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

115: @*/
116: PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg)
117: {

122:   PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));
123:   return(0);
124: }

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

132:     Logically Collective on PetscViewer

134:   Input Parameter:
135: .  viewer - the PetscViewer, must be of type HDF5

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

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

143:   Level: intermediate

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

147: @*/
148: PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg)
149: {
150:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

154:   *flg = hdf5->basedimension2;
155:   return(0);
156: }

160: PetscErrorCode  PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg)
161: {
162:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

165:   hdf5->spoutput = flg;
166:   return(0);
167: }

171: /*@C
172:      PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is
173:        compiled with double precision PetscReal.

175:     Logically Collective on PetscViewer

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

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


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

188:   Level: intermediate

190: .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
191:           PetscReal

193: @*/
194: PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg)
195: {

200:   PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));
201:   return(0);
202: }

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

210:     Logically Collective on PetscViewer

212:   Input Parameter:
213: .  viewer - the PetscViewer, must be of type HDF5

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

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

221:   Level: intermediate

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

226: @*/
227: PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg)
228: {
229:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

233:   *flg = hdf5->spoutput;
234:   return(0);
235: }

239: PetscErrorCode  PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[])
240: {
241:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
242: #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
243:   MPI_Info          info = MPI_INFO_NULL;
244: #endif
245:   hid_t             plist_id;
246:   PetscErrorCode    ierr;

249:   PetscStrallocpy(name, &hdf5->filename);
250:   /* Set up file access property list with parallel I/O access */
251:   PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS));
252: #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
253:   PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), info));
254: #endif
255:   /* Create or open the file collectively */
256:   switch (hdf5->btype) {
257:   case FILE_MODE_READ:
258:     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id));
259:     break;
260:   case FILE_MODE_APPEND:
261:     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id));
262:     break;
263:   case FILE_MODE_WRITE:
264:     PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id));
265:     break;
266:   default:
267:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
268:   }
269:   if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name);
270:   PetscStackCallHDF5(H5Pclose,(plist_id));
271:   return(0);
272: }

276: PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v)
277: {
278:   PetscViewer_HDF5 *hdf5;
279:   PetscErrorCode   ierr;

282:   PetscNewLog(v,&hdf5);

284:   v->data                = (void*) hdf5;
285:   v->ops->destroy        = PetscViewerDestroy_HDF5;
286:   v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5;
287:   v->ops->flush          = 0;
288:   hdf5->btype            = (PetscFileMode) -1;
289:   hdf5->filename         = 0;
290:   hdf5->timestep         = -1;
291:   hdf5->groups           = NULL;

293:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);
294:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);
295:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);
296:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);
297:   return(0);
298: }

302: /*@C
303:    PetscViewerHDF5Open - Opens a file for HDF5 input/output.

305:    Collective on MPI_Comm

307:    Input Parameters:
308: +  comm - MPI communicator
309: .  name - name of file
310: -  type - type of file
311: $    FILE_MODE_WRITE - create new file for binary output
312: $    FILE_MODE_READ - open existing file for binary input
313: $    FILE_MODE_APPEND - open existing file for binary output

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

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

322:    Level: beginner

324:    Note:
325:    This PetscViewer should be destroyed with PetscViewerDestroy().

327:    Concepts: HDF5 files
328:    Concepts: PetscViewerHDF5^creating

330: .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(),
331:           PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(),
332:           MatLoad(), PetscFileMode, PetscViewer
333: @*/
334: PetscErrorCode  PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
335: {

339:   PetscViewerCreate(comm, hdf5v);
340:   PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);
341:   PetscViewerFileSetMode(*hdf5v, type);
342:   PetscViewerFileSetName(*hdf5v, name);
343:   return(0);
344: }

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

351:   Not collective

353:   Input Parameter:
354: . viewer - the PetscViewer

356:   Output Parameter:
357: . file_id - The file id

359:   Level: intermediate

361: .seealso: PetscViewerHDF5Open()
362: @*/
363: PetscErrorCode  PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
364: {
365:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

369:   if (file_id) *file_id = hdf5->file_id;
370:   return(0);
371: }

375: /*@C
376:   PetscViewerHDF5PushGroup - Set the current HDF5 group for output

378:   Not collective

380:   Input Parameters:
381: + viewer - the PetscViewer
382: - name - The group name

384:   Level: intermediate

386: .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
387: @*/
388: PetscErrorCode  PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name)
389: {
390:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
391:   GroupList        *groupNode;
392:   PetscErrorCode   ierr;

397:   PetscMalloc(sizeof(GroupList), &groupNode);
398:   PetscStrallocpy(name, (char**) &groupNode->name);

400:   groupNode->next = hdf5->groups;
401:   hdf5->groups    = groupNode;
402:   return(0);
403: }

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

410:   Not collective

412:   Input Parameter:
413: . viewer - the PetscViewer

415:   Level: intermediate

417: .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup()
418: @*/
419: PetscErrorCode  PetscViewerHDF5PopGroup(PetscViewer viewer)
420: {
421:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
422:   GroupList        *groupNode;
423:   PetscErrorCode   ierr;

427:   if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
428:   groupNode    = hdf5->groups;
429:   hdf5->groups = hdf5->groups->next;
430:   PetscFree(groupNode->name);
431:   PetscFree(groupNode);
432:   return(0);
433: }

437: /*@C
438:   PetscViewerHDF5GetGroup - Get the current HDF5 group for output. If none has been assigned, returns NULL.

440:   Not collective

442:   Input Parameter:
443: . viewer - the PetscViewer

445:   Output Parameter:
446: . name - The group name

448:   Level: intermediate

450: .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup()
451: @*/
452: PetscErrorCode  PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name)
453: {
454:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data;

459:   if (hdf5->groups) *name = hdf5->groups->name;
460:   else *name = NULL;
461:   return(0);
462: }

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

469:   Not collective

471:   Input Parameter:
472: . viewer - the PetscViewer

474:   Level: intermediate

476: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep()
477: @*/
478: PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
479: {
480:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

484:   ++hdf5->timestep;
485:   return(0);
486: }

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

494:   Not collective

496:   Input Parameters:
497: + viewer - the PetscViewer
498: - timestep - The timestep number

500:   Level: intermediate

502: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
503: @*/
504: PetscErrorCode  PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
505: {
506:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

510:   hdf5->timestep = timestep;
511:   return(0);
512: }

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

519:   Not collective

521:   Input Parameter:
522: . viewer - the PetscViewer

524:   Output Parameter:
525: . timestep - The timestep number

527:   Level: intermediate

529: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep()
530: @*/
531: PetscErrorCode  PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
532: {
533:   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;

538:   *timestep = hdf5->timestep;
539:   return(0);
540: }

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

547:   Not collective

549:   Input Parameter:
550: . ptype - the PETSc datatype name (for example PETSC_DOUBLE)

552:   Output Parameter:
553: . mtype - the MPI datatype (for example MPI_DOUBLE, ...)

555:   Level: advanced

557: .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
558: @*/
559: PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
560: {
562:   if (ptype == PETSC_INT)
563: #if defined(PETSC_USE_64BIT_INDICES)
564:                                        *htype = H5T_NATIVE_LLONG;
565: #else
566:                                        *htype = H5T_NATIVE_INT;
567: #endif
568:   else if (ptype == PETSC_DOUBLE)      *htype = H5T_NATIVE_DOUBLE;
569:   else if (ptype == PETSC_LONG)        *htype = H5T_NATIVE_LONG;
570:   else if (ptype == PETSC_SHORT)       *htype = H5T_NATIVE_SHORT;
571:   else if (ptype == PETSC_ENUM)        *htype = H5T_NATIVE_DOUBLE;
572:   else if (ptype == PETSC_BOOL)        *htype = H5T_NATIVE_DOUBLE;
573:   else if (ptype == PETSC_FLOAT)       *htype = H5T_NATIVE_FLOAT;
574:   else if (ptype == PETSC_CHAR)        *htype = H5T_NATIVE_CHAR;
575:   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
576:   else if (ptype == PETSC_STRING)      *htype = H5Tcopy(H5T_C_S1);
577:   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
578:   return(0);
579: }

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

586:   Not collective

588:   Input Parameter:
589: . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...)

591:   Output Parameter:
592: . ptype - the PETSc datatype name (for example PETSC_DOUBLE)

594:   Level: advanced

596: .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
597: @*/
598: PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
599: {
601: #if defined(PETSC_USE_64BIT_INDICES)
602:   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_LONG;
603:   else if (htype == H5T_NATIVE_LLONG)  *ptype = PETSC_INT;
604: #else
605:   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_INT;
606: #endif
607:   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
608:   else if (htype == H5T_NATIVE_LONG)   *ptype = PETSC_LONG;
609:   else if (htype == H5T_NATIVE_SHORT)  *ptype = PETSC_SHORT;
610:   else if (htype == H5T_NATIVE_FLOAT)  *ptype = PETSC_FLOAT;
611:   else if (htype == H5T_NATIVE_CHAR)   *ptype = PETSC_CHAR;
612:   else if (htype == H5T_NATIVE_UCHAR)  *ptype = PETSC_CHAR;
613:   else if (htype == H5T_C_S1)          *ptype = PETSC_STRING;
614:   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
615:   return(0);
616: }

620: /*@
621:  PetscViewerHDF5WriteAttribute - Write a scalar attribute

623:   Input Parameters:
624: + viewer - The HDF5 viewer
625: . parent - The parent name
626: . name   - The attribute name
627: . datatype - The attribute type
628: - value    - The attribute value

630:   Level: advanced

632: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute()
633: @*/
634: PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value)
635: {
636:   hid_t          h5, dataspace, dataset, attribute, dtype;

644:   PetscDataTypeToHDF5DataType(datatype, &dtype);
645:   if (datatype == PETSC_STRING) {
646:     size_t len;
647:     PetscStrlen((const char *) value, &len);
648:     PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
649:   }
650:   PetscViewerHDF5GetFileId(viewer, &h5);
651:   PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR));
652: #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
653:   PetscStackCallHDF5Return(dataset,H5Dopen2,(h5, parent, H5P_DEFAULT));
654:   PetscStackCallHDF5Return(attribute,H5Acreate2,(dataset, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
655: #else
656:   PetscStackCallHDF5Return(dataset,H5Dopen,(h5, parent));
657:   PetscStackCallHDF5Return(attribute,H5Acreate,(dataset, name, dtype, dataspace, H5P_DEFAULT));
658: #endif
659:   PetscStackCallHDF5(H5Awrite,(attribute, dtype, value));
660:   if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype));
661:   PetscStackCallHDF5(H5Aclose,(attribute));
662:   PetscStackCallHDF5(H5Dclose,(dataset));
663:   PetscStackCallHDF5(H5Sclose,(dataspace));
664:   return(0);
665: }

669: /*@
670:  PetscViewerHDF5ReadAttribute - Read a scalar attribute

672:   Input Parameters:
673: + viewer - The HDF5 viewer
674: . parent - The parent name
675: . name   - The attribute name
676: - datatype - The attribute type

678:   Output Parameter:
679: . value    - The attribute value

681:   Level: advanced

683: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute()
684: @*/
685: PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, void *value)
686: {
687:   hid_t          h5, dataspace, dataset, attribute, dtype;

695:   PetscDataTypeToHDF5DataType(datatype, &dtype);
696:   PetscViewerHDF5GetFileId(viewer, &h5);
697:   PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR));
698: #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
699:   PetscStackCallHDF5Return(dataset,H5Dopen2,(h5, parent, H5P_DEFAULT));
700: #else
701:   PetscStackCallHDF5Return(dataset,H5Dopen,(h5, parent));
702: #endif
703:   PetscStackCallHDF5Return(attribute,H5Aopen_name,(dataset, name));
704:   PetscStackCallHDF5(H5Aread,(attribute, dtype, value));
705:   PetscStackCallHDF5(H5Aclose,(attribute));
706:   PetscStackCallHDF5(H5Dclose,(dataset));
707:   PetscStackCallHDF5(H5Sclose,(dataspace));
708:   return(0);
709: }

713: static PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, const char name[], H5O_type_t otype, PetscBool *has)
714: {
715:   hid_t          h5;

722:   *has = PETSC_FALSE;
723:   PetscViewerHDF5GetFileId(viewer, &h5);
724:   if (H5Lexists(h5, name, H5P_DEFAULT)) {
725:     H5O_info_t info;
726:     hid_t      obj;

728:     PetscStackCallHDF5Return(obj,H5Oopen,(h5, name, H5P_DEFAULT));
729:     PetscStackCallHDF5(H5Oget_info,(obj, &info));
730:     if (otype == info.type) *has = PETSC_TRUE;
731:     PetscStackCallHDF5(H5Oclose,(obj));
732:   }
733:   return(0);
734: }

738: /*@
739:  PetscViewerHDF5HasAttribute - Check whether a scalar attribute exists

741:   Input Parameters:
742: + viewer - The HDF5 viewer
743: . parent - The parent name
744: - name   - The attribute name

746:   Output Parameter:
747: . has    - Flag for attribute existence

749:   Level: advanced

751: .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute()
752: @*/
753: PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
754: {
755:   hid_t          h5, dataset;
756:   htri_t         hhas;
757:   PetscBool      exists;

765:   *has = PETSC_FALSE;
766:   PetscViewerHDF5GetFileId(viewer, &h5);
767:   PetscViewerHDF5HasObject(viewer, parent, H5O_TYPE_DATASET, &exists);
768:   if (exists) {
769: #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
770:     PetscStackCall("H5Dopen2",dataset = H5Dopen2(h5, parent, H5P_DEFAULT));
771: #else
772:     PetscStackCall("H5Dopen",dataset = H5Dopen(h5, parent));
773: #endif
774:     if (dataset < 0) return(0);
775:     PetscStackCall("H5Aexists",hhas = H5Aexists(dataset, name));
776:     if (hhas < 0) {
777:       PetscStackCallHDF5(H5Dclose,(dataset));
778:       return(0);
779:     }
780:     PetscStackCallHDF5(H5Dclose,(dataset));
781:     *has = hhas ? PETSC_TRUE : PETSC_FALSE;
782:   }
783:   return(0);
784: }

786: /*
787:   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
788:   is attached to a communicator, in this case the attribute is a PetscViewer.
789: */
790: static int Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;

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

797:   Collective on MPI_Comm

799:   Input Parameter:
800: . comm - the MPI communicator to share the HDF5 PetscViewer

802:   Level: intermediate

804:   Options Database Keys:
805: . -viewer_hdf5_filename <name>

807:   Environmental variables:
808: . PETSC_VIEWER_HDF5_FILENAME

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

815: .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy()
816: @*/
817: PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
818: {
820:   PetscBool      flg;
821:   PetscViewer    viewer;
822:   char           fname[PETSC_MAX_PATH_LEN];
823:   MPI_Comm       ncomm;

826:   PetscCommDuplicate(comm,&ncomm,NULL);if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
827:   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
828:     MPI_Keyval_create(MPI_NULL_COPY_FN,MPI_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0);
829:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
830:   }
831:   MPI_Attr_get(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg);
832:   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
833:   if (!flg) { /* PetscViewer not yet created */
834:     PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg);
835:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
836:     if (!flg) {
837:       PetscStrcpy(fname,"output.h5");
838:       if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
839:     }
840:     PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer);
841:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
842:     PetscObjectRegisterDestroy((PetscObject)viewer);
843:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
844:     MPI_Attr_put(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer);
845:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
846:   }
847:   PetscCommDestroy(&ncomm);
848:   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");return(0);}
849:   PetscFunctionReturn(viewer);
850: }