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);} 5482ea9b62SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 5582ea9b62SBarry Smith PetscFunctionReturn(0); 5682ea9b62SBarry Smith } 5782ea9b62SBarry Smith 581b793a25SVaclav Hapla static PetscErrorCode PetscViewerView_HDF5(PetscViewer v,PetscViewer viewer) 591b793a25SVaclav Hapla { 601b793a25SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data; 6103fa8834SVaclav Hapla PetscBool flg; 621b793a25SVaclav Hapla PetscErrorCode ierr; 631b793a25SVaclav Hapla 641b793a25SVaclav Hapla PetscFunctionBegin; 651b793a25SVaclav Hapla if (hdf5->filename) { 661b793a25SVaclav Hapla ierr = PetscViewerASCIIPrintf(viewer,"Filename: %s\n",hdf5->filename);CHKERRQ(ierr); 671b793a25SVaclav Hapla } 681b793a25SVaclav Hapla ierr = PetscViewerASCIIPrintf(viewer,"Vectors with blocksize 1 saved as 2D datasets: %s\n",PetscBools[hdf5->basedimension2]);CHKERRQ(ierr); 691b793a25SVaclav Hapla ierr = PetscViewerASCIIPrintf(viewer,"Enforce single precision storage: %s\n",PetscBools[hdf5->spoutput]);CHKERRQ(ierr); 7003fa8834SVaclav Hapla ierr = PetscViewerHDF5GetCollective(v,&flg);CHKERRQ(ierr); 7103fa8834SVaclav Hapla ierr = PetscViewerASCIIPrintf(viewer,"MPI-IO transfer mode: %s\n",flg ? "collective" : "independent");CHKERRQ(ierr); 721b793a25SVaclav Hapla PetscFunctionReturn(0); 731b793a25SVaclav Hapla } 741b793a25SVaclav Hapla 755c6c1daeSBarry Smith static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer) 765c6c1daeSBarry Smith { 775c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data; 785c6c1daeSBarry Smith PetscErrorCode ierr; 795c6c1daeSBarry Smith 805c6c1daeSBarry Smith PetscFunctionBegin; 815c6c1daeSBarry Smith ierr = PetscFree(hdf5->filename);CHKERRQ(ierr); 82729ab6d9SBarry Smith if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 835c6c1daeSBarry Smith PetscFunctionReturn(0); 845c6c1daeSBarry Smith } 855c6c1daeSBarry Smith 86*6226335aSVaclav Hapla static PetscErrorCode PetscViewerFlush_HDF5(PetscViewer viewer) 87*6226335aSVaclav Hapla { 88*6226335aSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data; 89*6226335aSVaclav Hapla 90*6226335aSVaclav Hapla PetscFunctionBegin; 91*6226335aSVaclav Hapla if (hdf5->file_id) PetscStackCallHDF5(H5Fflush,(hdf5->file_id, H5F_SCOPE_LOCAL)); 92*6226335aSVaclav Hapla PetscFunctionReturn(0); 93*6226335aSVaclav Hapla } 94*6226335aSVaclav Hapla 959b1403beSVaclav Hapla static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer) 965c6c1daeSBarry Smith { 975c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 985c6c1daeSBarry Smith PetscErrorCode ierr; 995c6c1daeSBarry Smith 1005c6c1daeSBarry Smith PetscFunctionBegin; 1019c5072fbSVaclav Hapla PetscStackCallHDF5(H5Pclose,(hdf5->dxpl_id)); 1025c6c1daeSBarry Smith ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr); 1035c6c1daeSBarry Smith while (hdf5->groups) { 1047d964842SVaclav Hapla PetscViewerHDF5GroupList *tmp = hdf5->groups->next; 1055c6c1daeSBarry Smith 1065c6c1daeSBarry Smith ierr = PetscFree(hdf5->groups->name);CHKERRQ(ierr); 1075c6c1daeSBarry Smith ierr = PetscFree(hdf5->groups);CHKERRQ(ierr); 1085c6c1daeSBarry Smith hdf5->groups = tmp; 1095c6c1daeSBarry Smith } 1105c6c1daeSBarry Smith ierr = PetscFree(hdf5);CHKERRQ(ierr); 1110b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr); 112d1232d7fSToby Isaac ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr); 1130b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr); 114058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetBaseDimension2_C",NULL);CHKERRQ(ierr); 115058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetSPOutput_C",NULL);CHKERRQ(ierr); 1165c6c1daeSBarry Smith PetscFunctionReturn(0); 1175c6c1daeSBarry Smith } 1185c6c1daeSBarry Smith 1199b1403beSVaclav Hapla static PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type) 1205c6c1daeSBarry Smith { 1215c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 1225c6c1daeSBarry Smith 1235c6c1daeSBarry Smith PetscFunctionBegin; 1245c6c1daeSBarry Smith hdf5->btype = type; 1255c6c1daeSBarry Smith PetscFunctionReturn(0); 1265c6c1daeSBarry Smith } 1275c6c1daeSBarry Smith 1280b62783dSVaclav Hapla static PetscErrorCode PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type) 1290b62783dSVaclav Hapla { 1300b62783dSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 1310b62783dSVaclav Hapla 1320b62783dSVaclav Hapla PetscFunctionBegin; 1330b62783dSVaclav Hapla *type = hdf5->btype; 1340b62783dSVaclav Hapla PetscFunctionReturn(0); 1350b62783dSVaclav Hapla } 1360b62783dSVaclav Hapla 1379b1403beSVaclav Hapla static PetscErrorCode PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg) 13882ea9b62SBarry Smith { 13982ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 14082ea9b62SBarry Smith 14182ea9b62SBarry Smith PetscFunctionBegin; 14282ea9b62SBarry Smith hdf5->basedimension2 = flg; 14382ea9b62SBarry Smith PetscFunctionReturn(0); 14482ea9b62SBarry Smith } 14582ea9b62SBarry Smith 146df863907SAlex Fikl /*@ 14782ea9b62SBarry Smith PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 14882ea9b62SBarry Smith dimension of 2. 14982ea9b62SBarry Smith 15082ea9b62SBarry Smith Logically Collective on PetscViewer 15182ea9b62SBarry Smith 15282ea9b62SBarry Smith Input Parameters: 15382ea9b62SBarry Smith + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 15482ea9b62SBarry 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 15582ea9b62SBarry Smith 15682ea9b62SBarry Smith Options Database: 15782ea9b62SBarry 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 15882ea9b62SBarry Smith 15995452b02SPatrick Sanan Notes: 16095452b02SPatrick 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 16182ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 16282ea9b62SBarry Smith 16382ea9b62SBarry Smith Level: intermediate 16482ea9b62SBarry Smith 16582ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 16682ea9b62SBarry Smith 16782ea9b62SBarry Smith @*/ 16882ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg) 16982ea9b62SBarry Smith { 17082ea9b62SBarry Smith PetscErrorCode ierr; 17182ea9b62SBarry Smith 17282ea9b62SBarry Smith PetscFunctionBegin; 17382ea9b62SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 17482ea9b62SBarry Smith ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 17582ea9b62SBarry Smith PetscFunctionReturn(0); 17682ea9b62SBarry Smith } 17782ea9b62SBarry Smith 178df863907SAlex Fikl /*@ 17982ea9b62SBarry Smith PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 18082ea9b62SBarry Smith dimension of 2. 18182ea9b62SBarry Smith 18282ea9b62SBarry Smith Logically Collective on PetscViewer 18382ea9b62SBarry Smith 18482ea9b62SBarry Smith Input Parameter: 18582ea9b62SBarry Smith . viewer - the PetscViewer, must be of type HDF5 18682ea9b62SBarry Smith 18782ea9b62SBarry Smith Output Parameter: 18882ea9b62SBarry 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 18982ea9b62SBarry Smith 19095452b02SPatrick Sanan Notes: 19195452b02SPatrick 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 19282ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 19382ea9b62SBarry Smith 19482ea9b62SBarry Smith Level: intermediate 19582ea9b62SBarry Smith 19682ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 19782ea9b62SBarry Smith 19882ea9b62SBarry Smith @*/ 19982ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg) 20082ea9b62SBarry Smith { 20182ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 20282ea9b62SBarry Smith 20382ea9b62SBarry Smith PetscFunctionBegin; 20482ea9b62SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 20582ea9b62SBarry Smith *flg = hdf5->basedimension2; 20682ea9b62SBarry Smith PetscFunctionReturn(0); 20782ea9b62SBarry Smith } 20882ea9b62SBarry Smith 2099b1403beSVaclav Hapla static PetscErrorCode PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg) 2109a0502c6SHåkon Strandenes { 2119a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 2129a0502c6SHåkon Strandenes 2139a0502c6SHåkon Strandenes PetscFunctionBegin; 2149a0502c6SHåkon Strandenes hdf5->spoutput = flg; 2159a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2169a0502c6SHåkon Strandenes } 2179a0502c6SHåkon Strandenes 218df863907SAlex Fikl /*@ 2199a0502c6SHåkon Strandenes PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is 2209a0502c6SHåkon Strandenes compiled with double precision PetscReal. 2219a0502c6SHåkon Strandenes 2229a0502c6SHåkon Strandenes Logically Collective on PetscViewer 2239a0502c6SHåkon Strandenes 2249a0502c6SHåkon Strandenes Input Parameters: 2259a0502c6SHåkon Strandenes + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 2269a0502c6SHåkon Strandenes - flg - if PETSC_TRUE the data will be written to disk with single precision 2279a0502c6SHåkon Strandenes 2289a0502c6SHåkon Strandenes Options Database: 2299a0502c6SHåkon Strandenes . -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision 2309a0502c6SHåkon Strandenes 23195452b02SPatrick Sanan Notes: 23295452b02SPatrick Sanan Setting this option does not make any difference if PETSc is compiled with single precision 2339a0502c6SHåkon Strandenes in the first place. It does not affect reading datasets (HDF5 handle this internally). 2349a0502c6SHåkon Strandenes 2359a0502c6SHåkon Strandenes Level: intermediate 2369a0502c6SHåkon Strandenes 2379a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 2389a0502c6SHåkon Strandenes PetscReal 2399a0502c6SHåkon Strandenes 2409a0502c6SHåkon Strandenes @*/ 2419a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg) 2429a0502c6SHåkon Strandenes { 2439a0502c6SHåkon Strandenes PetscErrorCode ierr; 2449a0502c6SHåkon Strandenes 2459a0502c6SHåkon Strandenes PetscFunctionBegin; 2469a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 2479a0502c6SHåkon Strandenes ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 2489a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2499a0502c6SHåkon Strandenes } 2509a0502c6SHåkon Strandenes 251df863907SAlex Fikl /*@ 2529a0502c6SHåkon Strandenes PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is 2539a0502c6SHåkon Strandenes compiled with double precision PetscReal. 2549a0502c6SHåkon Strandenes 2559a0502c6SHåkon Strandenes Logically Collective on PetscViewer 2569a0502c6SHåkon Strandenes 2579a0502c6SHåkon Strandenes Input Parameter: 2589a0502c6SHåkon Strandenes . viewer - the PetscViewer, must be of type HDF5 2599a0502c6SHåkon Strandenes 2609a0502c6SHåkon Strandenes Output Parameter: 2619a0502c6SHåkon Strandenes . flg - if PETSC_TRUE the data will be written to disk with single precision 2629a0502c6SHåkon Strandenes 26395452b02SPatrick Sanan Notes: 26495452b02SPatrick Sanan Setting this option does not make any difference if PETSc is compiled with single precision 2659a0502c6SHåkon Strandenes in the first place. It does not affect reading datasets (HDF5 handle this internally). 2669a0502c6SHåkon Strandenes 2679a0502c6SHåkon Strandenes Level: intermediate 2689a0502c6SHåkon Strandenes 2699a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 2709a0502c6SHåkon Strandenes PetscReal 2719a0502c6SHåkon Strandenes 2729a0502c6SHåkon Strandenes @*/ 2739a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg) 2749a0502c6SHåkon Strandenes { 2759a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 2769a0502c6SHåkon Strandenes 2779a0502c6SHåkon Strandenes PetscFunctionBegin; 2789a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 2799a0502c6SHåkon Strandenes *flg = hdf5->spoutput; 2809a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2819a0502c6SHåkon Strandenes } 2829a0502c6SHåkon Strandenes 283ee8b9732SVaclav Hapla static PetscErrorCode PetscViewerHDF5SetCollective_HDF5(PetscViewer viewer, PetscBool flg) 284ee8b9732SVaclav Hapla { 285ee8b9732SVaclav Hapla PetscFunctionBegin; 286ee8b9732SVaclav Hapla /* H5FD_MPIO_COLLECTIVE is wrong in hdf5 1.10.2, and is the same as H5FD_MPIO_INDEPENDENT in earlier versions 287ee8b9732SVaclav Hapla - see e.g. https://gitlab.cosma.dur.ac.uk/swift/swiftsim/issues/431 */ 288c796909eSBarry Smith #if H5_VERSION_GE(1,10,3) && defined(H5_HAVE_PARALLEL) 289ee8b9732SVaclav Hapla { 290ee8b9732SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 291ee8b9732SVaclav Hapla PetscStackCallHDF5(H5Pset_dxpl_mpio,(hdf5->dxpl_id, flg ? H5FD_MPIO_COLLECTIVE : H5FD_MPIO_INDEPENDENT)); 292ee8b9732SVaclav Hapla } 293ee8b9732SVaclav Hapla #else 294ee8b9732SVaclav Hapla if (flg) { 295ee8b9732SVaclav Hapla PetscErrorCode ierr; 296c796909eSBarry 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); 297ee8b9732SVaclav Hapla } 298ee8b9732SVaclav Hapla #endif 299ee8b9732SVaclav Hapla PetscFunctionReturn(0); 300ee8b9732SVaclav Hapla } 301ee8b9732SVaclav Hapla 302ee8b9732SVaclav Hapla /*@ 303ee8b9732SVaclav Hapla PetscViewerHDF5SetCollective - Use collective MPI-IO transfer mode for HDF5 reads and writes. 304ee8b9732SVaclav Hapla 305ee8b9732SVaclav Hapla Logically Collective; flg must contain common value 306ee8b9732SVaclav Hapla 307ee8b9732SVaclav Hapla Input Parameters: 308ee8b9732SVaclav Hapla + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 309ee8b9732SVaclav Hapla - flg - PETSC_TRUE for collective mode; PETSC_FALSE for independent mode (default) 310ee8b9732SVaclav Hapla 311ee8b9732SVaclav Hapla Options Database: 312ee8b9732SVaclav Hapla . -viewer_hdf5_collective - turns on (true) or off (false) collective transfers 313ee8b9732SVaclav Hapla 314ee8b9732SVaclav Hapla Notes: 315ee8b9732SVaclav Hapla Collective mode gives the MPI-IO layer underneath HDF5 a chance to do some additional collective optimizations and hence can perform better. 31653effed7SBarry 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. 317ee8b9732SVaclav Hapla 318ee8b9732SVaclav Hapla Developer notes: 319ee8b9732SVaclav Hapla In the HDF5 layer, PETSC_TRUE / PETSC_FALSE means H5Pset_dxpl_mpio() is called with H5FD_MPIO_COLLECTIVE / H5FD_MPIO_INDEPENDENT, respectively. 320ee8b9732SVaclav Hapla This in turn means use of MPI_File_{read,write}_all / MPI_File_{read,write} in the MPI-IO layer, respectively. 321ee8b9732SVaclav Hapla See HDF5 documentation and MPI-IO documentation for details. 322ee8b9732SVaclav Hapla 323ee8b9732SVaclav Hapla Level: intermediate 324ee8b9732SVaclav Hapla 325ee8b9732SVaclav Hapla .seealso: PetscViewerHDF5GetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open() 326ee8b9732SVaclav Hapla 327ee8b9732SVaclav Hapla @*/ 328ee8b9732SVaclav Hapla PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer viewer,PetscBool flg) 329ee8b9732SVaclav Hapla { 330ee8b9732SVaclav Hapla PetscErrorCode ierr; 331ee8b9732SVaclav Hapla 332ee8b9732SVaclav Hapla PetscFunctionBegin; 333ee8b9732SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 334ee8b9732SVaclav Hapla PetscValidLogicalCollectiveBool(viewer,flg,2); 335ee8b9732SVaclav Hapla ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetCollective_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 336ee8b9732SVaclav Hapla PetscFunctionReturn(0); 337ee8b9732SVaclav Hapla } 338ee8b9732SVaclav Hapla 339ee8b9732SVaclav Hapla static PetscErrorCode PetscViewerHDF5GetCollective_HDF5(PetscViewer viewer, PetscBool *flg) 340ee8b9732SVaclav Hapla { 341c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL) 342ee8b9732SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 343ee8b9732SVaclav Hapla H5FD_mpio_xfer_t mode; 344c796909eSBarry Smith #endif 345ee8b9732SVaclav Hapla 346ee8b9732SVaclav Hapla PetscFunctionBegin; 347c796909eSBarry Smith #if !defined(H5_HAVE_PARALLEL) 348c796909eSBarry Smith *flg = PETSC_FALSE; 349c796909eSBarry Smith #else 350ee8b9732SVaclav Hapla PetscStackCallHDF5(H5Pget_dxpl_mpio,(hdf5->dxpl_id, &mode)); 351ee8b9732SVaclav Hapla *flg = (mode == H5FD_MPIO_COLLECTIVE) ? PETSC_TRUE : PETSC_FALSE; 352c796909eSBarry Smith #endif 353ee8b9732SVaclav Hapla PetscFunctionReturn(0); 354ee8b9732SVaclav Hapla } 355ee8b9732SVaclav Hapla 356ee8b9732SVaclav Hapla /*@ 357ee8b9732SVaclav Hapla PetscViewerHDF5GetCollective - Return flag whether collective MPI-IO transfer mode is used for HDF5 reads and writes. 358ee8b9732SVaclav Hapla 359ee8b9732SVaclav Hapla Not Collective 360ee8b9732SVaclav Hapla 361ee8b9732SVaclav Hapla Input Parameters: 362ee8b9732SVaclav Hapla . viewer - the HDF5 PetscViewer 363ee8b9732SVaclav Hapla 364ee8b9732SVaclav Hapla Output Parameters: 365ee8b9732SVaclav Hapla . flg - the flag 366ee8b9732SVaclav Hapla 367ee8b9732SVaclav Hapla Level: intermediate 368ee8b9732SVaclav Hapla 369ee8b9732SVaclav Hapla Notes: 370c796909eSBarry 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. 371ee8b9732SVaclav Hapla For more details, see PetscViewerHDF5SetCollective(). 372ee8b9732SVaclav Hapla 373ee8b9732SVaclav Hapla .seealso: PetscViewerHDF5SetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open() 374ee8b9732SVaclav Hapla 375ee8b9732SVaclav Hapla @*/ 376ee8b9732SVaclav Hapla PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer viewer,PetscBool *flg) 377ee8b9732SVaclav Hapla { 378ee8b9732SVaclav Hapla PetscErrorCode ierr; 379ee8b9732SVaclav Hapla 380ee8b9732SVaclav Hapla PetscFunctionBegin; 381ee8b9732SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 382534a8f05SLisandro Dalcin PetscValidBoolPointer(flg,2); 383ee8b9732SVaclav Hapla 384ee8b9732SVaclav Hapla ierr = PetscUseMethod(viewer,"PetscViewerHDF5GetCollective_C",(PetscViewer,PetscBool*),(viewer,flg));CHKERRQ(ierr); 385ee8b9732SVaclav Hapla PetscFunctionReturn(0); 386ee8b9732SVaclav Hapla } 387ee8b9732SVaclav Hapla 3889b1403beSVaclav Hapla static PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[]) 3895c6c1daeSBarry Smith { 3905c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 3915c6c1daeSBarry Smith hid_t plist_id; 3925c6c1daeSBarry Smith PetscErrorCode ierr; 3935c6c1daeSBarry Smith 3945c6c1daeSBarry Smith PetscFunctionBegin; 395f683cc58SToby Isaac if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 396f683cc58SToby Isaac if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);} 3975c6c1daeSBarry Smith ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr); 3985c6c1daeSBarry Smith /* Set up file access property list with parallel I/O access */ 399729ab6d9SBarry Smith PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS)); 400c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL) 401d28bfa55SVaclav Hapla PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), MPI_INFO_NULL)); 402c796909eSBarry Smith #endif 4035c6c1daeSBarry Smith /* Create or open the file collectively */ 4045c6c1daeSBarry Smith switch (hdf5->btype) { 4055c6c1daeSBarry Smith case FILE_MODE_READ: 406729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id)); 4075c6c1daeSBarry Smith break; 4085c6c1daeSBarry Smith case FILE_MODE_APPEND: 4097e4fd573SVaclav Hapla case FILE_MODE_UPDATE: 4107e4fd573SVaclav Hapla { 4117e4fd573SVaclav Hapla PetscBool flg; 4127e4fd573SVaclav Hapla ierr = PetscTestFile(hdf5->filename, 'r', &flg);CHKERRQ(ierr); 4137e4fd573SVaclav Hapla if (flg) PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id)); 4147e4fd573SVaclav Hapla else PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_EXCL, H5P_DEFAULT, plist_id)); 4155c6c1daeSBarry Smith break; 4167e4fd573SVaclav Hapla } 4175c6c1daeSBarry Smith case FILE_MODE_WRITE: 418729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id)); 4195c6c1daeSBarry Smith break; 4207e4fd573SVaclav Hapla case FILE_MODE_UNDEFINED: 4217e4fd573SVaclav Hapla SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 4225c6c1daeSBarry Smith default: 4237e4fd573SVaclav Hapla SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP, "Unsupported file mode %s",PetscFileModes[hdf5->btype]); 4245c6c1daeSBarry Smith } 4255c6c1daeSBarry Smith if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name); 426729ab6d9SBarry Smith PetscStackCallHDF5(H5Pclose,(plist_id)); 4275c6c1daeSBarry Smith PetscFunctionReturn(0); 4285c6c1daeSBarry Smith } 4295c6c1daeSBarry Smith 430d1232d7fSToby Isaac static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name) 431d1232d7fSToby Isaac { 432d1232d7fSToby Isaac PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data; 433d1232d7fSToby Isaac 434d1232d7fSToby Isaac PetscFunctionBegin; 435d1232d7fSToby Isaac *name = vhdf5->filename; 436d1232d7fSToby Isaac PetscFunctionReturn(0); 437d1232d7fSToby Isaac } 438d1232d7fSToby Isaac 439b723ab35SVaclav Hapla static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer) 440b723ab35SVaclav Hapla { 4415dc64a97SVaclav Hapla /* 442b723ab35SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 443b723ab35SVaclav Hapla PetscErrorCode ierr; 4445dc64a97SVaclav Hapla */ 445b723ab35SVaclav Hapla 446b723ab35SVaclav Hapla PetscFunctionBegin; 447b723ab35SVaclav Hapla PetscFunctionReturn(0); 448b723ab35SVaclav Hapla } 449b723ab35SVaclav Hapla 4508556b5ebSBarry Smith /*MC 4518556b5ebSBarry Smith PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file 4528556b5ebSBarry Smith 4538556b5ebSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET, 4548556b5ebSBarry Smith PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING, 4558556b5ebSBarry Smith PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB, 4568556b5ebSBarry Smith PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType() 4578556b5ebSBarry Smith 4581b266c99SBarry Smith Level: beginner 4598556b5ebSBarry Smith M*/ 460d1232d7fSToby Isaac 4618cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v) 4625c6c1daeSBarry Smith { 4635c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5; 4645c6c1daeSBarry Smith PetscErrorCode ierr; 4655c6c1daeSBarry Smith 4665c6c1daeSBarry Smith PetscFunctionBegin; 46799335c34SVaclav Hapla #if !defined(H5_HAVE_PARALLEL) 46899335c34SVaclav Hapla { 46999335c34SVaclav Hapla PetscMPIInt size; 47099335c34SVaclav Hapla ierr = MPI_Comm_size(PetscObjectComm((PetscObject)v), &size);CHKERRMPI(ierr); 47199335c34SVaclav 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)"); 47299335c34SVaclav Hapla } 47399335c34SVaclav Hapla #endif 47499335c34SVaclav Hapla 475b00a9115SJed Brown ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr); 4765c6c1daeSBarry Smith 4775c6c1daeSBarry Smith v->data = (void*) hdf5; 4785c6c1daeSBarry Smith v->ops->destroy = PetscViewerDestroy_HDF5; 47982ea9b62SBarry Smith v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5; 480b723ab35SVaclav Hapla v->ops->setup = PetscViewerSetUp_HDF5; 4811b793a25SVaclav Hapla v->ops->view = PetscViewerView_HDF5; 482*6226335aSVaclav Hapla v->ops->flush = PetscViewerFlush_HDF5; 4837e4fd573SVaclav Hapla hdf5->btype = FILE_MODE_UNDEFINED; 484908793a3SLisandro Dalcin hdf5->filename = NULL; 4855c6c1daeSBarry Smith hdf5->timestep = -1; 4860298fd71SBarry Smith hdf5->groups = NULL; 4875c6c1daeSBarry Smith 4889c5072fbSVaclav Hapla PetscStackCallHDF5Return(hdf5->dxpl_id,H5Pcreate,(H5P_DATASET_XFER)); 4899c5072fbSVaclav Hapla 4900b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr); 491d1232d7fSToby Isaac ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr); 4920b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr); 4930b62783dSVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_HDF5);CHKERRQ(ierr); 49482ea9b62SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr); 4959a0502c6SHåkon Strandenes ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr); 496ee8b9732SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetCollective_C",PetscViewerHDF5SetCollective_HDF5);CHKERRQ(ierr); 497ee8b9732SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetCollective_C",PetscViewerHDF5GetCollective_HDF5);CHKERRQ(ierr); 4985c6c1daeSBarry Smith PetscFunctionReturn(0); 4995c6c1daeSBarry Smith } 5005c6c1daeSBarry Smith 5015c6c1daeSBarry Smith /*@C 5025c6c1daeSBarry Smith PetscViewerHDF5Open - Opens a file for HDF5 input/output. 5035c6c1daeSBarry Smith 504d083f849SBarry Smith Collective 5055c6c1daeSBarry Smith 5065c6c1daeSBarry Smith Input Parameters: 5075c6c1daeSBarry Smith + comm - MPI communicator 5085c6c1daeSBarry Smith . name - name of file 5095c6c1daeSBarry Smith - type - type of file 5105c6c1daeSBarry Smith 5115c6c1daeSBarry Smith Output Parameter: 5125c6c1daeSBarry Smith . hdf5v - PetscViewer for HDF5 input/output to use with the specified file 5135c6c1daeSBarry Smith 51482ea9b62SBarry Smith Options Database: 515a2b725a8SWilliam 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 516a2b725a8SWilliam Gropp - -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal 51782ea9b62SBarry Smith 5185c6c1daeSBarry Smith Level: beginner 5195c6c1daeSBarry Smith 5207e4fd573SVaclav Hapla Notes: 5217e4fd573SVaclav Hapla Reading is always available, regardless of the mode. Available modes are 5227e4fd573SVaclav Hapla + FILE_MODE_READ - open existing HDF5 file for read only access, fail if file does not exist [H5Fopen() with H5F_ACC_RDONLY] 5237e4fd573SVaclav Hapla . FILE_MODE_WRITE - if file exists, fully overwrite it, else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_TRUNC] 5247e4fd573SVaclav 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] 5257e4fd573SVaclav Hapla - FILE_MODE_UPDATE - same as FILE_MODE_APPEND 5267e4fd573SVaclav Hapla 5277e4fd573SVaclav 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. 5287e4fd573SVaclav Hapla 5295c6c1daeSBarry Smith This PetscViewer should be destroyed with PetscViewerDestroy(). 5305c6c1daeSBarry Smith 5316a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(), 5329a0502c6SHåkon Strandenes PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(), 533a56f64adSBarry Smith MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName() 5345c6c1daeSBarry Smith @*/ 5355c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v) 5365c6c1daeSBarry Smith { 5375c6c1daeSBarry Smith PetscErrorCode ierr; 5385c6c1daeSBarry Smith 5395c6c1daeSBarry Smith PetscFunctionBegin; 5405c6c1daeSBarry Smith ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr); 5415c6c1daeSBarry Smith ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr); 5425c6c1daeSBarry Smith ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr); 5435c6c1daeSBarry Smith ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr); 544b663d2d2SVaclav Hapla ierr = PetscViewerSetFromOptions(*hdf5v);CHKERRQ(ierr); 5455c6c1daeSBarry Smith PetscFunctionReturn(0); 5465c6c1daeSBarry Smith } 5475c6c1daeSBarry Smith 5485c6c1daeSBarry Smith /*@C 5495c6c1daeSBarry Smith PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls 5505c6c1daeSBarry Smith 5515c6c1daeSBarry Smith Not collective 5525c6c1daeSBarry Smith 5535c6c1daeSBarry Smith Input Parameter: 5545c6c1daeSBarry Smith . viewer - the PetscViewer 5555c6c1daeSBarry Smith 5565c6c1daeSBarry Smith Output Parameter: 5575c6c1daeSBarry Smith . file_id - The file id 5585c6c1daeSBarry Smith 5595c6c1daeSBarry Smith Level: intermediate 5605c6c1daeSBarry Smith 5615c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open() 5625c6c1daeSBarry Smith @*/ 5635c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id) 5645c6c1daeSBarry Smith { 5655c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 5665c6c1daeSBarry Smith 5675c6c1daeSBarry Smith PetscFunctionBegin; 5685c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 5695c6c1daeSBarry Smith if (file_id) *file_id = hdf5->file_id; 5705c6c1daeSBarry Smith PetscFunctionReturn(0); 5715c6c1daeSBarry Smith } 5725c6c1daeSBarry Smith 5735c6c1daeSBarry Smith /*@C 5745c6c1daeSBarry Smith PetscViewerHDF5PushGroup - Set the current HDF5 group for output 5755c6c1daeSBarry Smith 5765c6c1daeSBarry Smith Not collective 5775c6c1daeSBarry Smith 5785c6c1daeSBarry Smith Input Parameters: 5795c6c1daeSBarry Smith + viewer - the PetscViewer 5805c6c1daeSBarry Smith - name - The group name 5815c6c1daeSBarry Smith 5825c6c1daeSBarry Smith Level: intermediate 5835c6c1daeSBarry Smith 5842e361344SVaclav Hapla Notes: 5852e361344SVaclav Hapla This is designed to mnemonically resemble the Unix cd command. 5862e361344SVaclav 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. 5872e361344SVaclav Hapla . NULL, empty string, or any sequence of all slashes (e.g. "///") is interpreted as the root group "/". 5882e361344SVaclav Hapla - "." means the current group is pushed again. 5892e361344SVaclav Hapla 5902e361344SVaclav Hapla Example: 5912e361344SVaclav Hapla Suppose the current group is "/a". 5922e361344SVaclav Hapla + If name is NULL, empty string, or a sequence of all slashes (e.g. "///"), then the new group will be "/". 5932e361344SVaclav Hapla . If name is ".", then the new group will be "/a". 5942e361344SVaclav Hapla . If name is "b", then the new group will be "/a/b". 5952e361344SVaclav Hapla - If name is "/b", then the new group will be "/b". 5962e361344SVaclav Hapla 5972e361344SVaclav Hapla Developer Notes: 5982e361344SVaclav Hapla The root group "/" is internally stored as NULL. 599820fc2d1SVaclav Hapla 600874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 6015c6c1daeSBarry Smith @*/ 602be7aa430SVaclav Hapla PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char name[]) 6035c6c1daeSBarry Smith { 6045c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 6057d964842SVaclav Hapla PetscViewerHDF5GroupList *groupNode; 6062e361344SVaclav Hapla size_t i,len; 6072e361344SVaclav Hapla char buf[PETSC_MAX_PATH_LEN]; 6082e361344SVaclav Hapla const char *gname; 6095c6c1daeSBarry Smith PetscErrorCode ierr; 6105c6c1daeSBarry Smith 6115c6c1daeSBarry Smith PetscFunctionBegin; 6125c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 613820fc2d1SVaclav Hapla if (name) PetscValidCharPointer(name,2); 614820fc2d1SVaclav Hapla ierr = PetscStrlen(name, &len);CHKERRQ(ierr); 6152e361344SVaclav Hapla gname = NULL; 6162e361344SVaclav Hapla if (len) { 6172e361344SVaclav Hapla if (len == 1 && name[0] == '.') { 6182e361344SVaclav Hapla /* use current name */ 6192e361344SVaclav Hapla gname = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : NULL; 6202e361344SVaclav Hapla } else if (name[0] == '/') { 6212e361344SVaclav Hapla /* absolute */ 6222e361344SVaclav Hapla for (i=1; i<len; i++) { 6232e361344SVaclav Hapla if (name[i] != '/') { 6242e361344SVaclav Hapla gname = name; 6252e361344SVaclav Hapla break; 6262e361344SVaclav Hapla } 6272e361344SVaclav Hapla } 6282e361344SVaclav Hapla } else { 6292e361344SVaclav Hapla /* relative */ 6302e361344SVaclav Hapla const char *parent = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : ""; 6312e361344SVaclav Hapla ierr = PetscSNPrintf(buf, sizeof(buf), "%s/%s", parent, name);CHKERRQ(ierr); 6322e361344SVaclav Hapla gname = buf; 6332e361344SVaclav Hapla } 6342e361344SVaclav Hapla } 63595dccacaSBarry Smith ierr = PetscNew(&groupNode);CHKERRQ(ierr); 6362e361344SVaclav Hapla ierr = PetscStrallocpy(gname, (char**) &groupNode->name);CHKERRQ(ierr); 6375c6c1daeSBarry Smith groupNode->next = hdf5->groups; 6385c6c1daeSBarry Smith hdf5->groups = groupNode; 6395c6c1daeSBarry Smith PetscFunctionReturn(0); 6405c6c1daeSBarry Smith } 6415c6c1daeSBarry Smith 6423ef9c667SSatish Balay /*@ 6435c6c1daeSBarry Smith PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value 6445c6c1daeSBarry Smith 6455c6c1daeSBarry Smith Not collective 6465c6c1daeSBarry Smith 6475c6c1daeSBarry Smith Input Parameter: 6485c6c1daeSBarry Smith . viewer - the PetscViewer 6495c6c1daeSBarry Smith 6505c6c1daeSBarry Smith Level: intermediate 6515c6c1daeSBarry Smith 652874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 6535c6c1daeSBarry Smith @*/ 6545c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer) 6555c6c1daeSBarry Smith { 6565c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 6577d964842SVaclav Hapla PetscViewerHDF5GroupList *groupNode; 6585c6c1daeSBarry Smith PetscErrorCode ierr; 6595c6c1daeSBarry Smith 6605c6c1daeSBarry Smith PetscFunctionBegin; 6615c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 66282f516ccSBarry Smith if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop"); 6635c6c1daeSBarry Smith groupNode = hdf5->groups; 6645c6c1daeSBarry Smith hdf5->groups = hdf5->groups->next; 6655c6c1daeSBarry Smith ierr = PetscFree(groupNode->name);CHKERRQ(ierr); 6665c6c1daeSBarry Smith ierr = PetscFree(groupNode);CHKERRQ(ierr); 6675c6c1daeSBarry Smith PetscFunctionReturn(0); 6685c6c1daeSBarry Smith } 6695c6c1daeSBarry Smith 6705c6c1daeSBarry Smith /*@C 671874270d9SVaclav Hapla PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup(). 672874270d9SVaclav Hapla If none has been assigned, returns NULL. 6735c6c1daeSBarry Smith 6745c6c1daeSBarry Smith Not collective 6755c6c1daeSBarry Smith 6765c6c1daeSBarry Smith Input Parameter: 6775c6c1daeSBarry Smith . viewer - the PetscViewer 6785c6c1daeSBarry Smith 6795c6c1daeSBarry Smith Output Parameter: 6805c6c1daeSBarry Smith . name - The group name 6815c6c1daeSBarry Smith 6825c6c1daeSBarry Smith Level: intermediate 6835c6c1daeSBarry Smith 684874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup() 6855c6c1daeSBarry Smith @*/ 686be7aa430SVaclav Hapla PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char *name[]) 6875c6c1daeSBarry Smith { 6885c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 6895c6c1daeSBarry Smith 6905c6c1daeSBarry Smith PetscFunctionBegin; 6915c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 692c959eef4SJed Brown PetscValidPointer(name,2); 693a297a907SKarl Rupp if (hdf5->groups) *name = hdf5->groups->name; 6940298fd71SBarry Smith else *name = NULL; 6955c6c1daeSBarry Smith PetscFunctionReturn(0); 6965c6c1daeSBarry Smith } 6975c6c1daeSBarry Smith 6988c557b5aSVaclav Hapla /*@ 699874270d9SVaclav Hapla PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(), 700874270d9SVaclav Hapla and return this group's ID and file ID. 701874270d9SVaclav Hapla If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID. 702874270d9SVaclav Hapla 703874270d9SVaclav Hapla Not collective 704874270d9SVaclav Hapla 705874270d9SVaclav Hapla Input Parameter: 706874270d9SVaclav Hapla . viewer - the PetscViewer 707874270d9SVaclav Hapla 708d8d19677SJose E. Roman Output Parameters: 709874270d9SVaclav Hapla + fileId - The HDF5 file ID 710874270d9SVaclav Hapla - groupId - The HDF5 group ID 711874270d9SVaclav Hapla 712e74616beSVaclav Hapla Notes: 713e74616beSVaclav Hapla If the viewer is writable, the group is created if it doesn't exist yet. 714e74616beSVaclav Hapla 715874270d9SVaclav Hapla Level: intermediate 716874270d9SVaclav Hapla 717874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 718874270d9SVaclav Hapla @*/ 71954dbf706SVaclav Hapla PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId) 72054dbf706SVaclav Hapla { 72190e03537SVaclav Hapla hid_t file_id; 72290e03537SVaclav Hapla H5O_type_t type; 72376d59af2SVaclav Hapla const char *groupName = NULL, *fileName = NULL; 72476d59af2SVaclav Hapla PetscBool writable, has; 72554dbf706SVaclav Hapla PetscErrorCode ierr; 72654dbf706SVaclav Hapla 72754dbf706SVaclav Hapla PetscFunctionBegin; 72876d59af2SVaclav Hapla ierr = PetscViewerWritable(viewer, &writable);CHKERRQ(ierr); 72954dbf706SVaclav Hapla ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr); 73076d59af2SVaclav Hapla ierr = PetscViewerFileGetName(viewer, &fileName);CHKERRQ(ierr); 73154dbf706SVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr); 73276d59af2SVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, groupName, writable, &has, &type);CHKERRQ(ierr); 73376d59af2SVaclav Hapla if (!has) { 73476d59af2SVaclav 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); 73576d59af2SVaclav Hapla else SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_LIB, "HDF5 failed to create group %s although file %s is open for writing", groupName, fileName); 73676d59af2SVaclav Hapla } 73776d59af2SVaclav 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); 73890e03537SVaclav Hapla PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT)); 73954dbf706SVaclav Hapla *fileId = file_id; 74054dbf706SVaclav Hapla PetscFunctionReturn(0); 74154dbf706SVaclav Hapla } 74254dbf706SVaclav Hapla 7433ef9c667SSatish Balay /*@ 744d7dd068bSVaclav Hapla PetscViewerHDF5PushTimestepping - Activate timestepping mode for subsequent HDF5 reading and writing. 7455c6c1daeSBarry Smith 7465c6c1daeSBarry Smith Not collective 7475c6c1daeSBarry Smith 7485c6c1daeSBarry Smith Input Parameter: 7495c6c1daeSBarry Smith . viewer - the PetscViewer 7505c6c1daeSBarry Smith 7515c6c1daeSBarry Smith Level: intermediate 7525c6c1daeSBarry Smith 753d7dd068bSVaclav Hapla Notes: 754d7dd068bSVaclav Hapla On first PetscViewerHDF5PushTimestepping(), the initial time step is set to 0. 755d7dd068bSVaclav Hapla Next timesteps can then be set using PetscViewerHDF5IncrementTimestep() or PetscViewerHDF5SetTimestep(). 756d7dd068bSVaclav 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()]. 757d7dd068bSVaclav Hapla Use PetscViewerHDF5PopTimestepping() to deactivate timestepping mode; calling it by the end of the program is NOT mandatory. 758d7dd068bSVaclav Hapla Current timestep is remembered between PetscViewerHDF5PopTimestepping() and the next PetscViewerHDF5PushTimestepping(). 759d7dd068bSVaclav Hapla 760d7dd068bSVaclav Hapla If a dataset was stored with timestepping, it can be loaded only in the timestepping mode again. 761d7dd068bSVaclav Hapla Loading a timestepped dataset with timestepping disabled, or vice-versa results in an error. 762d7dd068bSVaclav Hapla 763d7dd068bSVaclav Hapla Developer notes: 764d7dd068bSVaclav Hapla Timestepped HDF5 dataset has an extra dimension and attribute "timestepping" set to true. 765d7dd068bSVaclav Hapla 766d7dd068bSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PopTimestepping(), PetscViewerHDF5IsTimestepping(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 767d7dd068bSVaclav Hapla @*/ 768d7dd068bSVaclav Hapla PetscErrorCode PetscViewerHDF5PushTimestepping(PetscViewer viewer) 769d7dd068bSVaclav Hapla { 770d7dd068bSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 771d7dd068bSVaclav Hapla 772d7dd068bSVaclav Hapla PetscFunctionBegin; 773d7dd068bSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 774d7dd068bSVaclav Hapla if (hdf5->timestepping) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping is already pushed"); 775d7dd068bSVaclav Hapla hdf5->timestepping = PETSC_TRUE; 776d7dd068bSVaclav Hapla if (hdf5->timestep < 0) hdf5->timestep = 0; 777d7dd068bSVaclav Hapla PetscFunctionReturn(0); 778d7dd068bSVaclav Hapla } 779d7dd068bSVaclav Hapla 780d7dd068bSVaclav Hapla /*@ 781d7dd068bSVaclav Hapla PetscViewerHDF5PopTimestepping - Deactivate timestepping mode for subsequent HDF5 reading and writing. 782d7dd068bSVaclav Hapla 783d7dd068bSVaclav Hapla Not collective 784d7dd068bSVaclav Hapla 785d7dd068bSVaclav Hapla Input Parameter: 786d7dd068bSVaclav Hapla . viewer - the PetscViewer 787d7dd068bSVaclav Hapla 788d7dd068bSVaclav Hapla Level: intermediate 789d7dd068bSVaclav Hapla 790d7dd068bSVaclav Hapla Notes: 791d7dd068bSVaclav Hapla See PetscViewerHDF5PushTimestepping() for details. 792d7dd068bSVaclav Hapla 793d7dd068bSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5IsTimestepping(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 794d7dd068bSVaclav Hapla @*/ 795d7dd068bSVaclav Hapla PetscErrorCode PetscViewerHDF5PopTimestepping(PetscViewer viewer) 796d7dd068bSVaclav Hapla { 797d7dd068bSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 798d7dd068bSVaclav Hapla 799d7dd068bSVaclav Hapla PetscFunctionBegin; 800d7dd068bSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 801d7dd068bSVaclav Hapla if (!hdf5->timestepping) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 802d7dd068bSVaclav Hapla hdf5->timestepping = PETSC_FALSE; 803d7dd068bSVaclav Hapla PetscFunctionReturn(0); 804d7dd068bSVaclav Hapla } 805d7dd068bSVaclav Hapla 806d7dd068bSVaclav Hapla /*@ 807d7dd068bSVaclav Hapla PetscViewerHDF5IsTimestepping - Ask the viewer whether it is in timestepping mode currently. 808d7dd068bSVaclav Hapla 809d7dd068bSVaclav Hapla Not collective 810d7dd068bSVaclav Hapla 811d7dd068bSVaclav Hapla Input Parameter: 812d7dd068bSVaclav Hapla . viewer - the PetscViewer 813d7dd068bSVaclav Hapla 814d7dd068bSVaclav Hapla Output Parameter: 815d7dd068bSVaclav Hapla . flg - is timestepping active? 816d7dd068bSVaclav Hapla 817d7dd068bSVaclav Hapla Level: intermediate 818d7dd068bSVaclav Hapla 819d7dd068bSVaclav Hapla Notes: 820d7dd068bSVaclav Hapla See PetscViewerHDF5PushTimestepping() for details. 821d7dd068bSVaclav Hapla 822d7dd068bSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5PopTimestepping(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 823d7dd068bSVaclav Hapla @*/ 824d7dd068bSVaclav Hapla PetscErrorCode PetscViewerHDF5IsTimestepping(PetscViewer viewer, PetscBool *flg) 825d7dd068bSVaclav Hapla { 826d7dd068bSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 827d7dd068bSVaclav Hapla 828d7dd068bSVaclav Hapla PetscFunctionBegin; 829d7dd068bSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 830d7dd068bSVaclav Hapla *flg = hdf5->timestepping; 831d7dd068bSVaclav Hapla PetscFunctionReturn(0); 832d7dd068bSVaclav Hapla } 833d7dd068bSVaclav Hapla 834d7dd068bSVaclav Hapla /*@ 835d7dd068bSVaclav Hapla PetscViewerHDF5IncrementTimestep - Increments current timestep for the HDF5 output. Fields are stacked in time. 836d7dd068bSVaclav Hapla 837d7dd068bSVaclav Hapla Not collective 838d7dd068bSVaclav Hapla 839d7dd068bSVaclav Hapla Input Parameter: 840d7dd068bSVaclav Hapla . viewer - the PetscViewer 841d7dd068bSVaclav Hapla 842d7dd068bSVaclav Hapla Level: intermediate 843d7dd068bSVaclav Hapla 844d7dd068bSVaclav Hapla Notes: 845d7dd068bSVaclav Hapla This can be called only if the viewer is in timestepping mode. See PetscViewerHDF5PushTimestepping() for details. 846d7dd068bSVaclav Hapla 847d7dd068bSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep() 8485c6c1daeSBarry Smith @*/ 8495c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer) 8505c6c1daeSBarry Smith { 8515c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 8525c6c1daeSBarry Smith 8535c6c1daeSBarry Smith PetscFunctionBegin; 8545c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 855d7dd068bSVaclav Hapla if (!hdf5->timestepping) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 8565c6c1daeSBarry Smith ++hdf5->timestep; 8575c6c1daeSBarry Smith PetscFunctionReturn(0); 8585c6c1daeSBarry Smith } 8595c6c1daeSBarry Smith 8603ef9c667SSatish Balay /*@ 861d7dd068bSVaclav Hapla PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. 8625c6c1daeSBarry Smith 863d7dd068bSVaclav Hapla Logically collective 8645c6c1daeSBarry Smith 8655c6c1daeSBarry Smith Input Parameters: 8665c6c1daeSBarry Smith + viewer - the PetscViewer 867d7dd068bSVaclav Hapla - timestep - The timestep 8685c6c1daeSBarry Smith 8695c6c1daeSBarry Smith Level: intermediate 8705c6c1daeSBarry Smith 871d7dd068bSVaclav Hapla Notes: 872d7dd068bSVaclav Hapla This can be called only if the viewer is in timestepping mode. See PetscViewerHDF5PushTimestepping() for details. 873d7dd068bSVaclav Hapla 874d7dd068bSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 8755c6c1daeSBarry Smith @*/ 8765c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep) 8775c6c1daeSBarry Smith { 8785c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 8795c6c1daeSBarry Smith 8805c6c1daeSBarry Smith PetscFunctionBegin; 8815c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 882d7dd068bSVaclav Hapla PetscValidLogicalCollectiveInt(viewer, timestep, 2); 883d7dd068bSVaclav Hapla if (timestep < 0) SETERRQ1(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestep %D is negative", timestep); 884d7dd068bSVaclav Hapla if (!hdf5->timestepping) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 8855c6c1daeSBarry Smith hdf5->timestep = timestep; 8865c6c1daeSBarry Smith PetscFunctionReturn(0); 8875c6c1daeSBarry Smith } 8885c6c1daeSBarry Smith 8893ef9c667SSatish Balay /*@ 8905c6c1daeSBarry Smith PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time. 8915c6c1daeSBarry Smith 8925c6c1daeSBarry Smith Not collective 8935c6c1daeSBarry Smith 8945c6c1daeSBarry Smith Input Parameter: 8955c6c1daeSBarry Smith . viewer - the PetscViewer 8965c6c1daeSBarry Smith 8975c6c1daeSBarry Smith Output Parameter: 898d7dd068bSVaclav Hapla . timestep - The timestep 8995c6c1daeSBarry Smith 9005c6c1daeSBarry Smith Level: intermediate 9015c6c1daeSBarry Smith 902d7dd068bSVaclav Hapla Notes: 903d7dd068bSVaclav Hapla This can be called only if the viewer is in the timestepping mode. See PetscViewerHDF5PushTimestepping() for details. 904d7dd068bSVaclav Hapla 905d7dd068bSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep() 9065c6c1daeSBarry Smith @*/ 9075c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep) 9085c6c1daeSBarry Smith { 9095c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 9105c6c1daeSBarry Smith 9115c6c1daeSBarry Smith PetscFunctionBegin; 9125c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 913d7dd068bSVaclav Hapla PetscValidIntPointer(timestep,2); 914d7dd068bSVaclav Hapla if (!hdf5->timestepping) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 9155c6c1daeSBarry Smith *timestep = hdf5->timestep; 9165c6c1daeSBarry Smith PetscFunctionReturn(0); 9175c6c1daeSBarry Smith } 9185c6c1daeSBarry Smith 91936bce6e8SMatthew G. Knepley /*@C 92036bce6e8SMatthew G. Knepley PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name. 92136bce6e8SMatthew G. Knepley 92236bce6e8SMatthew G. Knepley Not collective 92336bce6e8SMatthew G. Knepley 92436bce6e8SMatthew G. Knepley Input Parameter: 92536bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 92636bce6e8SMatthew G. Knepley 92736bce6e8SMatthew G. Knepley Output Parameter: 92836bce6e8SMatthew G. Knepley . mtype - the MPI datatype (for example MPI_DOUBLE, ...) 92936bce6e8SMatthew G. Knepley 93036bce6e8SMatthew G. Knepley Level: advanced 93136bce6e8SMatthew G. Knepley 93236bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 93336bce6e8SMatthew G. Knepley @*/ 93436bce6e8SMatthew G. Knepley PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype) 93536bce6e8SMatthew G. Knepley { 93636bce6e8SMatthew G. Knepley PetscFunctionBegin; 93736bce6e8SMatthew G. Knepley if (ptype == PETSC_INT) 93836bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 93936bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_LLONG; 94036bce6e8SMatthew G. Knepley #else 94136bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_INT; 94236bce6e8SMatthew G. Knepley #endif 94336bce6e8SMatthew G. Knepley else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE; 94436bce6e8SMatthew G. Knepley else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG; 94536bce6e8SMatthew G. Knepley else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT; 946c9450970SBarry Smith else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_INT; 947de3ba262SVaclav Hapla else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_INT; 94836bce6e8SMatthew G. Knepley else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT; 94936bce6e8SMatthew G. Knepley else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR; 95036bce6e8SMatthew G. Knepley else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR; 9517e97c476SMatthew G. Knepley else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1); 95236bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype"); 95336bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 95436bce6e8SMatthew G. Knepley } 95536bce6e8SMatthew G. Knepley 95636bce6e8SMatthew G. Knepley /*@C 95736bce6e8SMatthew G. Knepley PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name 95836bce6e8SMatthew G. Knepley 95936bce6e8SMatthew G. Knepley Not collective 96036bce6e8SMatthew G. Knepley 96136bce6e8SMatthew G. Knepley Input Parameter: 96236bce6e8SMatthew G. Knepley . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...) 96336bce6e8SMatthew G. Knepley 96436bce6e8SMatthew G. Knepley Output Parameter: 96536bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 96636bce6e8SMatthew G. Knepley 96736bce6e8SMatthew G. Knepley Level: advanced 96836bce6e8SMatthew G. Knepley 96936bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 97036bce6e8SMatthew G. Knepley @*/ 97136bce6e8SMatthew G. Knepley PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype) 97236bce6e8SMatthew G. Knepley { 97336bce6e8SMatthew G. Knepley PetscFunctionBegin; 97436bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 97536bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG; 97636bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT; 97736bce6e8SMatthew G. Knepley #else 97836bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT; 97936bce6e8SMatthew G. Knepley #endif 98036bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE; 98136bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG; 98236bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT; 98336bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT; 98436bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR; 98536bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR; 9867e97c476SMatthew G. Knepley else if (htype == H5T_C_S1) *ptype = PETSC_STRING; 98736bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype"); 98836bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 98936bce6e8SMatthew G. Knepley } 99036bce6e8SMatthew G. Knepley 991df863907SAlex Fikl /*@C 992b8ee85a1SVaclav Hapla PetscViewerHDF5WriteAttribute - Write an attribute 99336bce6e8SMatthew G. Knepley 994343a20b2SVaclav Hapla Collective 995343a20b2SVaclav Hapla 99636bce6e8SMatthew G. Knepley Input Parameters: 99736bce6e8SMatthew G. Knepley + viewer - The HDF5 viewer 9984302210dSVaclav Hapla . parent - The parent dataset/group name 99936bce6e8SMatthew G. Knepley . name - The attribute name 100036bce6e8SMatthew G. Knepley . datatype - The attribute type 100136bce6e8SMatthew G. Knepley - value - The attribute value 100236bce6e8SMatthew G. Knepley 100336bce6e8SMatthew G. Knepley Level: advanced 100436bce6e8SMatthew G. Knepley 10054302210dSVaclav Hapla Notes: 1006343a20b2SVaclav 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. 10074302210dSVaclav Hapla 1008577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 100936bce6e8SMatthew G. Knepley @*/ 10104302210dSVaclav Hapla PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value) 101136bce6e8SMatthew G. Knepley { 10124302210dSVaclav Hapla char *parentAbsPath; 101360568592SMatthew G. Knepley hid_t h5, dataspace, obj, attribute, dtype; 1014080f144cSVaclav Hapla PetscBool has; 101536bce6e8SMatthew G. Knepley PetscErrorCode ierr; 101636bce6e8SMatthew G. Knepley 101736bce6e8SMatthew G. Knepley PetscFunctionBegin; 10185cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 10194302210dSVaclav Hapla if (parent) PetscValidCharPointer(parent, 2); 1020c57a1a9aSVaclav Hapla PetscValidCharPointer(name, 3); 10214302210dSVaclav Hapla PetscValidLogicalCollectiveEnum(viewer,datatype,4); 1022b7e81f0fSVaclav Hapla PetscValidPointer(value, 5); 10234302210dSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath);CHKERRQ(ierr); 10244302210dSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_TRUE, NULL, NULL);CHKERRQ(ierr); 10254302210dSVaclav Hapla ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has);CHKERRQ(ierr); 102636bce6e8SMatthew G. Knepley ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 10277e97c476SMatthew G. Knepley if (datatype == PETSC_STRING) { 10287e97c476SMatthew G. Knepley size_t len; 10297e97c476SMatthew G. Knepley ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr); 1030729ab6d9SBarry Smith PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 10317e97c476SMatthew G. Knepley } 103236bce6e8SMatthew G. Knepley ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 1033729ab6d9SBarry Smith PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR)); 10344302210dSVaclav Hapla PetscStackCallHDF5Return(obj,H5Oopen,(h5, parentAbsPath, H5P_DEFAULT)); 1035080f144cSVaclav Hapla if (has) { 1036080f144cSVaclav Hapla PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 1037080f144cSVaclav Hapla } else { 103860568592SMatthew G. Knepley PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT)); 1039080f144cSVaclav Hapla } 1040729ab6d9SBarry Smith PetscStackCallHDF5(H5Awrite,(attribute, dtype, value)); 1041729ab6d9SBarry Smith if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype)); 1042729ab6d9SBarry Smith PetscStackCallHDF5(H5Aclose,(attribute)); 104360568592SMatthew G. Knepley PetscStackCallHDF5(H5Oclose,(obj)); 1044729ab6d9SBarry Smith PetscStackCallHDF5(H5Sclose,(dataspace)); 10454302210dSVaclav Hapla ierr = PetscFree(parentAbsPath);CHKERRQ(ierr); 104636bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 104736bce6e8SMatthew G. Knepley } 104836bce6e8SMatthew G. Knepley 1049df863907SAlex Fikl /*@C 1050577e0e04SVaclav Hapla PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name 1051577e0e04SVaclav Hapla 1052343a20b2SVaclav Hapla Collective 1053343a20b2SVaclav Hapla 1054577e0e04SVaclav Hapla Input Parameters: 1055577e0e04SVaclav Hapla + viewer - The HDF5 viewer 1056577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1057577e0e04SVaclav Hapla . name - The attribute name 1058577e0e04SVaclav Hapla . datatype - The attribute type 1059577e0e04SVaclav Hapla - value - The attribute value 1060577e0e04SVaclav Hapla 1061577e0e04SVaclav Hapla Notes: 1062343a20b2SVaclav 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). 1063577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 1064577e0e04SVaclav Hapla 1065577e0e04SVaclav Hapla Level: advanced 1066577e0e04SVaclav Hapla 1067577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1068577e0e04SVaclav Hapla @*/ 1069577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value) 1070577e0e04SVaclav Hapla { 1071577e0e04SVaclav Hapla PetscErrorCode ierr; 1072577e0e04SVaclav Hapla 1073577e0e04SVaclav Hapla PetscFunctionBegin; 1074577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1075577e0e04SVaclav Hapla PetscValidHeader(obj,2); 1076577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 1077b7e81f0fSVaclav Hapla PetscValidPointer(value,5); 1078577e0e04SVaclav Hapla ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 1079577e0e04SVaclav Hapla ierr = PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr); 1080577e0e04SVaclav Hapla PetscFunctionReturn(0); 1081577e0e04SVaclav Hapla } 1082577e0e04SVaclav Hapla 1083577e0e04SVaclav Hapla /*@C 1084b8ee85a1SVaclav Hapla PetscViewerHDF5ReadAttribute - Read an attribute 1085ad0c61feSMatthew G. Knepley 1086343a20b2SVaclav Hapla Collective 1087343a20b2SVaclav Hapla 1088ad0c61feSMatthew G. Knepley Input Parameters: 1089ad0c61feSMatthew G. Knepley + viewer - The HDF5 viewer 10904302210dSVaclav Hapla . parent - The parent dataset/group name 1091ad0c61feSMatthew G. Knepley . name - The attribute name 1092a2d6be1bSVaclav Hapla . datatype - The attribute type 1093a2d6be1bSVaclav Hapla - defaultValue - The pointer to the default value 1094ad0c61feSMatthew G. Knepley 1095ad0c61feSMatthew G. Knepley Output Parameter: 1096a2d6be1bSVaclav Hapla . value - The pointer to the read HDF5 attribute value 1097ad0c61feSMatthew G. Knepley 1098a2d6be1bSVaclav Hapla Notes: 1099a2d6be1bSVaclav Hapla If defaultValue is NULL and the attribute is not found, an error occurs. 1100a2d6be1bSVaclav Hapla If defaultValue is not NULL and the attribute is not found, defaultValue is copied to value. 1101a2d6be1bSVaclav Hapla The pointers defaultValue and value can be the same; for instance 1102a2d6be1bSVaclav Hapla $ flg = PETSC_FALSE; 1103a2d6be1bSVaclav Hapla $ ierr = PetscViewerHDF5ReadAttribute(viewer,name,"attr",PETSC_BOOL,&flg,&flg);CHKERRQ(ierr); 1104a2d6be1bSVaclav Hapla is valid, but make sure the default value is initialized. 1105a2d6be1bSVaclav Hapla 1106a2d6be1bSVaclav Hapla If the datatype is PETSC_STRING, the output string is newly allocated so one must PetscFree() it when no longer needed. 1107708d4cb3SBarry Smith 1108343a20b2SVaclav 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. 1109ad0c61feSMatthew G. Knepley 1110343a20b2SVaclav Hapla Level: advanced 11114302210dSVaclav Hapla 1112577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1113ad0c61feSMatthew G. Knepley @*/ 1114d70ec8faSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *defaultValue, void *value) 1115ad0c61feSMatthew G. Knepley { 11164302210dSVaclav Hapla char *parentAbsPath; 1117a2d6be1bSVaclav Hapla hid_t h5, obj, attribute, dtype; 1118969748fdSVaclav Hapla PetscBool has; 1119ad0c61feSMatthew G. Knepley PetscErrorCode ierr; 1120ad0c61feSMatthew G. Knepley 1121ad0c61feSMatthew G. Knepley PetscFunctionBegin; 11225cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 11234302210dSVaclav Hapla if (parent) PetscValidCharPointer(parent, 2); 1124c57a1a9aSVaclav Hapla PetscValidCharPointer(name, 3); 1125a2d6be1bSVaclav Hapla if (defaultValue) PetscValidPointer(defaultValue, 5); 1126a2d6be1bSVaclav Hapla PetscValidPointer(value, 6); 1127a2d6be1bSVaclav Hapla ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 11284302210dSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath);CHKERRQ(ierr); 11294302210dSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, &has, NULL);CHKERRQ(ierr); 11304302210dSVaclav Hapla if (has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has);CHKERRQ(ierr);} 1131a2d6be1bSVaclav Hapla if (!has) { 1132a2d6be1bSVaclav Hapla if (defaultValue) { 1133a2d6be1bSVaclav Hapla if (defaultValue != value) { 1134a2d6be1bSVaclav Hapla if (datatype == PETSC_STRING) { 1135a2d6be1bSVaclav Hapla ierr = PetscStrallocpy(*(char**)defaultValue, (char**)value);CHKERRQ(ierr); 1136a2d6be1bSVaclav Hapla } else { 1137a2d6be1bSVaclav Hapla size_t len; 1138a2d6be1bSVaclav Hapla PetscStackCallHDF5Return(len,H5Tget_size,(dtype)); 1139a2d6be1bSVaclav Hapla ierr = PetscMemcpy(value, defaultValue, len);CHKERRQ(ierr); 1140a2d6be1bSVaclav Hapla } 1141a2d6be1bSVaclav Hapla } 1142a2d6be1bSVaclav Hapla ierr = PetscFree(parentAbsPath);CHKERRQ(ierr); 1143a2d6be1bSVaclav Hapla PetscFunctionReturn(0); 1144a2d6be1bSVaclav Hapla } else SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist and default value not provided", parentAbsPath, name); 1145a2d6be1bSVaclav Hapla } 1146ad0c61feSMatthew G. Knepley ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 11474302210dSVaclav Hapla PetscStackCallHDF5Return(obj,H5Oopen,(h5, parentAbsPath, H5P_DEFAULT)); 114860568592SMatthew G. Knepley PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 1149f0b58479SMatthew G. Knepley if (datatype == PETSC_STRING) { 1150f0b58479SMatthew G. Knepley size_t len; 1151a2d6be1bSVaclav Hapla hid_t atype; 1152e90c5075SMarek Pecha PetscStackCallHDF5Return(atype,H5Aget_type,(attribute)); 1153a2d6be1bSVaclav Hapla PetscStackCallHDF5Return(len,H5Tget_size,(atype)); 1154708d4cb3SBarry Smith ierr = PetscMalloc((len+1) * sizeof(char), value);CHKERRQ(ierr); 1155708d4cb3SBarry Smith PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 1156708d4cb3SBarry Smith PetscStackCallHDF5(H5Aread,(attribute, dtype, *(char**)value)); 1157708d4cb3SBarry Smith } else { 115870efba33SBarry Smith PetscStackCallHDF5(H5Aread,(attribute, dtype, value)); 1159708d4cb3SBarry Smith } 1160729ab6d9SBarry Smith PetscStackCallHDF5(H5Aclose,(attribute)); 1161e90c5075SMarek Pecha /* H5Oclose can be used to close groups, datasets, or committed datatypes */ 1162e90c5075SMarek Pecha PetscStackCallHDF5(H5Oclose,(obj)); 11634302210dSVaclav Hapla ierr = PetscFree(parentAbsPath);CHKERRQ(ierr); 1164ad0c61feSMatthew G. Knepley PetscFunctionReturn(0); 1165ad0c61feSMatthew G. Knepley } 1166ad0c61feSMatthew G. Knepley 1167577e0e04SVaclav Hapla /*@C 1168577e0e04SVaclav Hapla PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name 1169577e0e04SVaclav Hapla 1170343a20b2SVaclav Hapla Collective 1171343a20b2SVaclav Hapla 1172577e0e04SVaclav Hapla Input Parameters: 1173577e0e04SVaclav Hapla + viewer - The HDF5 viewer 1174577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1175577e0e04SVaclav Hapla . name - The attribute name 1176577e0e04SVaclav Hapla - datatype - The attribute type 1177577e0e04SVaclav Hapla 1178577e0e04SVaclav Hapla Output Parameter: 1179577e0e04SVaclav Hapla . value - The attribute value 1180577e0e04SVaclav Hapla 1181577e0e04SVaclav Hapla Notes: 1182577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1183577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 1184577e0e04SVaclav Hapla 1185577e0e04SVaclav Hapla Level: advanced 1186577e0e04SVaclav Hapla 1187577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute() PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1188577e0e04SVaclav Hapla @*/ 1189a2d6be1bSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *defaultValue, void *value) 1190577e0e04SVaclav Hapla { 1191577e0e04SVaclav Hapla PetscErrorCode ierr; 1192577e0e04SVaclav Hapla 1193577e0e04SVaclav Hapla PetscFunctionBegin; 1194577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1195577e0e04SVaclav Hapla PetscValidHeader(obj,2); 1196577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 1197064a246eSJacob Faibussowitsch PetscValidPointer(value, 6); 1198577e0e04SVaclav Hapla ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 1199a2d6be1bSVaclav Hapla ierr = PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, defaultValue, value);CHKERRQ(ierr); 1200577e0e04SVaclav Hapla PetscFunctionReturn(0); 1201577e0e04SVaclav Hapla } 1202577e0e04SVaclav Hapla 1203a75c3fd4SVaclav Hapla PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_) 1204a75c3fd4SVaclav Hapla { 1205a75c3fd4SVaclav Hapla htri_t exists; 1206a75c3fd4SVaclav Hapla hid_t group; 1207a75c3fd4SVaclav Hapla 1208a75c3fd4SVaclav Hapla PetscFunctionBegin; 1209a75c3fd4SVaclav Hapla PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT)); 1210a75c3fd4SVaclav Hapla if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT)); 1211a75c3fd4SVaclav Hapla if (!exists && createGroup) { 1212a75c3fd4SVaclav Hapla PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT)); 1213a75c3fd4SVaclav Hapla PetscStackCallHDF5(H5Gclose,(group)); 1214a75c3fd4SVaclav Hapla exists = PETSC_TRUE; 1215a75c3fd4SVaclav Hapla } 1216a75c3fd4SVaclav Hapla *exists_ = (PetscBool) exists; 1217a75c3fd4SVaclav Hapla PetscFunctionReturn(0); 1218a75c3fd4SVaclav Hapla } 1219a75c3fd4SVaclav Hapla 1220bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype) 12215cdeb108SMatthew G. Knepley { 122290e03537SVaclav Hapla const char rootGroupName[] = "/"; 12235cdeb108SMatthew G. Knepley hid_t h5; 1224e5bf9ebcSVaclav Hapla PetscBool exists=PETSC_FALSE; 122515b861d2SVaclav Hapla PetscInt i; 122615b861d2SVaclav Hapla int n; 122785688ae2SVaclav Hapla char **hierarchy; 122885688ae2SVaclav Hapla char buf[PETSC_MAX_PATH_LEN]=""; 12295cdeb108SMatthew G. Knepley PetscErrorCode ierr; 12305cdeb108SMatthew G. Knepley 12315cdeb108SMatthew G. Knepley PetscFunctionBegin; 12325cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 123390e03537SVaclav Hapla if (name) PetscValidCharPointer(name, 2); 123490e03537SVaclav Hapla else name = rootGroupName; 1235ccd66a83SVaclav Hapla if (has) { 1236064a246eSJacob Faibussowitsch PetscValidBoolPointer(has, 4); 12375cdeb108SMatthew G. Knepley *has = PETSC_FALSE; 1238ccd66a83SVaclav Hapla } 1239ccd66a83SVaclav Hapla if (otype) { 1240064a246eSJacob Faibussowitsch PetscValidIntPointer(otype, 5); 124156cc0592SVaclav Hapla *otype = H5O_TYPE_UNKNOWN; 1242ccd66a83SVaclav Hapla } 1243c0cd0301SJed Brown ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 124485688ae2SVaclav Hapla 124585688ae2SVaclav Hapla /* 124685688ae2SVaclav Hapla Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing. 124785688ae2SVaclav Hapla Hence, each of them needs to be tested separately: 124885688ae2SVaclav Hapla 1) whether it's a valid link 124985688ae2SVaclav Hapla 2) whether this link resolves to an object 125085688ae2SVaclav Hapla See H5Oexists_by_name() documentation. 125185688ae2SVaclav Hapla */ 125285688ae2SVaclav Hapla ierr = PetscStrToArray(name,'/',&n,&hierarchy);CHKERRQ(ierr); 125385688ae2SVaclav Hapla if (!n) { 125485688ae2SVaclav Hapla /* Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */ 1255ccd66a83SVaclav Hapla if (has) *has = PETSC_TRUE; 1256ccd66a83SVaclav Hapla if (otype) *otype = H5O_TYPE_GROUP; 125785688ae2SVaclav Hapla ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 125885688ae2SVaclav Hapla PetscFunctionReturn(0); 125985688ae2SVaclav Hapla } 126085688ae2SVaclav Hapla for (i=0; i<n; i++) { 126185688ae2SVaclav Hapla ierr = PetscStrcat(buf,"/");CHKERRQ(ierr); 126285688ae2SVaclav Hapla ierr = PetscStrcat(buf,hierarchy[i]);CHKERRQ(ierr); 1263a75c3fd4SVaclav Hapla ierr = PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);CHKERRQ(ierr); 1264a75c3fd4SVaclav Hapla if (!exists) break; 126585688ae2SVaclav Hapla } 126685688ae2SVaclav Hapla ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 126785688ae2SVaclav Hapla 126885688ae2SVaclav Hapla /* If the object exists, get its type */ 1269ccd66a83SVaclav Hapla if (exists && otype) { 12705cdeb108SMatthew G. Knepley H5O_info_t info; 12715cdeb108SMatthew G. Knepley 127276276748SVaclav Hapla /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */ 127304633f7fSVaclav Hapla PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT)); 127456cc0592SVaclav Hapla *otype = info.type; 12755cdeb108SMatthew G. Knepley } 1276ccd66a83SVaclav Hapla if (has) *has = exists; 12775cdeb108SMatthew G. Knepley PetscFunctionReturn(0); 12785cdeb108SMatthew G. Knepley } 12795cdeb108SMatthew G. Knepley 12804302210dSVaclav Hapla /*@C 12810a9f38efSVaclav Hapla PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file 12820a9f38efSVaclav Hapla 1283343a20b2SVaclav Hapla Collective 1284343a20b2SVaclav Hapla 12850a9f38efSVaclav Hapla Input Parameters: 12864302210dSVaclav Hapla + viewer - The HDF5 viewer 12874302210dSVaclav Hapla - path - The group path 12880a9f38efSVaclav Hapla 12890a9f38efSVaclav Hapla Output Parameter: 12900a9f38efSVaclav Hapla . has - Flag for group existence 12910a9f38efSVaclav Hapla 12920a9f38efSVaclav Hapla Level: advanced 12930a9f38efSVaclav Hapla 12944302210dSVaclav Hapla Notes: 1295785c443eSVaclav 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. 1296785c443eSVaclav Hapla So NULL or empty path means the current pushed group. 1297343a20b2SVaclav Hapla If path exists but is not a group, PETSC_FALSE is returned. 12984302210dSVaclav Hapla 129989e0ef10SVaclav Hapla .seealso: PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasDataset(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5GetGroup(), PetscViewerHDF5OpenGroup() 13000a9f38efSVaclav Hapla @*/ 13014302210dSVaclav Hapla PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, const char path[], PetscBool *has) 13020a9f38efSVaclav Hapla { 13030a9f38efSVaclav Hapla H5O_type_t type; 13044302210dSVaclav Hapla char *abspath; 13050a9f38efSVaclav Hapla PetscErrorCode ierr; 13060a9f38efSVaclav Hapla 13070a9f38efSVaclav Hapla PetscFunctionBegin; 13080a9f38efSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 13094302210dSVaclav Hapla if (path) PetscValidCharPointer(path,2); 13104302210dSVaclav Hapla PetscValidBoolPointer(has,3); 13114302210dSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, path, &abspath);CHKERRQ(ierr); 13124302210dSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type);CHKERRQ(ierr); 13134302210dSVaclav Hapla *has = (PetscBool)(type == H5O_TYPE_GROUP); 13144302210dSVaclav Hapla ierr = PetscFree(abspath);CHKERRQ(ierr); 13150a9f38efSVaclav Hapla PetscFunctionReturn(0); 13160a9f38efSVaclav Hapla } 13170a9f38efSVaclav Hapla 131889e0ef10SVaclav Hapla /*@C 131989e0ef10SVaclav Hapla PetscViewerHDF5HasDataset - Check whether a given dataset exists in the HDF5 file 132089e0ef10SVaclav Hapla 1321343a20b2SVaclav Hapla Collective 1322343a20b2SVaclav Hapla 132389e0ef10SVaclav Hapla Input Parameters: 132489e0ef10SVaclav Hapla + viewer - The HDF5 viewer 132589e0ef10SVaclav Hapla - path - The dataset path 132689e0ef10SVaclav Hapla 132789e0ef10SVaclav Hapla Output Parameter: 132889e0ef10SVaclav Hapla . has - Flag whether dataset exists 132989e0ef10SVaclav Hapla 133089e0ef10SVaclav Hapla Level: advanced 133189e0ef10SVaclav Hapla 133289e0ef10SVaclav Hapla Notes: 1333343a20b2SVaclav 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. 133489e0ef10SVaclav Hapla If path is NULL or empty, has is set to PETSC_FALSE. 1335343a20b2SVaclav Hapla If path exists but is not a dataset, has is set to PETSC_FALSE as well. 133689e0ef10SVaclav Hapla 133789e0ef10SVaclav Hapla .seealso: PetscViewerHDF5HasObject(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasGroup(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5GetGroup() 133889e0ef10SVaclav Hapla @*/ 133989e0ef10SVaclav Hapla PetscErrorCode PetscViewerHDF5HasDataset(PetscViewer viewer, const char path[], PetscBool *has) 134089e0ef10SVaclav Hapla { 134189e0ef10SVaclav Hapla H5O_type_t type; 134289e0ef10SVaclav Hapla char *abspath; 134389e0ef10SVaclav Hapla PetscErrorCode ierr; 134489e0ef10SVaclav Hapla 134589e0ef10SVaclav Hapla PetscFunctionBegin; 134689e0ef10SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 134789e0ef10SVaclav Hapla if (path) PetscValidCharPointer(path,2); 134889e0ef10SVaclav Hapla PetscValidBoolPointer(has,3); 134989e0ef10SVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, path, &abspath);CHKERRQ(ierr); 135089e0ef10SVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type);CHKERRQ(ierr); 135189e0ef10SVaclav Hapla *has = (PetscBool)(type == H5O_TYPE_DATASET); 135289e0ef10SVaclav Hapla ierr = PetscFree(abspath);CHKERRQ(ierr); 135389e0ef10SVaclav Hapla PetscFunctionReturn(0); 135489e0ef10SVaclav Hapla } 135589e0ef10SVaclav Hapla 13560a9f38efSVaclav Hapla /*@ 1357e3f143f7SVaclav Hapla PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group 1358ecce1506SVaclav Hapla 1359343a20b2SVaclav Hapla Collective 1360343a20b2SVaclav Hapla 1361ecce1506SVaclav Hapla Input Parameters: 1362ecce1506SVaclav Hapla + viewer - The HDF5 viewer 1363ecce1506SVaclav Hapla - obj - The named object 1364ecce1506SVaclav Hapla 1365ecce1506SVaclav Hapla Output Parameter: 136689e0ef10SVaclav Hapla . has - Flag for dataset existence 1367ecce1506SVaclav Hapla 1368e3f143f7SVaclav Hapla Notes: 136989e0ef10SVaclav Hapla If the object is unnamed, an error occurs. 1370343a20b2SVaclav Hapla If the path current_group/object_name exists but is not a dataset, has is set to PETSC_FALSE as well. 1371e3f143f7SVaclav Hapla 1372ecce1506SVaclav Hapla Level: advanced 1373ecce1506SVaclav Hapla 137489e0ef10SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasDataset(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5GetGroup() 1375ecce1506SVaclav Hapla @*/ 1376ecce1506SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has) 1377ecce1506SVaclav Hapla { 137889e0ef10SVaclav Hapla size_t len; 1379ecce1506SVaclav Hapla PetscErrorCode ierr; 1380ecce1506SVaclav Hapla 1381ecce1506SVaclav Hapla PetscFunctionBegin; 1382c57a1a9aSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1383c57a1a9aSVaclav Hapla PetscValidHeader(obj,2); 13844302210dSVaclav Hapla PetscValidBoolPointer(has,3); 138589e0ef10SVaclav Hapla ierr = PetscStrlen(obj->name, &len);CHKERRQ(ierr); 138689e0ef10SVaclav Hapla if (!len) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named"); 138789e0ef10SVaclav Hapla ierr = PetscViewerHDF5HasDataset(viewer, obj->name, has);CHKERRQ(ierr); 1388ecce1506SVaclav Hapla PetscFunctionReturn(0); 1389ecce1506SVaclav Hapla } 1390ecce1506SVaclav Hapla 1391df863907SAlex Fikl /*@C 1392ece83b6cSVaclav Hapla PetscViewerHDF5HasAttribute - Check whether an attribute exists 1393e7bdbf8eSMatthew G. Knepley 1394343a20b2SVaclav Hapla Collective 1395343a20b2SVaclav Hapla 1396e7bdbf8eSMatthew G. Knepley Input Parameters: 1397e7bdbf8eSMatthew G. Knepley + viewer - The HDF5 viewer 13984302210dSVaclav Hapla . parent - The parent dataset/group name 1399e7bdbf8eSMatthew G. Knepley - name - The attribute name 1400e7bdbf8eSMatthew G. Knepley 1401e7bdbf8eSMatthew G. Knepley Output Parameter: 1402e7bdbf8eSMatthew G. Knepley . has - Flag for attribute existence 1403e7bdbf8eSMatthew G. Knepley 1404e7bdbf8eSMatthew G. Knepley Level: advanced 1405e7bdbf8eSMatthew G. Knepley 14064302210dSVaclav Hapla Notes: 1407343a20b2SVaclav 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. 14084302210dSVaclav Hapla 1409577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1410e7bdbf8eSMatthew G. Knepley @*/ 14114302210dSVaclav Hapla PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 1412e7bdbf8eSMatthew G. Knepley { 14134302210dSVaclav Hapla char *parentAbsPath; 1414e7bdbf8eSMatthew G. Knepley PetscErrorCode ierr; 1415e7bdbf8eSMatthew G. Knepley 1416e7bdbf8eSMatthew G. Knepley PetscFunctionBegin; 14175cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 14184302210dSVaclav Hapla if (parent) PetscValidCharPointer(parent,2); 1419c57a1a9aSVaclav Hapla PetscValidCharPointer(name,3); 14204302210dSVaclav Hapla PetscValidBoolPointer(has,4); 14214302210dSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath);CHKERRQ(ierr); 14224302210dSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, has, NULL);CHKERRQ(ierr); 14234302210dSVaclav Hapla if (*has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, has);CHKERRQ(ierr);} 14244302210dSVaclav Hapla ierr = PetscFree(parentAbsPath);CHKERRQ(ierr); 142506db490cSVaclav Hapla PetscFunctionReturn(0); 142606db490cSVaclav Hapla } 142706db490cSVaclav Hapla 1428577e0e04SVaclav Hapla /*@C 1429577e0e04SVaclav Hapla PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name 1430577e0e04SVaclav Hapla 1431343a20b2SVaclav Hapla Collective 1432343a20b2SVaclav Hapla 1433577e0e04SVaclav Hapla Input Parameters: 1434577e0e04SVaclav Hapla + viewer - The HDF5 viewer 1435577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1436577e0e04SVaclav Hapla - name - The attribute name 1437577e0e04SVaclav Hapla 1438577e0e04SVaclav Hapla Output Parameter: 1439577e0e04SVaclav Hapla . has - Flag for attribute existence 1440577e0e04SVaclav Hapla 1441577e0e04SVaclav Hapla Notes: 1442577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1443577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 1444577e0e04SVaclav Hapla 1445577e0e04SVaclav Hapla Level: advanced 1446577e0e04SVaclav Hapla 1447577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1448577e0e04SVaclav Hapla @*/ 1449577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has) 1450577e0e04SVaclav Hapla { 1451577e0e04SVaclav Hapla PetscErrorCode ierr; 1452577e0e04SVaclav Hapla 1453577e0e04SVaclav Hapla PetscFunctionBegin; 1454577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1455577e0e04SVaclav Hapla PetscValidHeader(obj,2); 1456577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 14574302210dSVaclav Hapla PetscValidBoolPointer(has,4); 1458577e0e04SVaclav Hapla ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 1459577e0e04SVaclav Hapla ierr = PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);CHKERRQ(ierr); 1460577e0e04SVaclav Hapla PetscFunctionReturn(0); 1461577e0e04SVaclav Hapla } 1462577e0e04SVaclav Hapla 146306db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 146406db490cSVaclav Hapla { 146506db490cSVaclav Hapla hid_t h5; 146606db490cSVaclav Hapla htri_t hhas; 146706db490cSVaclav Hapla PetscErrorCode ierr; 146806db490cSVaclav Hapla 146906db490cSVaclav Hapla PetscFunctionBegin; 147006db490cSVaclav Hapla ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 14712f936e54SVaclav Hapla PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT)); 14725cdeb108SMatthew G. Knepley *has = hhas ? PETSC_TRUE : PETSC_FALSE; 1473e7bdbf8eSMatthew G. Knepley PetscFunctionReturn(0); 1474e7bdbf8eSMatthew G. Knepley } 1475e7bdbf8eSMatthew G. Knepley 1476a75e6a4aSMatthew G. Knepley /* 1477a75e6a4aSMatthew G. Knepley The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that 1478a75e6a4aSMatthew G. Knepley is attached to a communicator, in this case the attribute is a PetscViewer. 1479a75e6a4aSMatthew G. Knepley */ 1480d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID; 1481a75e6a4aSMatthew G. Knepley 1482a75e6a4aSMatthew G. Knepley /*@C 1483a75e6a4aSMatthew G. Knepley PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator. 1484a75e6a4aSMatthew G. Knepley 1485d083f849SBarry Smith Collective 1486a75e6a4aSMatthew G. Knepley 1487a75e6a4aSMatthew G. Knepley Input Parameter: 1488a75e6a4aSMatthew G. Knepley . comm - the MPI communicator to share the HDF5 PetscViewer 1489a75e6a4aSMatthew G. Knepley 1490a75e6a4aSMatthew G. Knepley Level: intermediate 1491a75e6a4aSMatthew G. Knepley 1492a75e6a4aSMatthew G. Knepley Options Database Keys: 149310699b91SBarry Smith . -viewer_hdf5_filename <name> - name of the HDF5 file 1494a75e6a4aSMatthew G. Knepley 1495a75e6a4aSMatthew G. Knepley Environmental variables: 149610699b91SBarry Smith . PETSC_VIEWER_HDF5_FILENAME - name of the HDF5 file 1497a75e6a4aSMatthew G. Knepley 1498a75e6a4aSMatthew G. Knepley Notes: 1499a75e6a4aSMatthew G. Knepley Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return 1500a75e6a4aSMatthew G. Knepley an error code. The HDF5 PetscViewer is usually used in the form 1501a75e6a4aSMatthew G. Knepley $ XXXView(XXX object, PETSC_VIEWER_HDF5_(comm)); 1502a75e6a4aSMatthew G. Knepley 1503a75e6a4aSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy() 1504a75e6a4aSMatthew G. Knepley @*/ 1505a75e6a4aSMatthew G. Knepley PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm) 1506a75e6a4aSMatthew G. Knepley { 1507a75e6a4aSMatthew G. Knepley PetscErrorCode ierr; 1508a75e6a4aSMatthew G. Knepley PetscBool flg; 1509a75e6a4aSMatthew G. Knepley PetscViewer viewer; 1510a75e6a4aSMatthew G. Knepley char fname[PETSC_MAX_PATH_LEN]; 1511a75e6a4aSMatthew G. Knepley MPI_Comm ncomm; 1512a75e6a4aSMatthew G. Knepley 1513a75e6a4aSMatthew G. Knepley PetscFunctionBegin; 1514908793a3SLisandro 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);} 1515a75e6a4aSMatthew G. Knepley if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) { 1516908793a3SLisandro Dalcin ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,NULL); 1517908793a3SLisandro Dalcin if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 1518a75e6a4aSMatthew G. Knepley } 151947435625SJed Brown ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg); 1520908793a3SLisandro Dalcin if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 1521a75e6a4aSMatthew G. Knepley if (!flg) { /* PetscViewer not yet created */ 1522a75e6a4aSMatthew G. Knepley ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg); 15232cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1524a75e6a4aSMatthew G. Knepley if (!flg) { 1525a75e6a4aSMatthew G. Knepley ierr = PetscStrcpy(fname,"output.h5"); 15262cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1527a75e6a4aSMatthew G. Knepley } 1528a75e6a4aSMatthew G. Knepley ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer); 15292cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1530a75e6a4aSMatthew G. Knepley ierr = PetscObjectRegisterDestroy((PetscObject)viewer); 15312cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 153247435625SJed Brown ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer); 1533908793a3SLisandro Dalcin if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 1534a75e6a4aSMatthew G. Knepley } 1535a75e6a4aSMatthew G. Knepley ierr = PetscCommDestroy(&ncomm); 15362cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1537a75e6a4aSMatthew G. Knepley PetscFunctionReturn(viewer); 1538a75e6a4aSMatthew G. Knepley } 1539