16b9fdc0aSVaclav Hapla #include <petsc/private/viewerhdf5impl.h> /*I "petscviewerhdf5.h" I*/ 25c6c1daeSBarry Smith 3bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer, const char[], PetscBool, PetscBool*, H5O_type_t*); 406db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer, const char[], const char[], PetscBool*); 506db490cSVaclav Hapla 64302210dSVaclav Hapla static PetscErrorCode PetscViewerHDF5GetAbsolutePath_Internal(PetscViewer viewer, const char path[], char *abspath[]) 76c132bc1SVaclav Hapla { 84302210dSVaclav Hapla PetscBool relative = PETSC_FALSE; 96c132bc1SVaclav Hapla const char *group; 106c132bc1SVaclav Hapla char buf[PETSC_MAX_PATH_LEN] = ""; 116c132bc1SVaclav Hapla PetscErrorCode ierr; 126c132bc1SVaclav Hapla 136c132bc1SVaclav Hapla PetscFunctionBegin; 146c132bc1SVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr); 154302210dSVaclav Hapla ierr = PetscViewerHDF5PathIsRelative(path, PETSC_TRUE, &relative);CHKERRQ(ierr); 164302210dSVaclav Hapla if (relative) { 174302210dSVaclav Hapla ierr = PetscStrcpy(buf, group);CHKERRQ(ierr); 186c132bc1SVaclav Hapla ierr = PetscStrcat(buf, "/");CHKERRQ(ierr); 194302210dSVaclav Hapla ierr = PetscStrcat(buf, path);CHKERRQ(ierr); 204302210dSVaclav Hapla ierr = PetscStrallocpy(buf, abspath);CHKERRQ(ierr); 214302210dSVaclav Hapla } else { 224302210dSVaclav Hapla ierr = PetscStrallocpy(path, abspath);CHKERRQ(ierr); 234302210dSVaclav Hapla } 246c132bc1SVaclav Hapla PetscFunctionReturn(0); 256c132bc1SVaclav Hapla } 266c132bc1SVaclav Hapla 27577e0e04SVaclav Hapla static PetscErrorCode PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer, PetscObject obj) 28577e0e04SVaclav Hapla { 29577e0e04SVaclav Hapla PetscBool has; 30577e0e04SVaclav Hapla PetscErrorCode ierr; 31577e0e04SVaclav Hapla 32577e0e04SVaclav Hapla PetscFunctionBegin; 33577e0e04SVaclav Hapla ierr = PetscViewerHDF5HasObject(viewer, obj, &has);CHKERRQ(ierr); 34577e0e04SVaclav Hapla if (!has) { 3589e0ef10SVaclav Hapla const char *group; 36577e0e04SVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr); 3789e0ef10SVaclav Hapla SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) \"%s\" not stored in group %s", obj->name, group ? group : "/"); 38577e0e04SVaclav Hapla } 39577e0e04SVaclav Hapla PetscFunctionReturn(0); 40577e0e04SVaclav Hapla } 41577e0e04SVaclav Hapla 424416b707SBarry Smith static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptionItems *PetscOptionsObject,PetscViewer v) 4382ea9b62SBarry Smith { 4482ea9b62SBarry Smith PetscErrorCode ierr; 45ee8b9732SVaclav Hapla PetscBool flg = PETSC_FALSE, set; 4682ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data; 4782ea9b62SBarry Smith 4882ea9b62SBarry Smith PetscFunctionBegin; 4982ea9b62SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");CHKERRQ(ierr); 5082ea9b62SBarry Smith ierr = PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);CHKERRQ(ierr); 519a0502c6SHåkon Strandenes ierr = PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL);CHKERRQ(ierr); 52ee8b9732SVaclav Hapla ierr = PetscOptionsBool("-viewer_hdf5_collective","Enable collective transfer mode","PetscViewerHDF5SetCollective",flg,&flg,&set);CHKERRQ(ierr); 53ee8b9732SVaclav Hapla if (set) {ierr = PetscViewerHDF5SetCollective(v,flg);CHKERRQ(ierr);} 5419a20e4cSMatthew G. Knepley flg = PETSC_FALSE; 5519a20e4cSMatthew G. Knepley ierr = PetscOptionsBool("-viewer_hdf5_default_timestepping","Set default timestepping state","PetscViewerHDF5SetDefaultTimestepping",flg,&flg,&set);CHKERRQ(ierr); 5619a20e4cSMatthew G. Knepley if (set) {ierr = PetscViewerHDF5SetDefaultTimestepping(v,flg);CHKERRQ(ierr);} 5782ea9b62SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 5882ea9b62SBarry Smith PetscFunctionReturn(0); 5982ea9b62SBarry Smith } 6082ea9b62SBarry Smith 611b793a25SVaclav Hapla static PetscErrorCode PetscViewerView_HDF5(PetscViewer v,PetscViewer viewer) 621b793a25SVaclav Hapla { 631b793a25SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data; 6403fa8834SVaclav Hapla PetscBool flg; 651b793a25SVaclav Hapla PetscErrorCode ierr; 661b793a25SVaclav Hapla 671b793a25SVaclav Hapla PetscFunctionBegin; 681b793a25SVaclav Hapla if (hdf5->filename) { 691b793a25SVaclav Hapla ierr = PetscViewerASCIIPrintf(viewer,"Filename: %s\n",hdf5->filename);CHKERRQ(ierr); 701b793a25SVaclav Hapla } 711b793a25SVaclav Hapla ierr = PetscViewerASCIIPrintf(viewer,"Vectors with blocksize 1 saved as 2D datasets: %s\n",PetscBools[hdf5->basedimension2]);CHKERRQ(ierr); 721b793a25SVaclav Hapla ierr = PetscViewerASCIIPrintf(viewer,"Enforce single precision storage: %s\n",PetscBools[hdf5->spoutput]);CHKERRQ(ierr); 7303fa8834SVaclav Hapla ierr = PetscViewerHDF5GetCollective(v,&flg);CHKERRQ(ierr); 7403fa8834SVaclav Hapla ierr = PetscViewerASCIIPrintf(viewer,"MPI-IO transfer mode: %s\n",flg ? "collective" : "independent");CHKERRQ(ierr); 7519a20e4cSMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer,"Default timestepping: %s\n",PetscBools[hdf5->defTimestepping]);CHKERRQ(ierr); 761b793a25SVaclav Hapla PetscFunctionReturn(0); 771b793a25SVaclav Hapla } 781b793a25SVaclav Hapla 795c6c1daeSBarry Smith static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer) 805c6c1daeSBarry Smith { 815c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data; 825c6c1daeSBarry Smith PetscErrorCode ierr; 835c6c1daeSBarry Smith 845c6c1daeSBarry Smith PetscFunctionBegin; 855c6c1daeSBarry Smith ierr = PetscFree(hdf5->filename);CHKERRQ(ierr); 86729ab6d9SBarry Smith if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 875c6c1daeSBarry Smith PetscFunctionReturn(0); 885c6c1daeSBarry Smith } 895c6c1daeSBarry Smith 906226335aSVaclav Hapla static PetscErrorCode PetscViewerFlush_HDF5(PetscViewer viewer) 916226335aSVaclav Hapla { 926226335aSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data; 936226335aSVaclav Hapla 946226335aSVaclav Hapla PetscFunctionBegin; 956226335aSVaclav Hapla if (hdf5->file_id) PetscStackCallHDF5(H5Fflush,(hdf5->file_id, H5F_SCOPE_LOCAL)); 966226335aSVaclav Hapla PetscFunctionReturn(0); 976226335aSVaclav Hapla } 986226335aSVaclav Hapla 999b1403beSVaclav Hapla static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer) 1005c6c1daeSBarry Smith { 1015c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 1025c6c1daeSBarry Smith PetscErrorCode ierr; 1035c6c1daeSBarry Smith 1045c6c1daeSBarry Smith PetscFunctionBegin; 1059c5072fbSVaclav Hapla PetscStackCallHDF5(H5Pclose,(hdf5->dxpl_id)); 1065c6c1daeSBarry Smith ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr); 1075c6c1daeSBarry Smith while (hdf5->groups) { 1087d964842SVaclav Hapla PetscViewerHDF5GroupList *tmp = hdf5->groups->next; 1095c6c1daeSBarry Smith 1105c6c1daeSBarry Smith ierr = PetscFree(hdf5->groups->name);CHKERRQ(ierr); 1115c6c1daeSBarry Smith ierr = PetscFree(hdf5->groups);CHKERRQ(ierr); 1125c6c1daeSBarry Smith hdf5->groups = tmp; 1135c6c1daeSBarry Smith } 1145c6c1daeSBarry Smith ierr = PetscFree(hdf5);CHKERRQ(ierr); 1150b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr); 116d1232d7fSToby Isaac ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr); 1170b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr); 118058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetBaseDimension2_C",NULL);CHKERRQ(ierr); 119058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetSPOutput_C",NULL);CHKERRQ(ierr); 12019a20e4cSMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5GetDefaultTimestepping_C",NULL);CHKERRQ(ierr); 12119a20e4cSMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetDefaultTimestepping_C",NULL);CHKERRQ(ierr); 1225c6c1daeSBarry Smith PetscFunctionReturn(0); 1235c6c1daeSBarry Smith } 1245c6c1daeSBarry Smith 1259b1403beSVaclav Hapla static PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type) 1265c6c1daeSBarry Smith { 1275c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 1285c6c1daeSBarry Smith 1295c6c1daeSBarry Smith PetscFunctionBegin; 1305c6c1daeSBarry Smith hdf5->btype = type; 1315c6c1daeSBarry Smith PetscFunctionReturn(0); 1325c6c1daeSBarry Smith } 1335c6c1daeSBarry Smith 1340b62783dSVaclav Hapla static PetscErrorCode PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type) 1350b62783dSVaclav Hapla { 1360b62783dSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 1370b62783dSVaclav Hapla 1380b62783dSVaclav Hapla PetscFunctionBegin; 1390b62783dSVaclav Hapla *type = hdf5->btype; 1400b62783dSVaclav Hapla PetscFunctionReturn(0); 1410b62783dSVaclav Hapla } 1420b62783dSVaclav Hapla 1439b1403beSVaclav Hapla static PetscErrorCode PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg) 14482ea9b62SBarry Smith { 14582ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 14682ea9b62SBarry Smith 14782ea9b62SBarry Smith PetscFunctionBegin; 14882ea9b62SBarry Smith hdf5->basedimension2 = flg; 14982ea9b62SBarry Smith PetscFunctionReturn(0); 15082ea9b62SBarry Smith } 15182ea9b62SBarry Smith 152df863907SAlex Fikl /*@ 15382ea9b62SBarry Smith PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 15482ea9b62SBarry Smith dimension of 2. 15582ea9b62SBarry Smith 15682ea9b62SBarry Smith Logically Collective on PetscViewer 15782ea9b62SBarry Smith 15882ea9b62SBarry Smith Input Parameters: 15982ea9b62SBarry Smith + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 16082ea9b62SBarry Smith - flg - if PETSC_TRUE the vector will always have at least a dimension of 2 even if that first dimension is of size 1 16182ea9b62SBarry Smith 16282ea9b62SBarry Smith Options Database: 16382ea9b62SBarry Smith . -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 16482ea9b62SBarry Smith 16595452b02SPatrick Sanan Notes: 16695452b02SPatrick Sanan Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof 16782ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 16882ea9b62SBarry Smith 16982ea9b62SBarry Smith Level: intermediate 17082ea9b62SBarry Smith 17182ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 17282ea9b62SBarry Smith 17382ea9b62SBarry Smith @*/ 17482ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg) 17582ea9b62SBarry Smith { 17682ea9b62SBarry Smith PetscErrorCode ierr; 17782ea9b62SBarry Smith 17882ea9b62SBarry Smith PetscFunctionBegin; 17982ea9b62SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 18082ea9b62SBarry Smith ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 18182ea9b62SBarry Smith PetscFunctionReturn(0); 18282ea9b62SBarry Smith } 18382ea9b62SBarry Smith 184df863907SAlex Fikl /*@ 18582ea9b62SBarry Smith PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 18682ea9b62SBarry Smith dimension of 2. 18782ea9b62SBarry Smith 18882ea9b62SBarry Smith Logically Collective on PetscViewer 18982ea9b62SBarry Smith 19082ea9b62SBarry Smith Input Parameter: 19182ea9b62SBarry Smith . viewer - the PetscViewer, must be of type HDF5 19282ea9b62SBarry Smith 19382ea9b62SBarry Smith Output Parameter: 19482ea9b62SBarry Smith . flg - if PETSC_TRUE the vector will always have at least a dimension of 2 even if that first dimension is of size 1 19582ea9b62SBarry Smith 19695452b02SPatrick Sanan Notes: 19795452b02SPatrick Sanan Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof 19882ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 19982ea9b62SBarry Smith 20082ea9b62SBarry Smith Level: intermediate 20182ea9b62SBarry Smith 20282ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 20382ea9b62SBarry Smith 20482ea9b62SBarry Smith @*/ 20582ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg) 20682ea9b62SBarry Smith { 20782ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 20882ea9b62SBarry Smith 20982ea9b62SBarry Smith PetscFunctionBegin; 21082ea9b62SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 21182ea9b62SBarry Smith *flg = hdf5->basedimension2; 21282ea9b62SBarry Smith PetscFunctionReturn(0); 21382ea9b62SBarry Smith } 21482ea9b62SBarry Smith 2159b1403beSVaclav Hapla static PetscErrorCode PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg) 2169a0502c6SHåkon Strandenes { 2179a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 2189a0502c6SHåkon Strandenes 2199a0502c6SHåkon Strandenes PetscFunctionBegin; 2209a0502c6SHåkon Strandenes hdf5->spoutput = flg; 2219a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2229a0502c6SHåkon Strandenes } 2239a0502c6SHåkon Strandenes 224df863907SAlex Fikl /*@ 2259a0502c6SHåkon Strandenes PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is 2269a0502c6SHåkon Strandenes compiled with double precision PetscReal. 2279a0502c6SHåkon Strandenes 2289a0502c6SHåkon Strandenes Logically Collective on PetscViewer 2299a0502c6SHåkon Strandenes 2309a0502c6SHåkon Strandenes Input Parameters: 2319a0502c6SHåkon Strandenes + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 2329a0502c6SHåkon Strandenes - flg - if PETSC_TRUE the data will be written to disk with single precision 2339a0502c6SHåkon Strandenes 2349a0502c6SHåkon Strandenes Options Database: 2359a0502c6SHåkon Strandenes . -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision 2369a0502c6SHåkon Strandenes 23795452b02SPatrick Sanan Notes: 23895452b02SPatrick Sanan Setting this option does not make any difference if PETSc is compiled with single precision 2399a0502c6SHåkon Strandenes in the first place. It does not affect reading datasets (HDF5 handle this internally). 2409a0502c6SHåkon Strandenes 2419a0502c6SHåkon Strandenes Level: intermediate 2429a0502c6SHåkon Strandenes 2439a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 2449a0502c6SHåkon Strandenes PetscReal 2459a0502c6SHåkon Strandenes 2469a0502c6SHåkon Strandenes @*/ 2479a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg) 2489a0502c6SHåkon Strandenes { 2499a0502c6SHåkon Strandenes PetscErrorCode ierr; 2509a0502c6SHåkon Strandenes 2519a0502c6SHåkon Strandenes PetscFunctionBegin; 2529a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 2539a0502c6SHåkon Strandenes ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 2549a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2559a0502c6SHåkon Strandenes } 2569a0502c6SHåkon Strandenes 257df863907SAlex Fikl /*@ 2589a0502c6SHåkon Strandenes PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is 2599a0502c6SHåkon Strandenes compiled with double precision PetscReal. 2609a0502c6SHåkon Strandenes 2619a0502c6SHåkon Strandenes Logically Collective on PetscViewer 2629a0502c6SHåkon Strandenes 2639a0502c6SHåkon Strandenes Input Parameter: 2649a0502c6SHåkon Strandenes . viewer - the PetscViewer, must be of type HDF5 2659a0502c6SHåkon Strandenes 2669a0502c6SHåkon Strandenes Output Parameter: 2679a0502c6SHåkon Strandenes . flg - if PETSC_TRUE the data will be written to disk with single precision 2689a0502c6SHåkon Strandenes 26995452b02SPatrick Sanan Notes: 27095452b02SPatrick Sanan Setting this option does not make any difference if PETSc is compiled with single precision 2719a0502c6SHåkon Strandenes in the first place. It does not affect reading datasets (HDF5 handle this internally). 2729a0502c6SHåkon Strandenes 2739a0502c6SHåkon Strandenes Level: intermediate 2749a0502c6SHåkon Strandenes 2759a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 2769a0502c6SHåkon Strandenes PetscReal 2779a0502c6SHåkon Strandenes 2789a0502c6SHåkon Strandenes @*/ 2799a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg) 2809a0502c6SHåkon Strandenes { 2819a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 2829a0502c6SHåkon Strandenes 2839a0502c6SHåkon Strandenes PetscFunctionBegin; 2849a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 2859a0502c6SHåkon Strandenes *flg = hdf5->spoutput; 2869a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2879a0502c6SHåkon Strandenes } 2889a0502c6SHåkon Strandenes 289ee8b9732SVaclav Hapla static PetscErrorCode PetscViewerHDF5SetCollective_HDF5(PetscViewer viewer, PetscBool flg) 290ee8b9732SVaclav Hapla { 291ee8b9732SVaclav Hapla PetscFunctionBegin; 292ee8b9732SVaclav Hapla /* H5FD_MPIO_COLLECTIVE is wrong in hdf5 1.10.2, and is the same as H5FD_MPIO_INDEPENDENT in earlier versions 293ee8b9732SVaclav Hapla - see e.g. https://gitlab.cosma.dur.ac.uk/swift/swiftsim/issues/431 */ 294c796909eSBarry Smith #if H5_VERSION_GE(1,10,3) && defined(H5_HAVE_PARALLEL) 295ee8b9732SVaclav Hapla { 296ee8b9732SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 297ee8b9732SVaclav Hapla PetscStackCallHDF5(H5Pset_dxpl_mpio,(hdf5->dxpl_id, flg ? H5FD_MPIO_COLLECTIVE : H5FD_MPIO_INDEPENDENT)); 298ee8b9732SVaclav Hapla } 299ee8b9732SVaclav Hapla #else 300ee8b9732SVaclav Hapla if (flg) { 301ee8b9732SVaclav Hapla PetscErrorCode ierr; 302c796909eSBarry Smith ierr = 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");CHKERRQ(ierr); 303ee8b9732SVaclav Hapla } 304ee8b9732SVaclav Hapla #endif 305ee8b9732SVaclav Hapla PetscFunctionReturn(0); 306ee8b9732SVaclav Hapla } 307ee8b9732SVaclav Hapla 308ee8b9732SVaclav Hapla /*@ 309ee8b9732SVaclav Hapla PetscViewerHDF5SetCollective - Use collective MPI-IO transfer mode for HDF5 reads and writes. 310ee8b9732SVaclav Hapla 311ee8b9732SVaclav Hapla Logically Collective; flg must contain common value 312ee8b9732SVaclav Hapla 313ee8b9732SVaclav Hapla Input Parameters: 314ee8b9732SVaclav Hapla + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 315ee8b9732SVaclav Hapla - flg - PETSC_TRUE for collective mode; PETSC_FALSE for independent mode (default) 316ee8b9732SVaclav Hapla 317ee8b9732SVaclav Hapla Options Database: 318ee8b9732SVaclav Hapla . -viewer_hdf5_collective - turns on (true) or off (false) collective transfers 319ee8b9732SVaclav Hapla 320ee8b9732SVaclav Hapla Notes: 321ee8b9732SVaclav Hapla Collective mode gives the MPI-IO layer underneath HDF5 a chance to do some additional collective optimizations and hence can perform better. 32253effed7SBarry Smith 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. 323ee8b9732SVaclav Hapla 324ee8b9732SVaclav Hapla Developer notes: 325ee8b9732SVaclav Hapla In the HDF5 layer, PETSC_TRUE / PETSC_FALSE means H5Pset_dxpl_mpio() is called with H5FD_MPIO_COLLECTIVE / H5FD_MPIO_INDEPENDENT, respectively. 326ee8b9732SVaclav Hapla This in turn means use of MPI_File_{read,write}_all / MPI_File_{read,write} in the MPI-IO layer, respectively. 327ee8b9732SVaclav Hapla See HDF5 documentation and MPI-IO documentation for details. 328ee8b9732SVaclav Hapla 329ee8b9732SVaclav Hapla Level: intermediate 330ee8b9732SVaclav Hapla 331ee8b9732SVaclav Hapla .seealso: PetscViewerHDF5GetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open() 332ee8b9732SVaclav Hapla 333ee8b9732SVaclav Hapla @*/ 334ee8b9732SVaclav Hapla PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer viewer,PetscBool flg) 335ee8b9732SVaclav Hapla { 336ee8b9732SVaclav Hapla PetscErrorCode ierr; 337ee8b9732SVaclav Hapla 338ee8b9732SVaclav Hapla PetscFunctionBegin; 339ee8b9732SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 340ee8b9732SVaclav Hapla PetscValidLogicalCollectiveBool(viewer,flg,2); 341ee8b9732SVaclav Hapla ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetCollective_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 342ee8b9732SVaclav Hapla PetscFunctionReturn(0); 343ee8b9732SVaclav Hapla } 344ee8b9732SVaclav Hapla 345ee8b9732SVaclav Hapla static PetscErrorCode PetscViewerHDF5GetCollective_HDF5(PetscViewer viewer, PetscBool *flg) 346ee8b9732SVaclav Hapla { 347c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL) 348ee8b9732SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 349ee8b9732SVaclav Hapla H5FD_mpio_xfer_t mode; 350c796909eSBarry Smith #endif 351ee8b9732SVaclav Hapla 352ee8b9732SVaclav Hapla PetscFunctionBegin; 353c796909eSBarry Smith #if !defined(H5_HAVE_PARALLEL) 354c796909eSBarry Smith *flg = PETSC_FALSE; 355c796909eSBarry Smith #else 356ee8b9732SVaclav Hapla PetscStackCallHDF5(H5Pget_dxpl_mpio,(hdf5->dxpl_id, &mode)); 357ee8b9732SVaclav Hapla *flg = (mode == H5FD_MPIO_COLLECTIVE) ? PETSC_TRUE : PETSC_FALSE; 358c796909eSBarry Smith #endif 359ee8b9732SVaclav Hapla PetscFunctionReturn(0); 360ee8b9732SVaclav Hapla } 361ee8b9732SVaclav Hapla 362ee8b9732SVaclav Hapla /*@ 363ee8b9732SVaclav Hapla PetscViewerHDF5GetCollective - Return flag whether collective MPI-IO transfer mode is used for HDF5 reads and writes. 364ee8b9732SVaclav Hapla 365ee8b9732SVaclav Hapla Not Collective 366ee8b9732SVaclav Hapla 367ee8b9732SVaclav Hapla Input Parameters: 368ee8b9732SVaclav Hapla . viewer - the HDF5 PetscViewer 369ee8b9732SVaclav Hapla 370ee8b9732SVaclav Hapla Output Parameters: 371ee8b9732SVaclav Hapla . flg - the flag 372ee8b9732SVaclav Hapla 373ee8b9732SVaclav Hapla Level: intermediate 374ee8b9732SVaclav Hapla 375ee8b9732SVaclav Hapla Notes: 376c796909eSBarry Smith 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. 377ee8b9732SVaclav Hapla For more details, see PetscViewerHDF5SetCollective(). 378ee8b9732SVaclav Hapla 379ee8b9732SVaclav Hapla .seealso: PetscViewerHDF5SetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open() 380ee8b9732SVaclav Hapla 381ee8b9732SVaclav Hapla @*/ 382ee8b9732SVaclav Hapla PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer viewer,PetscBool *flg) 383ee8b9732SVaclav Hapla { 384ee8b9732SVaclav Hapla PetscErrorCode ierr; 385ee8b9732SVaclav Hapla 386ee8b9732SVaclav Hapla PetscFunctionBegin; 387ee8b9732SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 388534a8f05SLisandro Dalcin PetscValidBoolPointer(flg,2); 389ee8b9732SVaclav Hapla 390ee8b9732SVaclav Hapla ierr = PetscUseMethod(viewer,"PetscViewerHDF5GetCollective_C",(PetscViewer,PetscBool*),(viewer,flg));CHKERRQ(ierr); 391ee8b9732SVaclav Hapla PetscFunctionReturn(0); 392ee8b9732SVaclav Hapla } 393ee8b9732SVaclav Hapla 3949b1403beSVaclav Hapla static PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[]) 3955c6c1daeSBarry Smith { 3965c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 3975c6c1daeSBarry Smith hid_t plist_id; 3985c6c1daeSBarry Smith PetscErrorCode ierr; 3995c6c1daeSBarry Smith 4005c6c1daeSBarry Smith PetscFunctionBegin; 401f683cc58SToby Isaac if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 402f683cc58SToby Isaac if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);} 4035c6c1daeSBarry Smith ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr); 4045c6c1daeSBarry Smith /* Set up file access property list with parallel I/O access */ 405729ab6d9SBarry Smith PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS)); 406c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL) 407d28bfa55SVaclav Hapla PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), MPI_INFO_NULL)); 408c796909eSBarry Smith #endif 4095c6c1daeSBarry Smith /* Create or open the file collectively */ 4105c6c1daeSBarry Smith switch (hdf5->btype) { 4115c6c1daeSBarry Smith case FILE_MODE_READ: 412*8a2871f6SBarry Smith if (PetscDefined(USE_DEBUG)) { 413*8a2871f6SBarry Smith PetscMPIInt rank; 414*8a2871f6SBarry Smith PetscBool flg; 415*8a2871f6SBarry Smith 416*8a2871f6SBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)viewer),&rank);CHKERRMPI(ierr); 417*8a2871f6SBarry Smith if (rank == 0) { 418*8a2871f6SBarry Smith ierr = PetscTestFile(hdf5->filename, 'r', &flg);CHKERRQ(ierr); 419*8a2871f6SBarry Smith if (!flg) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"File %s requested for reading does not exist",hdf5->filename); 420*8a2871f6SBarry Smith } 421*8a2871f6SBarry Smith ierr = MPI_Barrier(PetscObjectComm((PetscObject)viewer));CHKERRMPI(ierr); 422*8a2871f6SBarry Smith } 423729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id)); 4245c6c1daeSBarry Smith break; 4255c6c1daeSBarry Smith case FILE_MODE_APPEND: 4267e4fd573SVaclav Hapla case FILE_MODE_UPDATE: 4277e4fd573SVaclav Hapla { 4287e4fd573SVaclav Hapla PetscBool flg; 4297e4fd573SVaclav Hapla ierr = PetscTestFile(hdf5->filename, 'r', &flg);CHKERRQ(ierr); 4307e4fd573SVaclav Hapla if (flg) PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id)); 4317e4fd573SVaclav Hapla else PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_EXCL, H5P_DEFAULT, plist_id)); 4325c6c1daeSBarry Smith break; 4337e4fd573SVaclav Hapla } 4345c6c1daeSBarry Smith case FILE_MODE_WRITE: 435729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id)); 4365c6c1daeSBarry Smith break; 4377e4fd573SVaclav Hapla case FILE_MODE_UNDEFINED: 4387e4fd573SVaclav Hapla SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 4395c6c1daeSBarry Smith default: 4407e4fd573SVaclav Hapla SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP, "Unsupported file mode %s",PetscFileModes[hdf5->btype]); 4415c6c1daeSBarry Smith } 4425c6c1daeSBarry Smith if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name); 443729ab6d9SBarry Smith PetscStackCallHDF5(H5Pclose,(plist_id)); 4445c6c1daeSBarry Smith PetscFunctionReturn(0); 4455c6c1daeSBarry Smith } 4465c6c1daeSBarry Smith 447d1232d7fSToby Isaac static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name) 448d1232d7fSToby Isaac { 449d1232d7fSToby Isaac PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data; 450d1232d7fSToby Isaac 451d1232d7fSToby Isaac PetscFunctionBegin; 452d1232d7fSToby Isaac *name = vhdf5->filename; 453d1232d7fSToby Isaac PetscFunctionReturn(0); 454d1232d7fSToby Isaac } 455d1232d7fSToby Isaac 456b723ab35SVaclav Hapla static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer) 457b723ab35SVaclav Hapla { 4585dc64a97SVaclav Hapla /* 459b723ab35SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 460b723ab35SVaclav Hapla PetscErrorCode ierr; 4615dc64a97SVaclav Hapla */ 462b723ab35SVaclav Hapla 463b723ab35SVaclav Hapla PetscFunctionBegin; 464b723ab35SVaclav Hapla PetscFunctionReturn(0); 465b723ab35SVaclav Hapla } 466b723ab35SVaclav Hapla 46719a20e4cSMatthew G. Knepley static PetscErrorCode PetscViewerHDF5SetDefaultTimestepping_HDF5(PetscViewer viewer, PetscBool flg) 46819a20e4cSMatthew G. Knepley { 46919a20e4cSMatthew G. Knepley PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 47019a20e4cSMatthew G. Knepley 47119a20e4cSMatthew G. Knepley PetscFunctionBegin; 47219a20e4cSMatthew G. Knepley hdf5->defTimestepping = flg; 47319a20e4cSMatthew G. Knepley PetscFunctionReturn(0); 47419a20e4cSMatthew G. Knepley } 47519a20e4cSMatthew G. Knepley 47619a20e4cSMatthew G. Knepley static PetscErrorCode PetscViewerHDF5GetDefaultTimestepping_HDF5(PetscViewer viewer, PetscBool *flg) 47719a20e4cSMatthew G. Knepley { 47819a20e4cSMatthew G. Knepley PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 47919a20e4cSMatthew G. Knepley 48019a20e4cSMatthew G. Knepley PetscFunctionBegin; 48119a20e4cSMatthew G. Knepley *flg = hdf5->defTimestepping; 48219a20e4cSMatthew G. Knepley PetscFunctionReturn(0); 48319a20e4cSMatthew G. Knepley } 48419a20e4cSMatthew G. Knepley 48519a20e4cSMatthew G. Knepley /*@ 48619a20e4cSMatthew G. Knepley PetscViewerHDF5SetDefaultTimestepping - Set the flag for default timestepping 48719a20e4cSMatthew G. Knepley 48819a20e4cSMatthew G. Knepley Logically Collective on PetscViewer 48919a20e4cSMatthew G. Knepley 49019a20e4cSMatthew G. Knepley Input Parameters: 49119a20e4cSMatthew G. Knepley + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 49219a20e4cSMatthew G. Knepley - flg - if PETSC_TRUE we will assume that timestepping is on 49319a20e4cSMatthew G. Knepley 49419a20e4cSMatthew G. Knepley Options Database: 49519a20e4cSMatthew G. Knepley . -viewer_hdf5_default_timestepping - turns on (true) or off (false) default timestepping 49619a20e4cSMatthew G. Knepley 49719a20e4cSMatthew G. Knepley Notes: 49819a20e4cSMatthew G. Knepley If the timestepping attribute is not found for an object, then the default timestepping is used 49919a20e4cSMatthew G. Knepley 50019a20e4cSMatthew G. Knepley Level: intermediate 50119a20e4cSMatthew G. Knepley 50219a20e4cSMatthew G. Knepley .seealso: PetscViewerHDF5GetDefaultTimestepping(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5GetTimestep() 50319a20e4cSMatthew G. Knepley @*/ 50419a20e4cSMatthew G. Knepley PetscErrorCode PetscViewerHDF5SetDefaultTimestepping(PetscViewer viewer, PetscBool flg) 50519a20e4cSMatthew G. Knepley { 50619a20e4cSMatthew G. Knepley PetscErrorCode ierr; 50719a20e4cSMatthew G. Knepley 50819a20e4cSMatthew G. Knepley PetscFunctionBegin; 50919a20e4cSMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 51019a20e4cSMatthew G. Knepley ierr = PetscTryMethod(viewer, "PetscViewerHDF5SetDefaultTimestepping_C", (PetscViewer,PetscBool), (viewer,flg));CHKERRQ(ierr); 51119a20e4cSMatthew G. Knepley PetscFunctionReturn(0); 51219a20e4cSMatthew G. Knepley } 51319a20e4cSMatthew G. Knepley 51419a20e4cSMatthew G. Knepley /*@ 51519a20e4cSMatthew G. Knepley PetscViewerHDF5GetDefaultTimestepping - Get the flag for default timestepping 51619a20e4cSMatthew G. Knepley 51719a20e4cSMatthew G. Knepley Not collective 51819a20e4cSMatthew G. Knepley 51919a20e4cSMatthew G. Knepley Input Parameter: 52019a20e4cSMatthew G. Knepley . viewer - the PetscViewer 52119a20e4cSMatthew G. Knepley 52219a20e4cSMatthew G. Knepley Output Parameter: 52319a20e4cSMatthew G. Knepley . flg - if PETSC_TRUE we will assume that timestepping is on 52419a20e4cSMatthew G. Knepley 52519a20e4cSMatthew G. Knepley Notes: 52619a20e4cSMatthew G. Knepley If the timestepping attribute is not found for an object, then the default timestepping is used 52719a20e4cSMatthew G. Knepley 52819a20e4cSMatthew G. Knepley Level: intermediate 52919a20e4cSMatthew G. Knepley 53019a20e4cSMatthew G. Knepley .seealso: PetscViewerHDF5SetDefaultTimestepping(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5GetTimestep() 53119a20e4cSMatthew G. Knepley @*/ 53219a20e4cSMatthew G. Knepley PetscErrorCode PetscViewerHDF5GetDefaultTimestepping(PetscViewer viewer, PetscBool *flg) 53319a20e4cSMatthew G. Knepley { 53419a20e4cSMatthew G. Knepley PetscErrorCode ierr; 53519a20e4cSMatthew G. Knepley 53619a20e4cSMatthew G. Knepley PetscFunctionBegin; 53719a20e4cSMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 53819a20e4cSMatthew G. Knepley ierr = PetscUseMethod(viewer, "PetscViewerHDF5GetDefaultTimestepping_C", (PetscViewer,PetscBool*), (viewer,flg));CHKERRQ(ierr); 53919a20e4cSMatthew G. Knepley PetscFunctionReturn(0); 54019a20e4cSMatthew G. Knepley } 54119a20e4cSMatthew G. Knepley 5428556b5ebSBarry Smith /*MC 5438556b5ebSBarry Smith PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file 5448556b5ebSBarry Smith 5458556b5ebSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET, 5468556b5ebSBarry Smith PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING, 5478556b5ebSBarry Smith PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB, 5488556b5ebSBarry Smith PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType() 5498556b5ebSBarry Smith 5501b266c99SBarry Smith Level: beginner 5518556b5ebSBarry Smith M*/ 552d1232d7fSToby Isaac 5538cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v) 5545c6c1daeSBarry Smith { 5555c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5; 5565c6c1daeSBarry Smith PetscErrorCode ierr; 5575c6c1daeSBarry Smith 5585c6c1daeSBarry Smith PetscFunctionBegin; 55999335c34SVaclav Hapla #if !defined(H5_HAVE_PARALLEL) 56099335c34SVaclav Hapla { 56199335c34SVaclav Hapla PetscMPIInt size; 56299335c34SVaclav Hapla ierr = MPI_Comm_size(PetscObjectComm((PetscObject)v), &size);CHKERRMPI(ierr); 56399335c34SVaclav Hapla 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)"); 56499335c34SVaclav Hapla } 56599335c34SVaclav Hapla #endif 56699335c34SVaclav Hapla 567b00a9115SJed Brown ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr); 5685c6c1daeSBarry Smith 5695c6c1daeSBarry Smith v->data = (void*) hdf5; 5705c6c1daeSBarry Smith v->ops->destroy = PetscViewerDestroy_HDF5; 57182ea9b62SBarry Smith v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5; 572b723ab35SVaclav Hapla v->ops->setup = PetscViewerSetUp_HDF5; 5731b793a25SVaclav Hapla v->ops->view = PetscViewerView_HDF5; 5746226335aSVaclav Hapla v->ops->flush = PetscViewerFlush_HDF5; 5757e4fd573SVaclav Hapla hdf5->btype = FILE_MODE_UNDEFINED; 576908793a3SLisandro Dalcin hdf5->filename = NULL; 5775c6c1daeSBarry Smith hdf5->timestep = -1; 5780298fd71SBarry Smith hdf5->groups = NULL; 5795c6c1daeSBarry Smith 5809c5072fbSVaclav Hapla PetscStackCallHDF5Return(hdf5->dxpl_id,H5Pcreate,(H5P_DATASET_XFER)); 5819c5072fbSVaclav Hapla 5820b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr); 583d1232d7fSToby Isaac ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr); 5840b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr); 5850b62783dSVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_HDF5);CHKERRQ(ierr); 58682ea9b62SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr); 5879a0502c6SHåkon Strandenes ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr); 588ee8b9732SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetCollective_C",PetscViewerHDF5SetCollective_HDF5);CHKERRQ(ierr); 589ee8b9732SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetCollective_C",PetscViewerHDF5GetCollective_HDF5);CHKERRQ(ierr); 59019a20e4cSMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetDefaultTimestepping_C",PetscViewerHDF5GetDefaultTimestepping_HDF5);CHKERRQ(ierr); 59119a20e4cSMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetDefaultTimestepping_C",PetscViewerHDF5SetDefaultTimestepping_HDF5);CHKERRQ(ierr); 5925c6c1daeSBarry Smith PetscFunctionReturn(0); 5935c6c1daeSBarry Smith } 5945c6c1daeSBarry Smith 5955c6c1daeSBarry Smith /*@C 5965c6c1daeSBarry Smith PetscViewerHDF5Open - Opens a file for HDF5 input/output. 5975c6c1daeSBarry Smith 598d083f849SBarry Smith Collective 5995c6c1daeSBarry Smith 6005c6c1daeSBarry Smith Input Parameters: 6015c6c1daeSBarry Smith + comm - MPI communicator 6025c6c1daeSBarry Smith . name - name of file 6035c6c1daeSBarry Smith - type - type of file 6045c6c1daeSBarry Smith 6055c6c1daeSBarry Smith Output Parameter: 6065c6c1daeSBarry Smith . hdf5v - PetscViewer for HDF5 input/output to use with the specified file 6075c6c1daeSBarry Smith 60882ea9b62SBarry Smith Options Database: 609a2b725a8SWilliam Gropp + -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 610a2b725a8SWilliam Gropp - -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal 61182ea9b62SBarry Smith 6125c6c1daeSBarry Smith Level: beginner 6135c6c1daeSBarry Smith 6147e4fd573SVaclav Hapla Notes: 6157e4fd573SVaclav Hapla Reading is always available, regardless of the mode. Available modes are 6167e4fd573SVaclav Hapla + FILE_MODE_READ - open existing HDF5 file for read only access, fail if file does not exist [H5Fopen() with H5F_ACC_RDONLY] 6177e4fd573SVaclav Hapla . FILE_MODE_WRITE - if file exists, fully overwrite it, else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_TRUNC] 6187e4fd573SVaclav Hapla . FILE_MODE_APPEND - if file exists, keep existing contents [H5Fopen() with H5F_ACC_RDWR], else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_EXCL] 6197e4fd573SVaclav Hapla - FILE_MODE_UPDATE - same as FILE_MODE_APPEND 6207e4fd573SVaclav Hapla 6217e4fd573SVaclav Hapla 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. 6227e4fd573SVaclav Hapla 6235c6c1daeSBarry Smith This PetscViewer should be destroyed with PetscViewerDestroy(). 6245c6c1daeSBarry Smith 6256a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(), 6269a0502c6SHåkon Strandenes PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(), 627a56f64adSBarry Smith MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName() 6285c6c1daeSBarry Smith @*/ 6295c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v) 6305c6c1daeSBarry Smith { 6315c6c1daeSBarry Smith PetscErrorCode ierr; 6325c6c1daeSBarry Smith 6335c6c1daeSBarry Smith PetscFunctionBegin; 6345c6c1daeSBarry Smith ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr); 6355c6c1daeSBarry Smith ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr); 6365c6c1daeSBarry Smith ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr); 6375c6c1daeSBarry Smith ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr); 638b663d2d2SVaclav Hapla ierr = PetscViewerSetFromOptions(*hdf5v);CHKERRQ(ierr); 6395c6c1daeSBarry Smith PetscFunctionReturn(0); 6405c6c1daeSBarry Smith } 6415c6c1daeSBarry Smith 6425c6c1daeSBarry Smith /*@C 6435c6c1daeSBarry Smith PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls 6445c6c1daeSBarry Smith 6455c6c1daeSBarry Smith Not collective 6465c6c1daeSBarry Smith 6475c6c1daeSBarry Smith Input Parameter: 6485c6c1daeSBarry Smith . viewer - the PetscViewer 6495c6c1daeSBarry Smith 6505c6c1daeSBarry Smith Output Parameter: 6515c6c1daeSBarry Smith . file_id - The file id 6525c6c1daeSBarry Smith 6535c6c1daeSBarry Smith Level: intermediate 6545c6c1daeSBarry Smith 6555c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open() 6565c6c1daeSBarry Smith @*/ 6575c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id) 6585c6c1daeSBarry Smith { 6595c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 6605c6c1daeSBarry Smith 6615c6c1daeSBarry Smith PetscFunctionBegin; 6625c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 6635c6c1daeSBarry Smith if (file_id) *file_id = hdf5->file_id; 6645c6c1daeSBarry Smith PetscFunctionReturn(0); 6655c6c1daeSBarry Smith } 6665c6c1daeSBarry Smith 6675c6c1daeSBarry Smith /*@C 6685c6c1daeSBarry Smith PetscViewerHDF5PushGroup - Set the current HDF5 group for output 6695c6c1daeSBarry Smith 6705c6c1daeSBarry Smith Not collective 6715c6c1daeSBarry Smith 6725c6c1daeSBarry Smith Input Parameters: 6735c6c1daeSBarry Smith + viewer - the PetscViewer 6745c6c1daeSBarry Smith - name - The group name 6755c6c1daeSBarry Smith 6765c6c1daeSBarry Smith Level: intermediate 6775c6c1daeSBarry Smith 6782e361344SVaclav Hapla Notes: 6792e361344SVaclav Hapla This is designed to mnemonically resemble the Unix cd command. 6802e361344SVaclav Hapla + 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. 6812e361344SVaclav Hapla . NULL, empty string, or any sequence of all slashes (e.g. "///") is interpreted as the root group "/". 6822e361344SVaclav Hapla - "." means the current group is pushed again. 6832e361344SVaclav Hapla 6842e361344SVaclav Hapla Example: 6852e361344SVaclav Hapla Suppose the current group is "/a". 6862e361344SVaclav Hapla + If name is NULL, empty string, or a sequence of all slashes (e.g. "///"), then the new group will be "/". 6872e361344SVaclav Hapla . If name is ".", then the new group will be "/a". 6882e361344SVaclav Hapla . If name is "b", then the new group will be "/a/b". 6892e361344SVaclav Hapla - If name is "/b", then the new group will be "/b". 6902e361344SVaclav Hapla 6912e361344SVaclav Hapla Developer Notes: 6922e361344SVaclav Hapla The root group "/" is internally stored as NULL. 693820fc2d1SVaclav Hapla 694874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 6955c6c1daeSBarry Smith @*/ 696be7aa430SVaclav Hapla PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char name[]) 6975c6c1daeSBarry Smith { 6985c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 6997d964842SVaclav Hapla PetscViewerHDF5GroupList *groupNode; 7002e361344SVaclav Hapla size_t i,len; 7012e361344SVaclav Hapla char buf[PETSC_MAX_PATH_LEN]; 7022e361344SVaclav Hapla const char *gname; 7035c6c1daeSBarry Smith PetscErrorCode ierr; 7045c6c1daeSBarry Smith 7055c6c1daeSBarry Smith PetscFunctionBegin; 7065c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 707820fc2d1SVaclav Hapla if (name) PetscValidCharPointer(name,2); 708820fc2d1SVaclav Hapla ierr = PetscStrlen(name, &len);CHKERRQ(ierr); 7092e361344SVaclav Hapla gname = NULL; 7102e361344SVaclav Hapla if (len) { 7112e361344SVaclav Hapla if (len == 1 && name[0] == '.') { 7122e361344SVaclav Hapla /* use current name */ 7132e361344SVaclav Hapla gname = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : NULL; 7142e361344SVaclav Hapla } else if (name[0] == '/') { 7152e361344SVaclav Hapla /* absolute */ 7162e361344SVaclav Hapla for (i=1; i<len; i++) { 7172e361344SVaclav Hapla if (name[i] != '/') { 7182e361344SVaclav Hapla gname = name; 7192e361344SVaclav Hapla break; 7202e361344SVaclav Hapla } 7212e361344SVaclav Hapla } 7222e361344SVaclav Hapla } else { 7232e361344SVaclav Hapla /* relative */ 7242e361344SVaclav Hapla const char *parent = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : ""; 7252e361344SVaclav Hapla ierr = PetscSNPrintf(buf, sizeof(buf), "%s/%s", parent, name);CHKERRQ(ierr); 7262e361344SVaclav Hapla gname = buf; 7272e361344SVaclav Hapla } 7282e361344SVaclav Hapla } 72995dccacaSBarry Smith ierr = PetscNew(&groupNode);CHKERRQ(ierr); 7302e361344SVaclav Hapla ierr = PetscStrallocpy(gname, (char**) &groupNode->name);CHKERRQ(ierr); 7315c6c1daeSBarry Smith groupNode->next = hdf5->groups; 7325c6c1daeSBarry Smith hdf5->groups = groupNode; 7335c6c1daeSBarry Smith PetscFunctionReturn(0); 7345c6c1daeSBarry Smith } 7355c6c1daeSBarry Smith 7363ef9c667SSatish Balay /*@ 7375c6c1daeSBarry Smith PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value 7385c6c1daeSBarry Smith 7395c6c1daeSBarry Smith Not collective 7405c6c1daeSBarry Smith 7415c6c1daeSBarry Smith Input Parameter: 7425c6c1daeSBarry Smith . viewer - the PetscViewer 7435c6c1daeSBarry Smith 7445c6c1daeSBarry Smith Level: intermediate 7455c6c1daeSBarry Smith 746874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 7475c6c1daeSBarry Smith @*/ 7485c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer) 7495c6c1daeSBarry Smith { 7505c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 7517d964842SVaclav Hapla PetscViewerHDF5GroupList *groupNode; 7525c6c1daeSBarry Smith PetscErrorCode ierr; 7535c6c1daeSBarry Smith 7545c6c1daeSBarry Smith PetscFunctionBegin; 7555c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 75682f516ccSBarry Smith if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop"); 7575c6c1daeSBarry Smith groupNode = hdf5->groups; 7585c6c1daeSBarry Smith hdf5->groups = hdf5->groups->next; 7595c6c1daeSBarry Smith ierr = PetscFree(groupNode->name);CHKERRQ(ierr); 7605c6c1daeSBarry Smith ierr = PetscFree(groupNode);CHKERRQ(ierr); 7615c6c1daeSBarry Smith PetscFunctionReturn(0); 7625c6c1daeSBarry Smith } 7635c6c1daeSBarry Smith 7645c6c1daeSBarry Smith /*@C 765874270d9SVaclav Hapla PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup(). 766874270d9SVaclav Hapla If none has been assigned, returns NULL. 7675c6c1daeSBarry Smith 7685c6c1daeSBarry Smith Not collective 7695c6c1daeSBarry Smith 7705c6c1daeSBarry Smith Input Parameter: 7715c6c1daeSBarry Smith . viewer - the PetscViewer 7725c6c1daeSBarry Smith 7735c6c1daeSBarry Smith Output Parameter: 7745c6c1daeSBarry Smith . name - The group name 7755c6c1daeSBarry Smith 7765c6c1daeSBarry Smith Level: intermediate 7775c6c1daeSBarry Smith 778874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup() 7795c6c1daeSBarry Smith @*/ 780be7aa430SVaclav Hapla PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char *name[]) 7815c6c1daeSBarry Smith { 7825c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 7835c6c1daeSBarry Smith 7845c6c1daeSBarry Smith PetscFunctionBegin; 7855c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 786c959eef4SJed Brown PetscValidPointer(name,2); 787a297a907SKarl Rupp if (hdf5->groups) *name = hdf5->groups->name; 7880298fd71SBarry Smith else *name = NULL; 7895c6c1daeSBarry Smith PetscFunctionReturn(0); 7905c6c1daeSBarry Smith } 7915c6c1daeSBarry Smith 7928c557b5aSVaclav Hapla /*@ 793874270d9SVaclav Hapla PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(), 794874270d9SVaclav Hapla and return this group's ID and file ID. 795874270d9SVaclav Hapla If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID. 796874270d9SVaclav Hapla 797874270d9SVaclav Hapla Not collective 798874270d9SVaclav Hapla 799874270d9SVaclav Hapla Input Parameter: 800874270d9SVaclav Hapla . viewer - the PetscViewer 801874270d9SVaclav Hapla 802d8d19677SJose E. Roman Output Parameters: 803874270d9SVaclav Hapla + fileId - The HDF5 file ID 804874270d9SVaclav Hapla - groupId - The HDF5 group ID 805874270d9SVaclav Hapla 806e74616beSVaclav Hapla Notes: 807e74616beSVaclav Hapla If the viewer is writable, the group is created if it doesn't exist yet. 808e74616beSVaclav Hapla 809874270d9SVaclav Hapla Level: intermediate 810874270d9SVaclav Hapla 811874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 812874270d9SVaclav Hapla @*/ 81354dbf706SVaclav Hapla PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId) 81454dbf706SVaclav Hapla { 81590e03537SVaclav Hapla hid_t file_id; 81690e03537SVaclav Hapla H5O_type_t type; 81776d59af2SVaclav Hapla const char *groupName = NULL, *fileName = NULL; 81876d59af2SVaclav Hapla PetscBool writable, has; 81954dbf706SVaclav Hapla PetscErrorCode ierr; 82054dbf706SVaclav Hapla 82154dbf706SVaclav Hapla PetscFunctionBegin; 82276d59af2SVaclav Hapla ierr = PetscViewerWritable(viewer, &writable);CHKERRQ(ierr); 82354dbf706SVaclav Hapla ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr); 82476d59af2SVaclav Hapla ierr = PetscViewerFileGetName(viewer, &fileName);CHKERRQ(ierr); 82554dbf706SVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr); 82676d59af2SVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, groupName, writable, &has, &type);CHKERRQ(ierr); 82776d59af2SVaclav Hapla if (!has) { 82876d59af2SVaclav Hapla 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); 82976d59af2SVaclav Hapla else SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_LIB, "HDF5 failed to create group %s although file %s is open for writing", groupName, fileName); 83076d59af2SVaclav Hapla } 83176d59af2SVaclav Hapla 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); 83290e03537SVaclav Hapla PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT)); 83354dbf706SVaclav Hapla *fileId = file_id; 83454dbf706SVaclav Hapla PetscFunctionReturn(0); 83554dbf706SVaclav Hapla } 83654dbf706SVaclav Hapla 8373ef9c667SSatish Balay /*@ 838d7dd068bSVaclav Hapla PetscViewerHDF5PushTimestepping - Activate timestepping mode for subsequent HDF5 reading and writing. 8395c6c1daeSBarry Smith 8405c6c1daeSBarry Smith Not collective 8415c6c1daeSBarry Smith 8425c6c1daeSBarry Smith Input Parameter: 8435c6c1daeSBarry Smith . viewer - the PetscViewer 8445c6c1daeSBarry Smith 8455c6c1daeSBarry Smith Level: intermediate 8465c6c1daeSBarry Smith 847d7dd068bSVaclav Hapla Notes: 848d7dd068bSVaclav Hapla On first PetscViewerHDF5PushTimestepping(), the initial time step is set to 0. 849d7dd068bSVaclav Hapla Next timesteps can then be set using PetscViewerHDF5IncrementTimestep() or PetscViewerHDF5SetTimestep(). 850d7dd068bSVaclav Hapla Current timestep value determines which timestep is read from or written to any dataset on the next HDF5 I/O operation [e.g. VecView()]. 851d7dd068bSVaclav Hapla Use PetscViewerHDF5PopTimestepping() to deactivate timestepping mode; calling it by the end of the program is NOT mandatory. 852d7dd068bSVaclav Hapla Current timestep is remembered between PetscViewerHDF5PopTimestepping() and the next PetscViewerHDF5PushTimestepping(). 853d7dd068bSVaclav Hapla 854d7dd068bSVaclav Hapla If a dataset was stored with timestepping, it can be loaded only in the timestepping mode again. 855d7dd068bSVaclav Hapla Loading a timestepped dataset with timestepping disabled, or vice-versa results in an error. 856d7dd068bSVaclav Hapla 857d7dd068bSVaclav Hapla Developer notes: 858d7dd068bSVaclav Hapla Timestepped HDF5 dataset has an extra dimension and attribute "timestepping" set to true. 859d7dd068bSVaclav Hapla 860d7dd068bSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PopTimestepping(), PetscViewerHDF5IsTimestepping(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 861d7dd068bSVaclav Hapla @*/ 862d7dd068bSVaclav Hapla PetscErrorCode PetscViewerHDF5PushTimestepping(PetscViewer viewer) 863d7dd068bSVaclav Hapla { 864d7dd068bSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 865d7dd068bSVaclav Hapla 866d7dd068bSVaclav Hapla PetscFunctionBegin; 867d7dd068bSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 868d7dd068bSVaclav Hapla if (hdf5->timestepping) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping is already pushed"); 869d7dd068bSVaclav Hapla hdf5->timestepping = PETSC_TRUE; 870d7dd068bSVaclav Hapla if (hdf5->timestep < 0) hdf5->timestep = 0; 871d7dd068bSVaclav Hapla PetscFunctionReturn(0); 872d7dd068bSVaclav Hapla } 873d7dd068bSVaclav Hapla 874d7dd068bSVaclav Hapla /*@ 875d7dd068bSVaclav Hapla PetscViewerHDF5PopTimestepping - Deactivate timestepping mode for subsequent HDF5 reading and writing. 876d7dd068bSVaclav Hapla 877d7dd068bSVaclav Hapla Not collective 878d7dd068bSVaclav Hapla 879d7dd068bSVaclav Hapla Input Parameter: 880d7dd068bSVaclav Hapla . viewer - the PetscViewer 881d7dd068bSVaclav Hapla 882d7dd068bSVaclav Hapla Level: intermediate 883d7dd068bSVaclav Hapla 884d7dd068bSVaclav Hapla Notes: 885d7dd068bSVaclav Hapla See PetscViewerHDF5PushTimestepping() for details. 886d7dd068bSVaclav Hapla 887d7dd068bSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5IsTimestepping(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 888d7dd068bSVaclav Hapla @*/ 889d7dd068bSVaclav Hapla PetscErrorCode PetscViewerHDF5PopTimestepping(PetscViewer viewer) 890d7dd068bSVaclav Hapla { 891d7dd068bSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 892d7dd068bSVaclav Hapla 893d7dd068bSVaclav Hapla PetscFunctionBegin; 894d7dd068bSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 895d7dd068bSVaclav Hapla if (!hdf5->timestepping) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 896d7dd068bSVaclav Hapla hdf5->timestepping = PETSC_FALSE; 897d7dd068bSVaclav Hapla PetscFunctionReturn(0); 898d7dd068bSVaclav Hapla } 899d7dd068bSVaclav Hapla 900d7dd068bSVaclav Hapla /*@ 901d7dd068bSVaclav Hapla PetscViewerHDF5IsTimestepping - Ask the viewer whether it is in timestepping mode currently. 902d7dd068bSVaclav Hapla 903d7dd068bSVaclav Hapla Not collective 904d7dd068bSVaclav Hapla 905d7dd068bSVaclav Hapla Input Parameter: 906d7dd068bSVaclav Hapla . viewer - the PetscViewer 907d7dd068bSVaclav Hapla 908d7dd068bSVaclav Hapla Output Parameter: 909d7dd068bSVaclav Hapla . flg - is timestepping active? 910d7dd068bSVaclav Hapla 911d7dd068bSVaclav Hapla Level: intermediate 912d7dd068bSVaclav Hapla 913d7dd068bSVaclav Hapla Notes: 914d7dd068bSVaclav Hapla See PetscViewerHDF5PushTimestepping() for details. 915d7dd068bSVaclav Hapla 916d7dd068bSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5PopTimestepping(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 917d7dd068bSVaclav Hapla @*/ 918d7dd068bSVaclav Hapla PetscErrorCode PetscViewerHDF5IsTimestepping(PetscViewer viewer, PetscBool *flg) 919d7dd068bSVaclav Hapla { 920d7dd068bSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 921d7dd068bSVaclav Hapla 922d7dd068bSVaclav Hapla PetscFunctionBegin; 923d7dd068bSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 924d7dd068bSVaclav Hapla *flg = hdf5->timestepping; 925d7dd068bSVaclav Hapla PetscFunctionReturn(0); 926d7dd068bSVaclav Hapla } 927d7dd068bSVaclav Hapla 928d7dd068bSVaclav Hapla /*@ 929d7dd068bSVaclav Hapla PetscViewerHDF5IncrementTimestep - Increments current timestep for the HDF5 output. Fields are stacked in time. 930d7dd068bSVaclav Hapla 931d7dd068bSVaclav Hapla Not collective 932d7dd068bSVaclav Hapla 933d7dd068bSVaclav Hapla Input Parameter: 934d7dd068bSVaclav Hapla . viewer - the PetscViewer 935d7dd068bSVaclav Hapla 936d7dd068bSVaclav Hapla Level: intermediate 937d7dd068bSVaclav Hapla 938d7dd068bSVaclav Hapla Notes: 939d7dd068bSVaclav Hapla This can be called only if the viewer is in timestepping mode. See PetscViewerHDF5PushTimestepping() for details. 940d7dd068bSVaclav Hapla 941d7dd068bSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep() 9425c6c1daeSBarry Smith @*/ 9435c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer) 9445c6c1daeSBarry Smith { 9455c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 9465c6c1daeSBarry Smith 9475c6c1daeSBarry Smith PetscFunctionBegin; 9485c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 949d7dd068bSVaclav Hapla if (!hdf5->timestepping) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 9505c6c1daeSBarry Smith ++hdf5->timestep; 9515c6c1daeSBarry Smith PetscFunctionReturn(0); 9525c6c1daeSBarry Smith } 9535c6c1daeSBarry Smith 9543ef9c667SSatish Balay /*@ 955d7dd068bSVaclav Hapla PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. 9565c6c1daeSBarry Smith 957d7dd068bSVaclav Hapla Logically collective 9585c6c1daeSBarry Smith 9595c6c1daeSBarry Smith Input Parameters: 9605c6c1daeSBarry Smith + viewer - the PetscViewer 961d7dd068bSVaclav Hapla - timestep - The timestep 9625c6c1daeSBarry Smith 9635c6c1daeSBarry Smith Level: intermediate 9645c6c1daeSBarry Smith 965d7dd068bSVaclav Hapla Notes: 966d7dd068bSVaclav Hapla This can be called only if the viewer is in timestepping mode. See PetscViewerHDF5PushTimestepping() for details. 967d7dd068bSVaclav Hapla 968d7dd068bSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 9695c6c1daeSBarry Smith @*/ 9705c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep) 9715c6c1daeSBarry Smith { 9725c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 9735c6c1daeSBarry Smith 9745c6c1daeSBarry Smith PetscFunctionBegin; 9755c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 976d7dd068bSVaclav Hapla PetscValidLogicalCollectiveInt(viewer, timestep, 2); 9773ca90d2dSJacob Faibussowitsch if (timestep < 0) SETERRQ1(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestep %" PetscInt_FMT " is negative", timestep); 978d7dd068bSVaclav Hapla if (!hdf5->timestepping) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 9795c6c1daeSBarry Smith hdf5->timestep = timestep; 9805c6c1daeSBarry Smith PetscFunctionReturn(0); 9815c6c1daeSBarry Smith } 9825c6c1daeSBarry Smith 9833ef9c667SSatish Balay /*@ 9845c6c1daeSBarry Smith PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time. 9855c6c1daeSBarry Smith 9865c6c1daeSBarry Smith Not collective 9875c6c1daeSBarry Smith 9885c6c1daeSBarry Smith Input Parameter: 9895c6c1daeSBarry Smith . viewer - the PetscViewer 9905c6c1daeSBarry Smith 9915c6c1daeSBarry Smith Output Parameter: 992d7dd068bSVaclav Hapla . timestep - The timestep 9935c6c1daeSBarry Smith 9945c6c1daeSBarry Smith Level: intermediate 9955c6c1daeSBarry Smith 996d7dd068bSVaclav Hapla Notes: 997d7dd068bSVaclav Hapla This can be called only if the viewer is in the timestepping mode. See PetscViewerHDF5PushTimestepping() for details. 998d7dd068bSVaclav Hapla 999d7dd068bSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep() 10005c6c1daeSBarry Smith @*/ 10015c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep) 10025c6c1daeSBarry Smith { 10035c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 10045c6c1daeSBarry Smith 10055c6c1daeSBarry Smith PetscFunctionBegin; 10065c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1007d7dd068bSVaclav Hapla PetscValidIntPointer(timestep,2); 1008d7dd068bSVaclav Hapla if (!hdf5->timestepping) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 10095c6c1daeSBarry Smith *timestep = hdf5->timestep; 10105c6c1daeSBarry Smith PetscFunctionReturn(0); 10115c6c1daeSBarry Smith } 10125c6c1daeSBarry Smith 101336bce6e8SMatthew G. Knepley /*@C 101436bce6e8SMatthew G. Knepley PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name. 101536bce6e8SMatthew G. Knepley 101636bce6e8SMatthew G. Knepley Not collective 101736bce6e8SMatthew G. Knepley 101836bce6e8SMatthew G. Knepley Input Parameter: 101936bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 102036bce6e8SMatthew G. Knepley 102136bce6e8SMatthew G. Knepley Output Parameter: 102236bce6e8SMatthew G. Knepley . mtype - the MPI datatype (for example MPI_DOUBLE, ...) 102336bce6e8SMatthew G. Knepley 102436bce6e8SMatthew G. Knepley Level: advanced 102536bce6e8SMatthew G. Knepley 102636bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 102736bce6e8SMatthew G. Knepley @*/ 102836bce6e8SMatthew G. Knepley PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype) 102936bce6e8SMatthew G. Knepley { 103036bce6e8SMatthew G. Knepley PetscFunctionBegin; 103136bce6e8SMatthew G. Knepley if (ptype == PETSC_INT) 103236bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 103336bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_LLONG; 103436bce6e8SMatthew G. Knepley #else 103536bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_INT; 103636bce6e8SMatthew G. Knepley #endif 103736bce6e8SMatthew G. Knepley else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE; 103836bce6e8SMatthew G. Knepley else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG; 103936bce6e8SMatthew G. Knepley else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT; 1040c9450970SBarry Smith else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_INT; 1041de3ba262SVaclav Hapla else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_INT; 104236bce6e8SMatthew G. Knepley else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT; 104336bce6e8SMatthew G. Knepley else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR; 104436bce6e8SMatthew G. Knepley else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR; 10457e97c476SMatthew G. Knepley else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1); 104636bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype"); 104736bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 104836bce6e8SMatthew G. Knepley } 104936bce6e8SMatthew G. Knepley 105036bce6e8SMatthew G. Knepley /*@C 105136bce6e8SMatthew G. Knepley PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name 105236bce6e8SMatthew G. Knepley 105336bce6e8SMatthew G. Knepley Not collective 105436bce6e8SMatthew G. Knepley 105536bce6e8SMatthew G. Knepley Input Parameter: 105636bce6e8SMatthew G. Knepley . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...) 105736bce6e8SMatthew G. Knepley 105836bce6e8SMatthew G. Knepley Output Parameter: 105936bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 106036bce6e8SMatthew G. Knepley 106136bce6e8SMatthew G. Knepley Level: advanced 106236bce6e8SMatthew G. Knepley 106336bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 106436bce6e8SMatthew G. Knepley @*/ 106536bce6e8SMatthew G. Knepley PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype) 106636bce6e8SMatthew G. Knepley { 106736bce6e8SMatthew G. Knepley PetscFunctionBegin; 106836bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 106936bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG; 107036bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT; 107136bce6e8SMatthew G. Knepley #else 107236bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT; 107336bce6e8SMatthew G. Knepley #endif 107436bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE; 107536bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG; 107636bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT; 107736bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT; 107836bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR; 107936bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR; 10807e97c476SMatthew G. Knepley else if (htype == H5T_C_S1) *ptype = PETSC_STRING; 108136bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype"); 108236bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 108336bce6e8SMatthew G. Knepley } 108436bce6e8SMatthew G. Knepley 1085df863907SAlex Fikl /*@C 1086b8ee85a1SVaclav Hapla PetscViewerHDF5WriteAttribute - Write an attribute 108736bce6e8SMatthew G. Knepley 1088343a20b2SVaclav Hapla Collective 1089343a20b2SVaclav Hapla 109036bce6e8SMatthew G. Knepley Input Parameters: 109136bce6e8SMatthew G. Knepley + viewer - The HDF5 viewer 10924302210dSVaclav Hapla . parent - The parent dataset/group name 109336bce6e8SMatthew G. Knepley . name - The attribute name 109436bce6e8SMatthew G. Knepley . datatype - The attribute type 109536bce6e8SMatthew G. Knepley - value - The attribute value 109636bce6e8SMatthew G. Knepley 109736bce6e8SMatthew G. Knepley Level: advanced 109836bce6e8SMatthew G. Knepley 10994302210dSVaclav Hapla Notes: 1100343a20b2SVaclav Hapla 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. 11014302210dSVaclav Hapla 1102577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 110336bce6e8SMatthew G. Knepley @*/ 11044302210dSVaclav Hapla PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value) 110536bce6e8SMatthew G. Knepley { 11064302210dSVaclav Hapla char *parentAbsPath; 110760568592SMatthew G. Knepley hid_t h5, dataspace, obj, attribute, dtype; 1108080f144cSVaclav Hapla PetscBool has; 110936bce6e8SMatthew G. Knepley PetscErrorCode ierr; 111036bce6e8SMatthew G. Knepley 111136bce6e8SMatthew G. Knepley PetscFunctionBegin; 11125cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 11134302210dSVaclav Hapla if (parent) PetscValidCharPointer(parent, 2); 1114c57a1a9aSVaclav Hapla PetscValidCharPointer(name, 3); 11154302210dSVaclav Hapla PetscValidLogicalCollectiveEnum(viewer,datatype,4); 1116b7e81f0fSVaclav Hapla PetscValidPointer(value, 5); 11174302210dSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath);CHKERRQ(ierr); 11184302210dSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_TRUE, NULL, NULL);CHKERRQ(ierr); 11194302210dSVaclav Hapla ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has);CHKERRQ(ierr); 112036bce6e8SMatthew G. Knepley ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 11217e97c476SMatthew G. Knepley if (datatype == PETSC_STRING) { 11227e97c476SMatthew G. Knepley size_t len; 11237e97c476SMatthew G. Knepley ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr); 1124729ab6d9SBarry Smith PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 11257e97c476SMatthew G. Knepley } 112636bce6e8SMatthew G. Knepley ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 1127729ab6d9SBarry Smith PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR)); 11284302210dSVaclav Hapla PetscStackCallHDF5Return(obj,H5Oopen,(h5, parentAbsPath, H5P_DEFAULT)); 1129080f144cSVaclav Hapla if (has) { 1130080f144cSVaclav Hapla PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 1131080f144cSVaclav Hapla } else { 113260568592SMatthew G. Knepley PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT)); 1133080f144cSVaclav Hapla } 1134729ab6d9SBarry Smith PetscStackCallHDF5(H5Awrite,(attribute, dtype, value)); 1135729ab6d9SBarry Smith if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype)); 1136729ab6d9SBarry Smith PetscStackCallHDF5(H5Aclose,(attribute)); 113760568592SMatthew G. Knepley PetscStackCallHDF5(H5Oclose,(obj)); 1138729ab6d9SBarry Smith PetscStackCallHDF5(H5Sclose,(dataspace)); 11394302210dSVaclav Hapla ierr = PetscFree(parentAbsPath);CHKERRQ(ierr); 114036bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 114136bce6e8SMatthew G. Knepley } 114236bce6e8SMatthew G. Knepley 1143df863907SAlex Fikl /*@C 1144577e0e04SVaclav Hapla PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name 1145577e0e04SVaclav Hapla 1146343a20b2SVaclav Hapla Collective 1147343a20b2SVaclav Hapla 1148577e0e04SVaclav Hapla Input Parameters: 1149577e0e04SVaclav Hapla + viewer - The HDF5 viewer 1150577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1151577e0e04SVaclav Hapla . name - The attribute name 1152577e0e04SVaclav Hapla . datatype - The attribute type 1153577e0e04SVaclav Hapla - value - The attribute value 1154577e0e04SVaclav Hapla 1155577e0e04SVaclav Hapla Notes: 1156343a20b2SVaclav Hapla 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). 1157577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 1158577e0e04SVaclav Hapla 1159577e0e04SVaclav Hapla Level: advanced 1160577e0e04SVaclav Hapla 1161577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1162577e0e04SVaclav Hapla @*/ 1163577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value) 1164577e0e04SVaclav Hapla { 1165577e0e04SVaclav Hapla PetscErrorCode ierr; 1166577e0e04SVaclav Hapla 1167577e0e04SVaclav Hapla PetscFunctionBegin; 1168577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1169577e0e04SVaclav Hapla PetscValidHeader(obj,2); 1170577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 1171b7e81f0fSVaclav Hapla PetscValidPointer(value,5); 1172577e0e04SVaclav Hapla ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 1173577e0e04SVaclav Hapla ierr = PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr); 1174577e0e04SVaclav Hapla PetscFunctionReturn(0); 1175577e0e04SVaclav Hapla } 1176577e0e04SVaclav Hapla 1177577e0e04SVaclav Hapla /*@C 1178b8ee85a1SVaclav Hapla PetscViewerHDF5ReadAttribute - Read an attribute 1179ad0c61feSMatthew G. Knepley 1180343a20b2SVaclav Hapla Collective 1181343a20b2SVaclav Hapla 1182ad0c61feSMatthew G. Knepley Input Parameters: 1183ad0c61feSMatthew G. Knepley + viewer - The HDF5 viewer 11844302210dSVaclav Hapla . parent - The parent dataset/group name 1185ad0c61feSMatthew G. Knepley . name - The attribute name 1186a2d6be1bSVaclav Hapla . datatype - The attribute type 1187a2d6be1bSVaclav Hapla - defaultValue - The pointer to the default value 1188ad0c61feSMatthew G. Knepley 1189ad0c61feSMatthew G. Knepley Output Parameter: 1190a2d6be1bSVaclav Hapla . value - The pointer to the read HDF5 attribute value 1191ad0c61feSMatthew G. Knepley 1192a2d6be1bSVaclav Hapla Notes: 1193a2d6be1bSVaclav Hapla If defaultValue is NULL and the attribute is not found, an error occurs. 1194a2d6be1bSVaclav Hapla If defaultValue is not NULL and the attribute is not found, defaultValue is copied to value. 1195a2d6be1bSVaclav Hapla The pointers defaultValue and value can be the same; for instance 1196a2d6be1bSVaclav Hapla $ flg = PETSC_FALSE; 1197a2d6be1bSVaclav Hapla $ ierr = PetscViewerHDF5ReadAttribute(viewer,name,"attr",PETSC_BOOL,&flg,&flg);CHKERRQ(ierr); 1198a2d6be1bSVaclav Hapla is valid, but make sure the default value is initialized. 1199a2d6be1bSVaclav Hapla 1200a2d6be1bSVaclav Hapla If the datatype is PETSC_STRING, the output string is newly allocated so one must PetscFree() it when no longer needed. 1201708d4cb3SBarry Smith 1202343a20b2SVaclav Hapla 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. 1203ad0c61feSMatthew G. Knepley 1204343a20b2SVaclav Hapla Level: advanced 12054302210dSVaclav Hapla 1206577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1207ad0c61feSMatthew G. Knepley @*/ 1208d70ec8faSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *defaultValue, void *value) 1209ad0c61feSMatthew G. Knepley { 12104302210dSVaclav Hapla char *parentAbsPath; 1211a2d6be1bSVaclav Hapla hid_t h5, obj, attribute, dtype; 1212969748fdSVaclav Hapla PetscBool has; 1213ad0c61feSMatthew G. Knepley PetscErrorCode ierr; 1214ad0c61feSMatthew G. Knepley 1215ad0c61feSMatthew G. Knepley PetscFunctionBegin; 12165cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 12174302210dSVaclav Hapla if (parent) PetscValidCharPointer(parent, 2); 1218c57a1a9aSVaclav Hapla PetscValidCharPointer(name, 3); 1219a2d6be1bSVaclav Hapla if (defaultValue) PetscValidPointer(defaultValue, 5); 1220a2d6be1bSVaclav Hapla PetscValidPointer(value, 6); 1221a2d6be1bSVaclav Hapla ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 12224302210dSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath);CHKERRQ(ierr); 12234302210dSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, &has, NULL);CHKERRQ(ierr); 12244302210dSVaclav Hapla if (has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has);CHKERRQ(ierr);} 1225a2d6be1bSVaclav Hapla if (!has) { 1226a2d6be1bSVaclav Hapla if (defaultValue) { 1227a2d6be1bSVaclav Hapla if (defaultValue != value) { 1228a2d6be1bSVaclav Hapla if (datatype == PETSC_STRING) { 1229a2d6be1bSVaclav Hapla ierr = PetscStrallocpy(*(char**)defaultValue, (char**)value);CHKERRQ(ierr); 1230a2d6be1bSVaclav Hapla } else { 1231a2d6be1bSVaclav Hapla size_t len; 1232a2d6be1bSVaclav Hapla PetscStackCallHDF5Return(len,H5Tget_size,(dtype)); 1233a2d6be1bSVaclav Hapla ierr = PetscMemcpy(value, defaultValue, len);CHKERRQ(ierr); 1234a2d6be1bSVaclav Hapla } 1235a2d6be1bSVaclav Hapla } 1236a2d6be1bSVaclav Hapla ierr = PetscFree(parentAbsPath);CHKERRQ(ierr); 1237a2d6be1bSVaclav Hapla PetscFunctionReturn(0); 1238a2d6be1bSVaclav Hapla } else SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist and default value not provided", parentAbsPath, name); 1239a2d6be1bSVaclav Hapla } 1240ad0c61feSMatthew G. Knepley ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 12414302210dSVaclav Hapla PetscStackCallHDF5Return(obj,H5Oopen,(h5, parentAbsPath, H5P_DEFAULT)); 124260568592SMatthew G. Knepley PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 1243f0b58479SMatthew G. Knepley if (datatype == PETSC_STRING) { 1244f0b58479SMatthew G. Knepley size_t len; 1245a2d6be1bSVaclav Hapla hid_t atype; 1246e90c5075SMarek Pecha PetscStackCallHDF5Return(atype,H5Aget_type,(attribute)); 1247a2d6be1bSVaclav Hapla PetscStackCallHDF5Return(len,H5Tget_size,(atype)); 1248708d4cb3SBarry Smith ierr = PetscMalloc((len+1) * sizeof(char), value);CHKERRQ(ierr); 1249708d4cb3SBarry Smith PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 1250708d4cb3SBarry Smith PetscStackCallHDF5(H5Aread,(attribute, dtype, *(char**)value)); 1251708d4cb3SBarry Smith } else { 125270efba33SBarry Smith PetscStackCallHDF5(H5Aread,(attribute, dtype, value)); 1253708d4cb3SBarry Smith } 1254729ab6d9SBarry Smith PetscStackCallHDF5(H5Aclose,(attribute)); 1255e90c5075SMarek Pecha /* H5Oclose can be used to close groups, datasets, or committed datatypes */ 1256e90c5075SMarek Pecha PetscStackCallHDF5(H5Oclose,(obj)); 12574302210dSVaclav Hapla ierr = PetscFree(parentAbsPath);CHKERRQ(ierr); 1258ad0c61feSMatthew G. Knepley PetscFunctionReturn(0); 1259ad0c61feSMatthew G. Knepley } 1260ad0c61feSMatthew G. Knepley 1261577e0e04SVaclav Hapla /*@C 1262577e0e04SVaclav Hapla PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name 1263577e0e04SVaclav Hapla 1264343a20b2SVaclav Hapla Collective 1265343a20b2SVaclav Hapla 1266577e0e04SVaclav Hapla Input Parameters: 1267577e0e04SVaclav Hapla + viewer - The HDF5 viewer 1268577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1269577e0e04SVaclav Hapla . name - The attribute name 1270577e0e04SVaclav Hapla - datatype - The attribute type 1271577e0e04SVaclav Hapla 1272577e0e04SVaclav Hapla Output Parameter: 1273577e0e04SVaclav Hapla . value - The attribute value 1274577e0e04SVaclav Hapla 1275577e0e04SVaclav Hapla Notes: 1276577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1277577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 1278577e0e04SVaclav Hapla 1279577e0e04SVaclav Hapla Level: advanced 1280577e0e04SVaclav Hapla 1281577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute() PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1282577e0e04SVaclav Hapla @*/ 1283a2d6be1bSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *defaultValue, void *value) 1284577e0e04SVaclav Hapla { 1285577e0e04SVaclav Hapla PetscErrorCode ierr; 1286577e0e04SVaclav Hapla 1287577e0e04SVaclav Hapla PetscFunctionBegin; 1288577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1289577e0e04SVaclav Hapla PetscValidHeader(obj,2); 1290577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 1291064a246eSJacob Faibussowitsch PetscValidPointer(value, 6); 1292577e0e04SVaclav Hapla ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 1293a2d6be1bSVaclav Hapla ierr = PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, defaultValue, value);CHKERRQ(ierr); 1294577e0e04SVaclav Hapla PetscFunctionReturn(0); 1295577e0e04SVaclav Hapla } 1296577e0e04SVaclav Hapla 1297a75c3fd4SVaclav Hapla PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_) 1298a75c3fd4SVaclav Hapla { 1299a75c3fd4SVaclav Hapla htri_t exists; 1300a75c3fd4SVaclav Hapla hid_t group; 1301a75c3fd4SVaclav Hapla 1302a75c3fd4SVaclav Hapla PetscFunctionBegin; 1303a75c3fd4SVaclav Hapla PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT)); 1304a75c3fd4SVaclav Hapla if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT)); 1305a75c3fd4SVaclav Hapla if (!exists && createGroup) { 1306a75c3fd4SVaclav Hapla PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT)); 1307a75c3fd4SVaclav Hapla PetscStackCallHDF5(H5Gclose,(group)); 1308a75c3fd4SVaclav Hapla exists = PETSC_TRUE; 1309a75c3fd4SVaclav Hapla } 1310a75c3fd4SVaclav Hapla *exists_ = (PetscBool) exists; 1311a75c3fd4SVaclav Hapla PetscFunctionReturn(0); 1312a75c3fd4SVaclav Hapla } 1313a75c3fd4SVaclav Hapla 1314bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype) 13155cdeb108SMatthew G. Knepley { 131690e03537SVaclav Hapla const char rootGroupName[] = "/"; 13175cdeb108SMatthew G. Knepley hid_t h5; 1318e5bf9ebcSVaclav Hapla PetscBool exists=PETSC_FALSE; 131915b861d2SVaclav Hapla PetscInt i; 132015b861d2SVaclav Hapla int n; 132185688ae2SVaclav Hapla char **hierarchy; 132285688ae2SVaclav Hapla char buf[PETSC_MAX_PATH_LEN]=""; 13235cdeb108SMatthew G. Knepley PetscErrorCode ierr; 13245cdeb108SMatthew G. Knepley 13255cdeb108SMatthew G. Knepley PetscFunctionBegin; 13265cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 132790e03537SVaclav Hapla if (name) PetscValidCharPointer(name, 2); 132890e03537SVaclav Hapla else name = rootGroupName; 1329ccd66a83SVaclav Hapla if (has) { 1330064a246eSJacob Faibussowitsch PetscValidBoolPointer(has, 4); 13315cdeb108SMatthew G. Knepley *has = PETSC_FALSE; 1332ccd66a83SVaclav Hapla } 1333ccd66a83SVaclav Hapla if (otype) { 1334064a246eSJacob Faibussowitsch PetscValidIntPointer(otype, 5); 133556cc0592SVaclav Hapla *otype = H5O_TYPE_UNKNOWN; 1336ccd66a83SVaclav Hapla } 1337c0cd0301SJed Brown ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 133885688ae2SVaclav Hapla 133985688ae2SVaclav Hapla /* 134085688ae2SVaclav Hapla Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing. 134185688ae2SVaclav Hapla Hence, each of them needs to be tested separately: 134285688ae2SVaclav Hapla 1) whether it's a valid link 134385688ae2SVaclav Hapla 2) whether this link resolves to an object 134485688ae2SVaclav Hapla See H5Oexists_by_name() documentation. 134585688ae2SVaclav Hapla */ 134685688ae2SVaclav Hapla ierr = PetscStrToArray(name,'/',&n,&hierarchy);CHKERRQ(ierr); 134785688ae2SVaclav Hapla if (!n) { 134885688ae2SVaclav Hapla /* Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */ 1349ccd66a83SVaclav Hapla if (has) *has = PETSC_TRUE; 1350ccd66a83SVaclav Hapla if (otype) *otype = H5O_TYPE_GROUP; 135185688ae2SVaclav Hapla ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 135285688ae2SVaclav Hapla PetscFunctionReturn(0); 135385688ae2SVaclav Hapla } 135485688ae2SVaclav Hapla for (i=0; i<n; i++) { 135585688ae2SVaclav Hapla ierr = PetscStrcat(buf,"/");CHKERRQ(ierr); 135685688ae2SVaclav Hapla ierr = PetscStrcat(buf,hierarchy[i]);CHKERRQ(ierr); 1357a75c3fd4SVaclav Hapla ierr = PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);CHKERRQ(ierr); 1358a75c3fd4SVaclav Hapla if (!exists) break; 135985688ae2SVaclav Hapla } 136085688ae2SVaclav Hapla ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 136185688ae2SVaclav Hapla 136285688ae2SVaclav Hapla /* If the object exists, get its type */ 1363ccd66a83SVaclav Hapla if (exists && otype) { 13645cdeb108SMatthew G. Knepley H5O_info_t info; 13655cdeb108SMatthew G. Knepley 136676276748SVaclav Hapla /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */ 136704633f7fSVaclav Hapla PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT)); 136856cc0592SVaclav Hapla *otype = info.type; 13695cdeb108SMatthew G. Knepley } 1370ccd66a83SVaclav Hapla if (has) *has = exists; 13715cdeb108SMatthew G. Knepley PetscFunctionReturn(0); 13725cdeb108SMatthew G. Knepley } 13735cdeb108SMatthew G. Knepley 13744302210dSVaclav Hapla /*@C 13750a9f38efSVaclav Hapla PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file 13760a9f38efSVaclav Hapla 1377343a20b2SVaclav Hapla Collective 1378343a20b2SVaclav Hapla 13790a9f38efSVaclav Hapla Input Parameters: 13804302210dSVaclav Hapla + viewer - The HDF5 viewer 13814302210dSVaclav Hapla - path - The group path 13820a9f38efSVaclav Hapla 13830a9f38efSVaclav Hapla Output Parameter: 13840a9f38efSVaclav Hapla . has - Flag for group existence 13850a9f38efSVaclav Hapla 13860a9f38efSVaclav Hapla Level: advanced 13870a9f38efSVaclav Hapla 13884302210dSVaclav Hapla Notes: 1389785c443eSVaclav Hapla If path starts with '/', it is taken as an absolute path overriding currently pushed group, else path is relative to the current pushed group. 1390785c443eSVaclav Hapla So NULL or empty path means the current pushed group. 1391343a20b2SVaclav Hapla If path exists but is not a group, PETSC_FALSE is returned. 13924302210dSVaclav Hapla 139389e0ef10SVaclav Hapla .seealso: PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasDataset(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5GetGroup(), PetscViewerHDF5OpenGroup() 13940a9f38efSVaclav Hapla @*/ 13954302210dSVaclav Hapla PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, const char path[], PetscBool *has) 13960a9f38efSVaclav Hapla { 13970a9f38efSVaclav Hapla H5O_type_t type; 13984302210dSVaclav Hapla char *abspath; 13990a9f38efSVaclav Hapla PetscErrorCode ierr; 14000a9f38efSVaclav Hapla 14010a9f38efSVaclav Hapla PetscFunctionBegin; 14020a9f38efSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 14034302210dSVaclav Hapla if (path) PetscValidCharPointer(path,2); 14044302210dSVaclav Hapla PetscValidBoolPointer(has,3); 14054302210dSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, path, &abspath);CHKERRQ(ierr); 14064302210dSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type);CHKERRQ(ierr); 14074302210dSVaclav Hapla *has = (PetscBool)(type == H5O_TYPE_GROUP); 14084302210dSVaclav Hapla ierr = PetscFree(abspath);CHKERRQ(ierr); 14090a9f38efSVaclav Hapla PetscFunctionReturn(0); 14100a9f38efSVaclav Hapla } 14110a9f38efSVaclav Hapla 141289e0ef10SVaclav Hapla /*@C 141389e0ef10SVaclav Hapla PetscViewerHDF5HasDataset - Check whether a given dataset exists in the HDF5 file 141489e0ef10SVaclav Hapla 1415343a20b2SVaclav Hapla Collective 1416343a20b2SVaclav Hapla 141789e0ef10SVaclav Hapla Input Parameters: 141889e0ef10SVaclav Hapla + viewer - The HDF5 viewer 141989e0ef10SVaclav Hapla - path - The dataset path 142089e0ef10SVaclav Hapla 142189e0ef10SVaclav Hapla Output Parameter: 142289e0ef10SVaclav Hapla . has - Flag whether dataset exists 142389e0ef10SVaclav Hapla 142489e0ef10SVaclav Hapla Level: advanced 142589e0ef10SVaclav Hapla 142689e0ef10SVaclav Hapla Notes: 1427343a20b2SVaclav Hapla If path starts with '/', it is taken as an absolute path overriding currently pushed group, else path is relative to the current pushed group. 142889e0ef10SVaclav Hapla If path is NULL or empty, has is set to PETSC_FALSE. 1429343a20b2SVaclav Hapla If path exists but is not a dataset, has is set to PETSC_FALSE as well. 143089e0ef10SVaclav Hapla 143189e0ef10SVaclav Hapla .seealso: PetscViewerHDF5HasObject(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasGroup(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5GetGroup() 143289e0ef10SVaclav Hapla @*/ 143389e0ef10SVaclav Hapla PetscErrorCode PetscViewerHDF5HasDataset(PetscViewer viewer, const char path[], PetscBool *has) 143489e0ef10SVaclav Hapla { 143589e0ef10SVaclav Hapla H5O_type_t type; 143689e0ef10SVaclav Hapla char *abspath; 143789e0ef10SVaclav Hapla PetscErrorCode ierr; 143889e0ef10SVaclav Hapla 143989e0ef10SVaclav Hapla PetscFunctionBegin; 144089e0ef10SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 144189e0ef10SVaclav Hapla if (path) PetscValidCharPointer(path,2); 144289e0ef10SVaclav Hapla PetscValidBoolPointer(has,3); 144389e0ef10SVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, path, &abspath);CHKERRQ(ierr); 144489e0ef10SVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type);CHKERRQ(ierr); 144589e0ef10SVaclav Hapla *has = (PetscBool)(type == H5O_TYPE_DATASET); 144689e0ef10SVaclav Hapla ierr = PetscFree(abspath);CHKERRQ(ierr); 144789e0ef10SVaclav Hapla PetscFunctionReturn(0); 144889e0ef10SVaclav Hapla } 144989e0ef10SVaclav Hapla 14500a9f38efSVaclav Hapla /*@ 1451e3f143f7SVaclav Hapla PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group 1452ecce1506SVaclav Hapla 1453343a20b2SVaclav Hapla Collective 1454343a20b2SVaclav Hapla 1455ecce1506SVaclav Hapla Input Parameters: 1456ecce1506SVaclav Hapla + viewer - The HDF5 viewer 1457ecce1506SVaclav Hapla - obj - The named object 1458ecce1506SVaclav Hapla 1459ecce1506SVaclav Hapla Output Parameter: 146089e0ef10SVaclav Hapla . has - Flag for dataset existence 1461ecce1506SVaclav Hapla 1462e3f143f7SVaclav Hapla Notes: 146389e0ef10SVaclav Hapla If the object is unnamed, an error occurs. 1464343a20b2SVaclav Hapla If the path current_group/object_name exists but is not a dataset, has is set to PETSC_FALSE as well. 1465e3f143f7SVaclav Hapla 1466ecce1506SVaclav Hapla Level: advanced 1467ecce1506SVaclav Hapla 146889e0ef10SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasDataset(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5GetGroup() 1469ecce1506SVaclav Hapla @*/ 1470ecce1506SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has) 1471ecce1506SVaclav Hapla { 147289e0ef10SVaclav Hapla size_t len; 1473ecce1506SVaclav Hapla PetscErrorCode ierr; 1474ecce1506SVaclav Hapla 1475ecce1506SVaclav Hapla PetscFunctionBegin; 1476c57a1a9aSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1477c57a1a9aSVaclav Hapla PetscValidHeader(obj,2); 14784302210dSVaclav Hapla PetscValidBoolPointer(has,3); 147989e0ef10SVaclav Hapla ierr = PetscStrlen(obj->name, &len);CHKERRQ(ierr); 148089e0ef10SVaclav Hapla if (!len) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named"); 148189e0ef10SVaclav Hapla ierr = PetscViewerHDF5HasDataset(viewer, obj->name, has);CHKERRQ(ierr); 1482ecce1506SVaclav Hapla PetscFunctionReturn(0); 1483ecce1506SVaclav Hapla } 1484ecce1506SVaclav Hapla 1485df863907SAlex Fikl /*@C 1486ece83b6cSVaclav Hapla PetscViewerHDF5HasAttribute - Check whether an attribute exists 1487e7bdbf8eSMatthew G. Knepley 1488343a20b2SVaclav Hapla Collective 1489343a20b2SVaclav Hapla 1490e7bdbf8eSMatthew G. Knepley Input Parameters: 1491e7bdbf8eSMatthew G. Knepley + viewer - The HDF5 viewer 14924302210dSVaclav Hapla . parent - The parent dataset/group name 1493e7bdbf8eSMatthew G. Knepley - name - The attribute name 1494e7bdbf8eSMatthew G. Knepley 1495e7bdbf8eSMatthew G. Knepley Output Parameter: 1496e7bdbf8eSMatthew G. Knepley . has - Flag for attribute existence 1497e7bdbf8eSMatthew G. Knepley 1498e7bdbf8eSMatthew G. Knepley Level: advanced 1499e7bdbf8eSMatthew G. Knepley 15004302210dSVaclav Hapla Notes: 1501343a20b2SVaclav Hapla 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. 15024302210dSVaclav Hapla 1503577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1504e7bdbf8eSMatthew G. Knepley @*/ 15054302210dSVaclav Hapla PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 1506e7bdbf8eSMatthew G. Knepley { 15074302210dSVaclav Hapla char *parentAbsPath; 1508e7bdbf8eSMatthew G. Knepley PetscErrorCode ierr; 1509e7bdbf8eSMatthew G. Knepley 1510e7bdbf8eSMatthew G. Knepley PetscFunctionBegin; 15115cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 15124302210dSVaclav Hapla if (parent) PetscValidCharPointer(parent,2); 1513c57a1a9aSVaclav Hapla PetscValidCharPointer(name,3); 15144302210dSVaclav Hapla PetscValidBoolPointer(has,4); 15154302210dSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath);CHKERRQ(ierr); 15164302210dSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, has, NULL);CHKERRQ(ierr); 15174302210dSVaclav Hapla if (*has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, has);CHKERRQ(ierr);} 15184302210dSVaclav Hapla ierr = PetscFree(parentAbsPath);CHKERRQ(ierr); 151906db490cSVaclav Hapla PetscFunctionReturn(0); 152006db490cSVaclav Hapla } 152106db490cSVaclav Hapla 1522577e0e04SVaclav Hapla /*@C 1523577e0e04SVaclav Hapla PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name 1524577e0e04SVaclav Hapla 1525343a20b2SVaclav Hapla Collective 1526343a20b2SVaclav Hapla 1527577e0e04SVaclav Hapla Input Parameters: 1528577e0e04SVaclav Hapla + viewer - The HDF5 viewer 1529577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1530577e0e04SVaclav Hapla - name - The attribute name 1531577e0e04SVaclav Hapla 1532577e0e04SVaclav Hapla Output Parameter: 1533577e0e04SVaclav Hapla . has - Flag for attribute existence 1534577e0e04SVaclav Hapla 1535577e0e04SVaclav Hapla Notes: 1536577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1537577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 1538577e0e04SVaclav Hapla 1539577e0e04SVaclav Hapla Level: advanced 1540577e0e04SVaclav Hapla 1541577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1542577e0e04SVaclav Hapla @*/ 1543577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has) 1544577e0e04SVaclav Hapla { 1545577e0e04SVaclav Hapla PetscErrorCode ierr; 1546577e0e04SVaclav Hapla 1547577e0e04SVaclav Hapla PetscFunctionBegin; 1548577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1549577e0e04SVaclav Hapla PetscValidHeader(obj,2); 1550577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 15514302210dSVaclav Hapla PetscValidBoolPointer(has,4); 1552577e0e04SVaclav Hapla ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 1553577e0e04SVaclav Hapla ierr = PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);CHKERRQ(ierr); 1554577e0e04SVaclav Hapla PetscFunctionReturn(0); 1555577e0e04SVaclav Hapla } 1556577e0e04SVaclav Hapla 155706db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 155806db490cSVaclav Hapla { 155906db490cSVaclav Hapla hid_t h5; 156006db490cSVaclav Hapla htri_t hhas; 156106db490cSVaclav Hapla PetscErrorCode ierr; 156206db490cSVaclav Hapla 156306db490cSVaclav Hapla PetscFunctionBegin; 156406db490cSVaclav Hapla ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 15652f936e54SVaclav Hapla PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT)); 15665cdeb108SMatthew G. Knepley *has = hhas ? PETSC_TRUE : PETSC_FALSE; 1567e7bdbf8eSMatthew G. Knepley PetscFunctionReturn(0); 1568e7bdbf8eSMatthew G. Knepley } 1569e7bdbf8eSMatthew G. Knepley 1570a75e6a4aSMatthew G. Knepley /* 1571a75e6a4aSMatthew G. Knepley The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that 1572a75e6a4aSMatthew G. Knepley is attached to a communicator, in this case the attribute is a PetscViewer. 1573a75e6a4aSMatthew G. Knepley */ 1574d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID; 1575a75e6a4aSMatthew G. Knepley 1576a75e6a4aSMatthew G. Knepley /*@C 1577a75e6a4aSMatthew G. Knepley PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator. 1578a75e6a4aSMatthew G. Knepley 1579d083f849SBarry Smith Collective 1580a75e6a4aSMatthew G. Knepley 1581a75e6a4aSMatthew G. Knepley Input Parameter: 1582a75e6a4aSMatthew G. Knepley . comm - the MPI communicator to share the HDF5 PetscViewer 1583a75e6a4aSMatthew G. Knepley 1584a75e6a4aSMatthew G. Knepley Level: intermediate 1585a75e6a4aSMatthew G. Knepley 1586a75e6a4aSMatthew G. Knepley Options Database Keys: 158710699b91SBarry Smith . -viewer_hdf5_filename <name> - name of the HDF5 file 1588a75e6a4aSMatthew G. Knepley 1589a75e6a4aSMatthew G. Knepley Environmental variables: 159010699b91SBarry Smith . PETSC_VIEWER_HDF5_FILENAME - name of the HDF5 file 1591a75e6a4aSMatthew G. Knepley 1592a75e6a4aSMatthew G. Knepley Notes: 1593a75e6a4aSMatthew G. Knepley Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return 1594a75e6a4aSMatthew G. Knepley an error code. The HDF5 PetscViewer is usually used in the form 1595a75e6a4aSMatthew G. Knepley $ XXXView(XXX object, PETSC_VIEWER_HDF5_(comm)); 1596a75e6a4aSMatthew G. Knepley 1597a75e6a4aSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy() 1598a75e6a4aSMatthew G. Knepley @*/ 1599a75e6a4aSMatthew G. Knepley PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm) 1600a75e6a4aSMatthew G. Knepley { 1601a75e6a4aSMatthew G. Knepley PetscErrorCode ierr; 1602a75e6a4aSMatthew G. Knepley PetscBool flg; 1603a75e6a4aSMatthew G. Knepley PetscViewer viewer; 1604a75e6a4aSMatthew G. Knepley char fname[PETSC_MAX_PATH_LEN]; 1605a75e6a4aSMatthew G. Knepley MPI_Comm ncomm; 1606a75e6a4aSMatthew G. Knepley 1607a75e6a4aSMatthew G. Knepley PetscFunctionBegin; 1608908793a3SLisandro Dalcin ierr = PetscCommDuplicate(comm,&ncomm,NULL);if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 1609a75e6a4aSMatthew G. Knepley if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) { 1610908793a3SLisandro Dalcin ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,NULL); 1611908793a3SLisandro Dalcin if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 1612a75e6a4aSMatthew G. Knepley } 161347435625SJed Brown ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg); 1614908793a3SLisandro Dalcin if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 1615a75e6a4aSMatthew G. Knepley if (!flg) { /* PetscViewer not yet created */ 1616a75e6a4aSMatthew G. Knepley ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg); 16172cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1618a75e6a4aSMatthew G. Knepley if (!flg) { 1619a75e6a4aSMatthew G. Knepley ierr = PetscStrcpy(fname,"output.h5"); 16202cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1621a75e6a4aSMatthew G. Knepley } 1622a75e6a4aSMatthew G. Knepley ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer); 16232cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1624a75e6a4aSMatthew G. Knepley ierr = PetscObjectRegisterDestroy((PetscObject)viewer); 16252cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 162647435625SJed Brown ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer); 1627908793a3SLisandro Dalcin if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 1628a75e6a4aSMatthew G. Knepley } 1629a75e6a4aSMatthew G. Knepley ierr = PetscCommDestroy(&ncomm); 16302cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1631a75e6a4aSMatthew G. Knepley PetscFunctionReturn(viewer); 1632a75e6a4aSMatthew G. Knepley } 1633