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 869b1403beSVaclav Hapla static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer) 875c6c1daeSBarry Smith { 885c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 895c6c1daeSBarry Smith PetscErrorCode ierr; 905c6c1daeSBarry Smith 915c6c1daeSBarry Smith PetscFunctionBegin; 929c5072fbSVaclav Hapla PetscStackCallHDF5(H5Pclose,(hdf5->dxpl_id)); 935c6c1daeSBarry Smith ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr); 945c6c1daeSBarry Smith while (hdf5->groups) { 957d964842SVaclav Hapla PetscViewerHDF5GroupList *tmp = hdf5->groups->next; 965c6c1daeSBarry Smith 975c6c1daeSBarry Smith ierr = PetscFree(hdf5->groups->name);CHKERRQ(ierr); 985c6c1daeSBarry Smith ierr = PetscFree(hdf5->groups);CHKERRQ(ierr); 995c6c1daeSBarry Smith hdf5->groups = tmp; 1005c6c1daeSBarry Smith } 1015c6c1daeSBarry Smith ierr = PetscFree(hdf5);CHKERRQ(ierr); 1020b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr); 103d1232d7fSToby Isaac ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr); 1040b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr); 105058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetBaseDimension2_C",NULL);CHKERRQ(ierr); 106058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetSPOutput_C",NULL);CHKERRQ(ierr); 1075c6c1daeSBarry Smith PetscFunctionReturn(0); 1085c6c1daeSBarry Smith } 1095c6c1daeSBarry Smith 1109b1403beSVaclav Hapla static PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type) 1115c6c1daeSBarry Smith { 1125c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 1135c6c1daeSBarry Smith 1145c6c1daeSBarry Smith PetscFunctionBegin; 1155c6c1daeSBarry Smith hdf5->btype = type; 1165c6c1daeSBarry Smith PetscFunctionReturn(0); 1175c6c1daeSBarry Smith } 1185c6c1daeSBarry Smith 1190b62783dSVaclav Hapla static PetscErrorCode PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type) 1200b62783dSVaclav Hapla { 1210b62783dSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 1220b62783dSVaclav Hapla 1230b62783dSVaclav Hapla PetscFunctionBegin; 1240b62783dSVaclav Hapla *type = hdf5->btype; 1250b62783dSVaclav Hapla PetscFunctionReturn(0); 1260b62783dSVaclav Hapla } 1270b62783dSVaclav Hapla 1289b1403beSVaclav Hapla static PetscErrorCode PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg) 12982ea9b62SBarry Smith { 13082ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 13182ea9b62SBarry Smith 13282ea9b62SBarry Smith PetscFunctionBegin; 13382ea9b62SBarry Smith hdf5->basedimension2 = flg; 13482ea9b62SBarry Smith PetscFunctionReturn(0); 13582ea9b62SBarry Smith } 13682ea9b62SBarry Smith 137df863907SAlex Fikl /*@ 13882ea9b62SBarry Smith PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 13982ea9b62SBarry Smith dimension of 2. 14082ea9b62SBarry Smith 14182ea9b62SBarry Smith Logically Collective on PetscViewer 14282ea9b62SBarry Smith 14382ea9b62SBarry Smith Input Parameters: 14482ea9b62SBarry Smith + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 14582ea9b62SBarry 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 14682ea9b62SBarry Smith 14782ea9b62SBarry Smith Options Database: 14882ea9b62SBarry 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 14982ea9b62SBarry Smith 15095452b02SPatrick Sanan Notes: 15195452b02SPatrick 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 15282ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 15382ea9b62SBarry Smith 15482ea9b62SBarry Smith Level: intermediate 15582ea9b62SBarry Smith 15682ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 15782ea9b62SBarry Smith 15882ea9b62SBarry Smith @*/ 15982ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg) 16082ea9b62SBarry Smith { 16182ea9b62SBarry Smith PetscErrorCode ierr; 16282ea9b62SBarry Smith 16382ea9b62SBarry Smith PetscFunctionBegin; 16482ea9b62SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 16582ea9b62SBarry Smith ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 16682ea9b62SBarry Smith PetscFunctionReturn(0); 16782ea9b62SBarry Smith } 16882ea9b62SBarry Smith 169df863907SAlex Fikl /*@ 17082ea9b62SBarry Smith PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 17182ea9b62SBarry Smith dimension of 2. 17282ea9b62SBarry Smith 17382ea9b62SBarry Smith Logically Collective on PetscViewer 17482ea9b62SBarry Smith 17582ea9b62SBarry Smith Input Parameter: 17682ea9b62SBarry Smith . viewer - the PetscViewer, must be of type HDF5 17782ea9b62SBarry Smith 17882ea9b62SBarry Smith Output Parameter: 17982ea9b62SBarry 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 18082ea9b62SBarry Smith 18195452b02SPatrick Sanan Notes: 18295452b02SPatrick 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 18382ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 18482ea9b62SBarry Smith 18582ea9b62SBarry Smith Level: intermediate 18682ea9b62SBarry Smith 18782ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 18882ea9b62SBarry Smith 18982ea9b62SBarry Smith @*/ 19082ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg) 19182ea9b62SBarry Smith { 19282ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 19382ea9b62SBarry Smith 19482ea9b62SBarry Smith PetscFunctionBegin; 19582ea9b62SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 19682ea9b62SBarry Smith *flg = hdf5->basedimension2; 19782ea9b62SBarry Smith PetscFunctionReturn(0); 19882ea9b62SBarry Smith } 19982ea9b62SBarry Smith 2009b1403beSVaclav Hapla static PetscErrorCode PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg) 2019a0502c6SHåkon Strandenes { 2029a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 2039a0502c6SHåkon Strandenes 2049a0502c6SHåkon Strandenes PetscFunctionBegin; 2059a0502c6SHåkon Strandenes hdf5->spoutput = flg; 2069a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2079a0502c6SHåkon Strandenes } 2089a0502c6SHåkon Strandenes 209df863907SAlex Fikl /*@ 2109a0502c6SHåkon Strandenes PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is 2119a0502c6SHåkon Strandenes compiled with double precision PetscReal. 2129a0502c6SHåkon Strandenes 2139a0502c6SHåkon Strandenes Logically Collective on PetscViewer 2149a0502c6SHåkon Strandenes 2159a0502c6SHåkon Strandenes Input Parameters: 2169a0502c6SHåkon Strandenes + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 2179a0502c6SHåkon Strandenes - flg - if PETSC_TRUE the data will be written to disk with single precision 2189a0502c6SHåkon Strandenes 2199a0502c6SHåkon Strandenes Options Database: 2209a0502c6SHåkon Strandenes . -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision 2219a0502c6SHåkon Strandenes 22295452b02SPatrick Sanan Notes: 22395452b02SPatrick Sanan Setting this option does not make any difference if PETSc is compiled with single precision 2249a0502c6SHåkon Strandenes in the first place. It does not affect reading datasets (HDF5 handle this internally). 2259a0502c6SHåkon Strandenes 2269a0502c6SHåkon Strandenes Level: intermediate 2279a0502c6SHåkon Strandenes 2289a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 2299a0502c6SHåkon Strandenes PetscReal 2309a0502c6SHåkon Strandenes 2319a0502c6SHåkon Strandenes @*/ 2329a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg) 2339a0502c6SHåkon Strandenes { 2349a0502c6SHåkon Strandenes PetscErrorCode ierr; 2359a0502c6SHåkon Strandenes 2369a0502c6SHåkon Strandenes PetscFunctionBegin; 2379a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 2389a0502c6SHåkon Strandenes ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 2399a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2409a0502c6SHåkon Strandenes } 2419a0502c6SHåkon Strandenes 242df863907SAlex Fikl /*@ 2439a0502c6SHåkon Strandenes PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is 2449a0502c6SHåkon Strandenes compiled with double precision PetscReal. 2459a0502c6SHåkon Strandenes 2469a0502c6SHåkon Strandenes Logically Collective on PetscViewer 2479a0502c6SHåkon Strandenes 2489a0502c6SHåkon Strandenes Input Parameter: 2499a0502c6SHåkon Strandenes . viewer - the PetscViewer, must be of type HDF5 2509a0502c6SHåkon Strandenes 2519a0502c6SHåkon Strandenes Output Parameter: 2529a0502c6SHåkon Strandenes . flg - if PETSC_TRUE the data will be written to disk with single precision 2539a0502c6SHåkon Strandenes 25495452b02SPatrick Sanan Notes: 25595452b02SPatrick Sanan Setting this option does not make any difference if PETSc is compiled with single precision 2569a0502c6SHåkon Strandenes in the first place. It does not affect reading datasets (HDF5 handle this internally). 2579a0502c6SHåkon Strandenes 2589a0502c6SHåkon Strandenes Level: intermediate 2599a0502c6SHåkon Strandenes 2609a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 2619a0502c6SHåkon Strandenes PetscReal 2629a0502c6SHåkon Strandenes 2639a0502c6SHåkon Strandenes @*/ 2649a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg) 2659a0502c6SHåkon Strandenes { 2669a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 2679a0502c6SHåkon Strandenes 2689a0502c6SHåkon Strandenes PetscFunctionBegin; 2699a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 2709a0502c6SHåkon Strandenes *flg = hdf5->spoutput; 2719a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2729a0502c6SHåkon Strandenes } 2739a0502c6SHåkon Strandenes 274ee8b9732SVaclav Hapla static PetscErrorCode PetscViewerHDF5SetCollective_HDF5(PetscViewer viewer, PetscBool flg) 275ee8b9732SVaclav Hapla { 276ee8b9732SVaclav Hapla PetscFunctionBegin; 277ee8b9732SVaclav Hapla /* H5FD_MPIO_COLLECTIVE is wrong in hdf5 1.10.2, and is the same as H5FD_MPIO_INDEPENDENT in earlier versions 278ee8b9732SVaclav Hapla - see e.g. https://gitlab.cosma.dur.ac.uk/swift/swiftsim/issues/431 */ 279c796909eSBarry Smith #if H5_VERSION_GE(1,10,3) && defined(H5_HAVE_PARALLEL) 280ee8b9732SVaclav Hapla { 281ee8b9732SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 282ee8b9732SVaclav Hapla PetscStackCallHDF5(H5Pset_dxpl_mpio,(hdf5->dxpl_id, flg ? H5FD_MPIO_COLLECTIVE : H5FD_MPIO_INDEPENDENT)); 283ee8b9732SVaclav Hapla } 284ee8b9732SVaclav Hapla #else 285ee8b9732SVaclav Hapla if (flg) { 286ee8b9732SVaclav Hapla PetscErrorCode ierr; 287c796909eSBarry 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); 288ee8b9732SVaclav Hapla } 289ee8b9732SVaclav Hapla #endif 290ee8b9732SVaclav Hapla PetscFunctionReturn(0); 291ee8b9732SVaclav Hapla } 292ee8b9732SVaclav Hapla 293ee8b9732SVaclav Hapla /*@ 294ee8b9732SVaclav Hapla PetscViewerHDF5SetCollective - Use collective MPI-IO transfer mode for HDF5 reads and writes. 295ee8b9732SVaclav Hapla 296ee8b9732SVaclav Hapla Logically Collective; flg must contain common value 297ee8b9732SVaclav Hapla 298ee8b9732SVaclav Hapla Input Parameters: 299ee8b9732SVaclav Hapla + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 300ee8b9732SVaclav Hapla - flg - PETSC_TRUE for collective mode; PETSC_FALSE for independent mode (default) 301ee8b9732SVaclav Hapla 302ee8b9732SVaclav Hapla Options Database: 303ee8b9732SVaclav Hapla . -viewer_hdf5_collective - turns on (true) or off (false) collective transfers 304ee8b9732SVaclav Hapla 305ee8b9732SVaclav Hapla Notes: 306ee8b9732SVaclav Hapla Collective mode gives the MPI-IO layer underneath HDF5 a chance to do some additional collective optimizations and hence can perform better. 30753effed7SBarry 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. 308ee8b9732SVaclav Hapla 309ee8b9732SVaclav Hapla Developer notes: 310ee8b9732SVaclav Hapla In the HDF5 layer, PETSC_TRUE / PETSC_FALSE means H5Pset_dxpl_mpio() is called with H5FD_MPIO_COLLECTIVE / H5FD_MPIO_INDEPENDENT, respectively. 311ee8b9732SVaclav Hapla This in turn means use of MPI_File_{read,write}_all / MPI_File_{read,write} in the MPI-IO layer, respectively. 312ee8b9732SVaclav Hapla See HDF5 documentation and MPI-IO documentation for details. 313ee8b9732SVaclav Hapla 314ee8b9732SVaclav Hapla Level: intermediate 315ee8b9732SVaclav Hapla 316ee8b9732SVaclav Hapla .seealso: PetscViewerHDF5GetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open() 317ee8b9732SVaclav Hapla 318ee8b9732SVaclav Hapla @*/ 319ee8b9732SVaclav Hapla PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer viewer,PetscBool flg) 320ee8b9732SVaclav Hapla { 321ee8b9732SVaclav Hapla PetscErrorCode ierr; 322ee8b9732SVaclav Hapla 323ee8b9732SVaclav Hapla PetscFunctionBegin; 324ee8b9732SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 325ee8b9732SVaclav Hapla PetscValidLogicalCollectiveBool(viewer,flg,2); 326ee8b9732SVaclav Hapla ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetCollective_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 327ee8b9732SVaclav Hapla PetscFunctionReturn(0); 328ee8b9732SVaclav Hapla } 329ee8b9732SVaclav Hapla 330ee8b9732SVaclav Hapla static PetscErrorCode PetscViewerHDF5GetCollective_HDF5(PetscViewer viewer, PetscBool *flg) 331ee8b9732SVaclav Hapla { 332c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL) 333ee8b9732SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 334ee8b9732SVaclav Hapla H5FD_mpio_xfer_t mode; 335c796909eSBarry Smith #endif 336ee8b9732SVaclav Hapla 337ee8b9732SVaclav Hapla PetscFunctionBegin; 338c796909eSBarry Smith #if !defined(H5_HAVE_PARALLEL) 339c796909eSBarry Smith *flg = PETSC_FALSE; 340c796909eSBarry Smith #else 341ee8b9732SVaclav Hapla PetscStackCallHDF5(H5Pget_dxpl_mpio,(hdf5->dxpl_id, &mode)); 342ee8b9732SVaclav Hapla *flg = (mode == H5FD_MPIO_COLLECTIVE) ? PETSC_TRUE : PETSC_FALSE; 343c796909eSBarry Smith #endif 344ee8b9732SVaclav Hapla PetscFunctionReturn(0); 345ee8b9732SVaclav Hapla } 346ee8b9732SVaclav Hapla 347ee8b9732SVaclav Hapla /*@ 348ee8b9732SVaclav Hapla PetscViewerHDF5GetCollective - Return flag whether collective MPI-IO transfer mode is used for HDF5 reads and writes. 349ee8b9732SVaclav Hapla 350ee8b9732SVaclav Hapla Not Collective 351ee8b9732SVaclav Hapla 352ee8b9732SVaclav Hapla Input Parameters: 353ee8b9732SVaclav Hapla . viewer - the HDF5 PetscViewer 354ee8b9732SVaclav Hapla 355ee8b9732SVaclav Hapla Output Parameters: 356ee8b9732SVaclav Hapla . flg - the flag 357ee8b9732SVaclav Hapla 358ee8b9732SVaclav Hapla Level: intermediate 359ee8b9732SVaclav Hapla 360ee8b9732SVaclav Hapla Notes: 361c796909eSBarry 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. 362ee8b9732SVaclav Hapla For more details, see PetscViewerHDF5SetCollective(). 363ee8b9732SVaclav Hapla 364ee8b9732SVaclav Hapla .seealso: PetscViewerHDF5SetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open() 365ee8b9732SVaclav Hapla 366ee8b9732SVaclav Hapla @*/ 367ee8b9732SVaclav Hapla PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer viewer,PetscBool *flg) 368ee8b9732SVaclav Hapla { 369ee8b9732SVaclav Hapla PetscErrorCode ierr; 370ee8b9732SVaclav Hapla 371ee8b9732SVaclav Hapla PetscFunctionBegin; 372ee8b9732SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 373534a8f05SLisandro Dalcin PetscValidBoolPointer(flg,2); 374ee8b9732SVaclav Hapla 375ee8b9732SVaclav Hapla ierr = PetscUseMethod(viewer,"PetscViewerHDF5GetCollective_C",(PetscViewer,PetscBool*),(viewer,flg));CHKERRQ(ierr); 376ee8b9732SVaclav Hapla PetscFunctionReturn(0); 377ee8b9732SVaclav Hapla } 378ee8b9732SVaclav Hapla 3799b1403beSVaclav Hapla static PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[]) 3805c6c1daeSBarry Smith { 3815c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 3825c6c1daeSBarry Smith hid_t plist_id; 3835c6c1daeSBarry Smith PetscErrorCode ierr; 3845c6c1daeSBarry Smith 3855c6c1daeSBarry Smith PetscFunctionBegin; 386f683cc58SToby Isaac if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 387f683cc58SToby Isaac if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);} 3885c6c1daeSBarry Smith ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr); 3895c6c1daeSBarry Smith /* Set up file access property list with parallel I/O access */ 390729ab6d9SBarry Smith PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS)); 391c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL) 392d28bfa55SVaclav Hapla PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), MPI_INFO_NULL)); 393c796909eSBarry Smith #endif 3945c6c1daeSBarry Smith /* Create or open the file collectively */ 3955c6c1daeSBarry Smith switch (hdf5->btype) { 3965c6c1daeSBarry Smith case FILE_MODE_READ: 397729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id)); 3985c6c1daeSBarry Smith break; 3995c6c1daeSBarry Smith case FILE_MODE_APPEND: 4007e4fd573SVaclav Hapla case FILE_MODE_UPDATE: 4017e4fd573SVaclav Hapla { 4027e4fd573SVaclav Hapla PetscBool flg; 4037e4fd573SVaclav Hapla ierr = PetscTestFile(hdf5->filename, 'r', &flg);CHKERRQ(ierr); 4047e4fd573SVaclav Hapla if (flg) PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id)); 4057e4fd573SVaclav Hapla else PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_EXCL, H5P_DEFAULT, plist_id)); 4065c6c1daeSBarry Smith break; 4077e4fd573SVaclav Hapla } 4085c6c1daeSBarry Smith case FILE_MODE_WRITE: 409729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id)); 4105c6c1daeSBarry Smith break; 4117e4fd573SVaclav Hapla case FILE_MODE_UNDEFINED: 4127e4fd573SVaclav Hapla SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 4135c6c1daeSBarry Smith default: 4147e4fd573SVaclav Hapla SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP, "Unsupported file mode %s",PetscFileModes[hdf5->btype]); 4155c6c1daeSBarry Smith } 4165c6c1daeSBarry Smith if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name); 417729ab6d9SBarry Smith PetscStackCallHDF5(H5Pclose,(plist_id)); 4185c6c1daeSBarry Smith PetscFunctionReturn(0); 4195c6c1daeSBarry Smith } 4205c6c1daeSBarry Smith 421d1232d7fSToby Isaac static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name) 422d1232d7fSToby Isaac { 423d1232d7fSToby Isaac PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data; 424d1232d7fSToby Isaac 425d1232d7fSToby Isaac PetscFunctionBegin; 426d1232d7fSToby Isaac *name = vhdf5->filename; 427d1232d7fSToby Isaac PetscFunctionReturn(0); 428d1232d7fSToby Isaac } 429d1232d7fSToby Isaac 430b723ab35SVaclav Hapla static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer) 431b723ab35SVaclav Hapla { 4325dc64a97SVaclav Hapla /* 433b723ab35SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 434b723ab35SVaclav Hapla PetscErrorCode ierr; 4355dc64a97SVaclav Hapla */ 436b723ab35SVaclav Hapla 437b723ab35SVaclav Hapla PetscFunctionBegin; 438b723ab35SVaclav Hapla PetscFunctionReturn(0); 439b723ab35SVaclav Hapla } 440b723ab35SVaclav Hapla 4418556b5ebSBarry Smith /*MC 4428556b5ebSBarry Smith PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file 4438556b5ebSBarry Smith 4448556b5ebSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET, 4458556b5ebSBarry Smith PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING, 4468556b5ebSBarry Smith PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB, 4478556b5ebSBarry Smith PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType() 4488556b5ebSBarry Smith 4491b266c99SBarry Smith Level: beginner 4508556b5ebSBarry Smith M*/ 451d1232d7fSToby Isaac 4528cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v) 4535c6c1daeSBarry Smith { 4545c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5; 4555c6c1daeSBarry Smith PetscErrorCode ierr; 4565c6c1daeSBarry Smith 4575c6c1daeSBarry Smith PetscFunctionBegin; 45899335c34SVaclav Hapla #if !defined(H5_HAVE_PARALLEL) 45999335c34SVaclav Hapla { 46099335c34SVaclav Hapla PetscMPIInt size; 46199335c34SVaclav Hapla ierr = MPI_Comm_size(PetscObjectComm((PetscObject)v), &size);CHKERRMPI(ierr); 46299335c34SVaclav 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)"); 46399335c34SVaclav Hapla } 46499335c34SVaclav Hapla #endif 46599335c34SVaclav Hapla 466b00a9115SJed Brown ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr); 4675c6c1daeSBarry Smith 4685c6c1daeSBarry Smith v->data = (void*) hdf5; 4695c6c1daeSBarry Smith v->ops->destroy = PetscViewerDestroy_HDF5; 47082ea9b62SBarry Smith v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5; 471b723ab35SVaclav Hapla v->ops->setup = PetscViewerSetUp_HDF5; 4721b793a25SVaclav Hapla v->ops->view = PetscViewerView_HDF5; 473908793a3SLisandro Dalcin v->ops->flush = NULL; 4747e4fd573SVaclav Hapla hdf5->btype = FILE_MODE_UNDEFINED; 475908793a3SLisandro Dalcin hdf5->filename = NULL; 4765c6c1daeSBarry Smith hdf5->timestep = -1; 4770298fd71SBarry Smith hdf5->groups = NULL; 4785c6c1daeSBarry Smith 4799c5072fbSVaclav Hapla PetscStackCallHDF5Return(hdf5->dxpl_id,H5Pcreate,(H5P_DATASET_XFER)); 4809c5072fbSVaclav Hapla 4810b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr); 482d1232d7fSToby Isaac ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr); 4830b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr); 4840b62783dSVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_HDF5);CHKERRQ(ierr); 48582ea9b62SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr); 4869a0502c6SHåkon Strandenes ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr); 487ee8b9732SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetCollective_C",PetscViewerHDF5SetCollective_HDF5);CHKERRQ(ierr); 488ee8b9732SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetCollective_C",PetscViewerHDF5GetCollective_HDF5);CHKERRQ(ierr); 4895c6c1daeSBarry Smith PetscFunctionReturn(0); 4905c6c1daeSBarry Smith } 4915c6c1daeSBarry Smith 4925c6c1daeSBarry Smith /*@C 4935c6c1daeSBarry Smith PetscViewerHDF5Open - Opens a file for HDF5 input/output. 4945c6c1daeSBarry Smith 495d083f849SBarry Smith Collective 4965c6c1daeSBarry Smith 4975c6c1daeSBarry Smith Input Parameters: 4985c6c1daeSBarry Smith + comm - MPI communicator 4995c6c1daeSBarry Smith . name - name of file 5005c6c1daeSBarry Smith - type - type of file 5015c6c1daeSBarry Smith 5025c6c1daeSBarry Smith Output Parameter: 5035c6c1daeSBarry Smith . hdf5v - PetscViewer for HDF5 input/output to use with the specified file 5045c6c1daeSBarry Smith 50582ea9b62SBarry Smith Options Database: 506a2b725a8SWilliam 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 507a2b725a8SWilliam Gropp - -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal 50882ea9b62SBarry Smith 5095c6c1daeSBarry Smith Level: beginner 5105c6c1daeSBarry Smith 5117e4fd573SVaclav Hapla Notes: 5127e4fd573SVaclav Hapla Reading is always available, regardless of the mode. Available modes are 5137e4fd573SVaclav Hapla + FILE_MODE_READ - open existing HDF5 file for read only access, fail if file does not exist [H5Fopen() with H5F_ACC_RDONLY] 5147e4fd573SVaclav Hapla . FILE_MODE_WRITE - if file exists, fully overwrite it, else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_TRUNC] 5157e4fd573SVaclav 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] 5167e4fd573SVaclav Hapla - FILE_MODE_UPDATE - same as FILE_MODE_APPEND 5177e4fd573SVaclav Hapla 5187e4fd573SVaclav 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. 5197e4fd573SVaclav Hapla 5205c6c1daeSBarry Smith This PetscViewer should be destroyed with PetscViewerDestroy(). 5215c6c1daeSBarry Smith 5226a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(), 5239a0502c6SHåkon Strandenes PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(), 524a56f64adSBarry Smith MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName() 5255c6c1daeSBarry Smith @*/ 5265c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v) 5275c6c1daeSBarry Smith { 5285c6c1daeSBarry Smith PetscErrorCode ierr; 5295c6c1daeSBarry Smith 5305c6c1daeSBarry Smith PetscFunctionBegin; 5315c6c1daeSBarry Smith ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr); 5325c6c1daeSBarry Smith ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr); 5335c6c1daeSBarry Smith ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr); 5345c6c1daeSBarry Smith ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr); 535b663d2d2SVaclav Hapla ierr = PetscViewerSetFromOptions(*hdf5v);CHKERRQ(ierr); 5365c6c1daeSBarry Smith PetscFunctionReturn(0); 5375c6c1daeSBarry Smith } 5385c6c1daeSBarry Smith 5395c6c1daeSBarry Smith /*@C 5405c6c1daeSBarry Smith PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls 5415c6c1daeSBarry Smith 5425c6c1daeSBarry Smith Not collective 5435c6c1daeSBarry Smith 5445c6c1daeSBarry Smith Input Parameter: 5455c6c1daeSBarry Smith . viewer - the PetscViewer 5465c6c1daeSBarry Smith 5475c6c1daeSBarry Smith Output Parameter: 5485c6c1daeSBarry Smith . file_id - The file id 5495c6c1daeSBarry Smith 5505c6c1daeSBarry Smith Level: intermediate 5515c6c1daeSBarry Smith 5525c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open() 5535c6c1daeSBarry Smith @*/ 5545c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id) 5555c6c1daeSBarry Smith { 5565c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 5575c6c1daeSBarry Smith 5585c6c1daeSBarry Smith PetscFunctionBegin; 5595c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 5605c6c1daeSBarry Smith if (file_id) *file_id = hdf5->file_id; 5615c6c1daeSBarry Smith PetscFunctionReturn(0); 5625c6c1daeSBarry Smith } 5635c6c1daeSBarry Smith 5645c6c1daeSBarry Smith /*@C 5655c6c1daeSBarry Smith PetscViewerHDF5PushGroup - Set the current HDF5 group for output 5665c6c1daeSBarry Smith 5675c6c1daeSBarry Smith Not collective 5685c6c1daeSBarry Smith 5695c6c1daeSBarry Smith Input Parameters: 5705c6c1daeSBarry Smith + viewer - the PetscViewer 5715c6c1daeSBarry Smith - name - The group name 5725c6c1daeSBarry Smith 5735c6c1daeSBarry Smith Level: intermediate 5745c6c1daeSBarry Smith 5752e361344SVaclav Hapla Notes: 5762e361344SVaclav Hapla This is designed to mnemonically resemble the Unix cd command. 5772e361344SVaclav 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. 5782e361344SVaclav Hapla . NULL, empty string, or any sequence of all slashes (e.g. "///") is interpreted as the root group "/". 5792e361344SVaclav Hapla - "." means the current group is pushed again. 5802e361344SVaclav Hapla 5812e361344SVaclav Hapla Example: 5822e361344SVaclav Hapla Suppose the current group is "/a". 5832e361344SVaclav Hapla + If name is NULL, empty string, or a sequence of all slashes (e.g. "///"), then the new group will be "/". 5842e361344SVaclav Hapla . If name is ".", then the new group will be "/a". 5852e361344SVaclav Hapla . If name is "b", then the new group will be "/a/b". 5862e361344SVaclav Hapla - If name is "/b", then the new group will be "/b". 5872e361344SVaclav Hapla 5882e361344SVaclav Hapla Developer Notes: 5892e361344SVaclav Hapla The root group "/" is internally stored as NULL. 590820fc2d1SVaclav Hapla 591874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 5925c6c1daeSBarry Smith @*/ 593be7aa430SVaclav Hapla PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char name[]) 5945c6c1daeSBarry Smith { 5955c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 5967d964842SVaclav Hapla PetscViewerHDF5GroupList *groupNode; 5972e361344SVaclav Hapla size_t i,len; 5982e361344SVaclav Hapla char buf[PETSC_MAX_PATH_LEN]; 5992e361344SVaclav Hapla const char *gname; 6005c6c1daeSBarry Smith PetscErrorCode ierr; 6015c6c1daeSBarry Smith 6025c6c1daeSBarry Smith PetscFunctionBegin; 6035c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 604820fc2d1SVaclav Hapla if (name) PetscValidCharPointer(name,2); 605820fc2d1SVaclav Hapla ierr = PetscStrlen(name, &len);CHKERRQ(ierr); 6062e361344SVaclav Hapla gname = NULL; 6072e361344SVaclav Hapla if (len) { 6082e361344SVaclav Hapla if (len == 1 && name[0] == '.') { 6092e361344SVaclav Hapla /* use current name */ 6102e361344SVaclav Hapla gname = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : NULL; 6112e361344SVaclav Hapla } else if (name[0] == '/') { 6122e361344SVaclav Hapla /* absolute */ 6132e361344SVaclav Hapla for (i=1; i<len; i++) { 6142e361344SVaclav Hapla if (name[i] != '/') { 6152e361344SVaclav Hapla gname = name; 6162e361344SVaclav Hapla break; 6172e361344SVaclav Hapla } 6182e361344SVaclav Hapla } 6192e361344SVaclav Hapla } else { 6202e361344SVaclav Hapla /* relative */ 6212e361344SVaclav Hapla const char *parent = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : ""; 6222e361344SVaclav Hapla ierr = PetscSNPrintf(buf, sizeof(buf), "%s/%s", parent, name);CHKERRQ(ierr); 6232e361344SVaclav Hapla gname = buf; 6242e361344SVaclav Hapla } 6252e361344SVaclav Hapla } 62695dccacaSBarry Smith ierr = PetscNew(&groupNode);CHKERRQ(ierr); 6272e361344SVaclav Hapla ierr = PetscStrallocpy(gname, (char**) &groupNode->name);CHKERRQ(ierr); 6285c6c1daeSBarry Smith groupNode->next = hdf5->groups; 6295c6c1daeSBarry Smith hdf5->groups = groupNode; 6305c6c1daeSBarry Smith PetscFunctionReturn(0); 6315c6c1daeSBarry Smith } 6325c6c1daeSBarry Smith 6333ef9c667SSatish Balay /*@ 6345c6c1daeSBarry Smith PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value 6355c6c1daeSBarry Smith 6365c6c1daeSBarry Smith Not collective 6375c6c1daeSBarry Smith 6385c6c1daeSBarry Smith Input Parameter: 6395c6c1daeSBarry Smith . viewer - the PetscViewer 6405c6c1daeSBarry Smith 6415c6c1daeSBarry Smith Level: intermediate 6425c6c1daeSBarry Smith 643874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 6445c6c1daeSBarry Smith @*/ 6455c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer) 6465c6c1daeSBarry Smith { 6475c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 6487d964842SVaclav Hapla PetscViewerHDF5GroupList *groupNode; 6495c6c1daeSBarry Smith PetscErrorCode ierr; 6505c6c1daeSBarry Smith 6515c6c1daeSBarry Smith PetscFunctionBegin; 6525c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 65382f516ccSBarry Smith if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop"); 6545c6c1daeSBarry Smith groupNode = hdf5->groups; 6555c6c1daeSBarry Smith hdf5->groups = hdf5->groups->next; 6565c6c1daeSBarry Smith ierr = PetscFree(groupNode->name);CHKERRQ(ierr); 6575c6c1daeSBarry Smith ierr = PetscFree(groupNode);CHKERRQ(ierr); 6585c6c1daeSBarry Smith PetscFunctionReturn(0); 6595c6c1daeSBarry Smith } 6605c6c1daeSBarry Smith 6615c6c1daeSBarry Smith /*@C 662874270d9SVaclav Hapla PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup(). 663874270d9SVaclav Hapla If none has been assigned, returns NULL. 6645c6c1daeSBarry Smith 6655c6c1daeSBarry Smith Not collective 6665c6c1daeSBarry Smith 6675c6c1daeSBarry Smith Input Parameter: 6685c6c1daeSBarry Smith . viewer - the PetscViewer 6695c6c1daeSBarry Smith 6705c6c1daeSBarry Smith Output Parameter: 6715c6c1daeSBarry Smith . name - The group name 6725c6c1daeSBarry Smith 6735c6c1daeSBarry Smith Level: intermediate 6745c6c1daeSBarry Smith 675874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup() 6765c6c1daeSBarry Smith @*/ 677be7aa430SVaclav Hapla PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char *name[]) 6785c6c1daeSBarry Smith { 6795c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 6805c6c1daeSBarry Smith 6815c6c1daeSBarry Smith PetscFunctionBegin; 6825c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 683c959eef4SJed Brown PetscValidPointer(name,2); 684a297a907SKarl Rupp if (hdf5->groups) *name = hdf5->groups->name; 6850298fd71SBarry Smith else *name = NULL; 6865c6c1daeSBarry Smith PetscFunctionReturn(0); 6875c6c1daeSBarry Smith } 6885c6c1daeSBarry Smith 6898c557b5aSVaclav Hapla /*@ 690874270d9SVaclav Hapla PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(), 691874270d9SVaclav Hapla and return this group's ID and file ID. 692874270d9SVaclav Hapla If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID. 693874270d9SVaclav Hapla 694874270d9SVaclav Hapla Not collective 695874270d9SVaclav Hapla 696874270d9SVaclav Hapla Input Parameter: 697874270d9SVaclav Hapla . viewer - the PetscViewer 698874270d9SVaclav Hapla 699*d8d19677SJose E. Roman Output Parameters: 700874270d9SVaclav Hapla + fileId - The HDF5 file ID 701874270d9SVaclav Hapla - groupId - The HDF5 group ID 702874270d9SVaclav Hapla 703e74616beSVaclav Hapla Notes: 704e74616beSVaclav Hapla If the viewer is writable, the group is created if it doesn't exist yet. 705e74616beSVaclav Hapla 706874270d9SVaclav Hapla Level: intermediate 707874270d9SVaclav Hapla 708874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 709874270d9SVaclav Hapla @*/ 71054dbf706SVaclav Hapla PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId) 71154dbf706SVaclav Hapla { 71290e03537SVaclav Hapla hid_t file_id; 71390e03537SVaclav Hapla H5O_type_t type; 71454dbf706SVaclav Hapla const char *groupName = NULL; 715e74616beSVaclav Hapla PetscBool create; 71654dbf706SVaclav Hapla PetscErrorCode ierr; 71754dbf706SVaclav Hapla 71854dbf706SVaclav Hapla PetscFunctionBegin; 719e74616beSVaclav Hapla ierr = PetscViewerWritable(viewer, &create);CHKERRQ(ierr); 72054dbf706SVaclav Hapla ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr); 72154dbf706SVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr); 722e74616beSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, groupName, create, NULL, &type);CHKERRQ(ierr); 72390e03537SVaclav Hapla if (type != H5O_TYPE_GROUP) SETERRQ1(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Path %s resolves to something which is not a group", groupName); 72490e03537SVaclav Hapla PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT)); 72554dbf706SVaclav Hapla *fileId = file_id; 72654dbf706SVaclav Hapla PetscFunctionReturn(0); 72754dbf706SVaclav Hapla } 72854dbf706SVaclav Hapla 7293ef9c667SSatish Balay /*@ 730d7dd068bSVaclav Hapla PetscViewerHDF5PushTimestepping - Activate timestepping mode for subsequent HDF5 reading and writing. 7315c6c1daeSBarry Smith 7325c6c1daeSBarry Smith Not collective 7335c6c1daeSBarry Smith 7345c6c1daeSBarry Smith Input Parameter: 7355c6c1daeSBarry Smith . viewer - the PetscViewer 7365c6c1daeSBarry Smith 7375c6c1daeSBarry Smith Level: intermediate 7385c6c1daeSBarry Smith 739d7dd068bSVaclav Hapla Notes: 740d7dd068bSVaclav Hapla On first PetscViewerHDF5PushTimestepping(), the initial time step is set to 0. 741d7dd068bSVaclav Hapla Next timesteps can then be set using PetscViewerHDF5IncrementTimestep() or PetscViewerHDF5SetTimestep(). 742d7dd068bSVaclav 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()]. 743d7dd068bSVaclav Hapla Use PetscViewerHDF5PopTimestepping() to deactivate timestepping mode; calling it by the end of the program is NOT mandatory. 744d7dd068bSVaclav Hapla Current timestep is remembered between PetscViewerHDF5PopTimestepping() and the next PetscViewerHDF5PushTimestepping(). 745d7dd068bSVaclav Hapla 746d7dd068bSVaclav Hapla If a dataset was stored with timestepping, it can be loaded only in the timestepping mode again. 747d7dd068bSVaclav Hapla Loading a timestepped dataset with timestepping disabled, or vice-versa results in an error. 748d7dd068bSVaclav Hapla 749d7dd068bSVaclav Hapla Developer notes: 750d7dd068bSVaclav Hapla Timestepped HDF5 dataset has an extra dimension and attribute "timestepping" set to true. 751d7dd068bSVaclav Hapla 752d7dd068bSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PopTimestepping(), PetscViewerHDF5IsTimestepping(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 753d7dd068bSVaclav Hapla @*/ 754d7dd068bSVaclav Hapla PetscErrorCode PetscViewerHDF5PushTimestepping(PetscViewer viewer) 755d7dd068bSVaclav Hapla { 756d7dd068bSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 757d7dd068bSVaclav Hapla 758d7dd068bSVaclav Hapla PetscFunctionBegin; 759d7dd068bSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 760d7dd068bSVaclav Hapla if (hdf5->timestepping) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping is already pushed"); 761d7dd068bSVaclav Hapla hdf5->timestepping = PETSC_TRUE; 762d7dd068bSVaclav Hapla if (hdf5->timestep < 0) hdf5->timestep = 0; 763d7dd068bSVaclav Hapla PetscFunctionReturn(0); 764d7dd068bSVaclav Hapla } 765d7dd068bSVaclav Hapla 766d7dd068bSVaclav Hapla /*@ 767d7dd068bSVaclav Hapla PetscViewerHDF5PopTimestepping - Deactivate timestepping mode for subsequent HDF5 reading and writing. 768d7dd068bSVaclav Hapla 769d7dd068bSVaclav Hapla Not collective 770d7dd068bSVaclav Hapla 771d7dd068bSVaclav Hapla Input Parameter: 772d7dd068bSVaclav Hapla . viewer - the PetscViewer 773d7dd068bSVaclav Hapla 774d7dd068bSVaclav Hapla Level: intermediate 775d7dd068bSVaclav Hapla 776d7dd068bSVaclav Hapla Notes: 777d7dd068bSVaclav Hapla See PetscViewerHDF5PushTimestepping() for details. 778d7dd068bSVaclav Hapla 779d7dd068bSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5IsTimestepping(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 780d7dd068bSVaclav Hapla @*/ 781d7dd068bSVaclav Hapla PetscErrorCode PetscViewerHDF5PopTimestepping(PetscViewer viewer) 782d7dd068bSVaclav Hapla { 783d7dd068bSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 784d7dd068bSVaclav Hapla 785d7dd068bSVaclav Hapla PetscFunctionBegin; 786d7dd068bSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 787d7dd068bSVaclav Hapla if (!hdf5->timestepping) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 788d7dd068bSVaclav Hapla hdf5->timestepping = PETSC_FALSE; 789d7dd068bSVaclav Hapla PetscFunctionReturn(0); 790d7dd068bSVaclav Hapla } 791d7dd068bSVaclav Hapla 792d7dd068bSVaclav Hapla /*@ 793d7dd068bSVaclav Hapla PetscViewerHDF5IsTimestepping - Ask the viewer whether it is in timestepping mode currently. 794d7dd068bSVaclav Hapla 795d7dd068bSVaclav Hapla Not collective 796d7dd068bSVaclav Hapla 797d7dd068bSVaclav Hapla Input Parameter: 798d7dd068bSVaclav Hapla . viewer - the PetscViewer 799d7dd068bSVaclav Hapla 800d7dd068bSVaclav Hapla Output Parameter: 801d7dd068bSVaclav Hapla . flg - is timestepping active? 802d7dd068bSVaclav Hapla 803d7dd068bSVaclav Hapla Level: intermediate 804d7dd068bSVaclav Hapla 805d7dd068bSVaclav Hapla Notes: 806d7dd068bSVaclav Hapla See PetscViewerHDF5PushTimestepping() for details. 807d7dd068bSVaclav Hapla 808d7dd068bSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5PopTimestepping(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 809d7dd068bSVaclav Hapla @*/ 810d7dd068bSVaclav Hapla PetscErrorCode PetscViewerHDF5IsTimestepping(PetscViewer viewer, PetscBool *flg) 811d7dd068bSVaclav Hapla { 812d7dd068bSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 813d7dd068bSVaclav Hapla 814d7dd068bSVaclav Hapla PetscFunctionBegin; 815d7dd068bSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 816d7dd068bSVaclav Hapla *flg = hdf5->timestepping; 817d7dd068bSVaclav Hapla PetscFunctionReturn(0); 818d7dd068bSVaclav Hapla } 819d7dd068bSVaclav Hapla 820d7dd068bSVaclav Hapla /*@ 821d7dd068bSVaclav Hapla PetscViewerHDF5IncrementTimestep - Increments current timestep for the HDF5 output. Fields are stacked in time. 822d7dd068bSVaclav Hapla 823d7dd068bSVaclav Hapla Not collective 824d7dd068bSVaclav Hapla 825d7dd068bSVaclav Hapla Input Parameter: 826d7dd068bSVaclav Hapla . viewer - the PetscViewer 827d7dd068bSVaclav Hapla 828d7dd068bSVaclav Hapla Level: intermediate 829d7dd068bSVaclav Hapla 830d7dd068bSVaclav Hapla Notes: 831d7dd068bSVaclav Hapla This can be called only if the viewer is in timestepping mode. See PetscViewerHDF5PushTimestepping() for details. 832d7dd068bSVaclav Hapla 833d7dd068bSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep() 8345c6c1daeSBarry Smith @*/ 8355c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer) 8365c6c1daeSBarry Smith { 8375c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 8385c6c1daeSBarry Smith 8395c6c1daeSBarry Smith PetscFunctionBegin; 8405c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 841d7dd068bSVaclav Hapla if (!hdf5->timestepping) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 8425c6c1daeSBarry Smith ++hdf5->timestep; 8435c6c1daeSBarry Smith PetscFunctionReturn(0); 8445c6c1daeSBarry Smith } 8455c6c1daeSBarry Smith 8463ef9c667SSatish Balay /*@ 847d7dd068bSVaclav Hapla PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. 8485c6c1daeSBarry Smith 849d7dd068bSVaclav Hapla Logically collective 8505c6c1daeSBarry Smith 8515c6c1daeSBarry Smith Input Parameters: 8525c6c1daeSBarry Smith + viewer - the PetscViewer 853d7dd068bSVaclav Hapla - timestep - The timestep 8545c6c1daeSBarry Smith 8555c6c1daeSBarry Smith Level: intermediate 8565c6c1daeSBarry Smith 857d7dd068bSVaclav Hapla Notes: 858d7dd068bSVaclav Hapla This can be called only if the viewer is in timestepping mode. See PetscViewerHDF5PushTimestepping() for details. 859d7dd068bSVaclav Hapla 860d7dd068bSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 8615c6c1daeSBarry Smith @*/ 8625c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep) 8635c6c1daeSBarry Smith { 8645c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 8655c6c1daeSBarry Smith 8665c6c1daeSBarry Smith PetscFunctionBegin; 8675c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 868d7dd068bSVaclav Hapla PetscValidLogicalCollectiveInt(viewer, timestep, 2); 869d7dd068bSVaclav Hapla if (timestep < 0) SETERRQ1(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestep %D is negative", timestep); 870d7dd068bSVaclav Hapla if (!hdf5->timestepping) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 8715c6c1daeSBarry Smith hdf5->timestep = timestep; 8725c6c1daeSBarry Smith PetscFunctionReturn(0); 8735c6c1daeSBarry Smith } 8745c6c1daeSBarry Smith 8753ef9c667SSatish Balay /*@ 8765c6c1daeSBarry Smith PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time. 8775c6c1daeSBarry Smith 8785c6c1daeSBarry Smith Not collective 8795c6c1daeSBarry Smith 8805c6c1daeSBarry Smith Input Parameter: 8815c6c1daeSBarry Smith . viewer - the PetscViewer 8825c6c1daeSBarry Smith 8835c6c1daeSBarry Smith Output Parameter: 884d7dd068bSVaclav Hapla . timestep - The timestep 8855c6c1daeSBarry Smith 8865c6c1daeSBarry Smith Level: intermediate 8875c6c1daeSBarry Smith 888d7dd068bSVaclav Hapla Notes: 889d7dd068bSVaclav Hapla This can be called only if the viewer is in the timestepping mode. See PetscViewerHDF5PushTimestepping() for details. 890d7dd068bSVaclav Hapla 891d7dd068bSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep() 8925c6c1daeSBarry Smith @*/ 8935c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep) 8945c6c1daeSBarry Smith { 8955c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 8965c6c1daeSBarry Smith 8975c6c1daeSBarry Smith PetscFunctionBegin; 8985c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 899d7dd068bSVaclav Hapla PetscValidIntPointer(timestep,2); 900d7dd068bSVaclav Hapla if (!hdf5->timestepping) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 9015c6c1daeSBarry Smith *timestep = hdf5->timestep; 9025c6c1daeSBarry Smith PetscFunctionReturn(0); 9035c6c1daeSBarry Smith } 9045c6c1daeSBarry Smith 90536bce6e8SMatthew G. Knepley /*@C 90636bce6e8SMatthew G. Knepley PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name. 90736bce6e8SMatthew G. Knepley 90836bce6e8SMatthew G. Knepley Not collective 90936bce6e8SMatthew G. Knepley 91036bce6e8SMatthew G. Knepley Input Parameter: 91136bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 91236bce6e8SMatthew G. Knepley 91336bce6e8SMatthew G. Knepley Output Parameter: 91436bce6e8SMatthew G. Knepley . mtype - the MPI datatype (for example MPI_DOUBLE, ...) 91536bce6e8SMatthew G. Knepley 91636bce6e8SMatthew G. Knepley Level: advanced 91736bce6e8SMatthew G. Knepley 91836bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 91936bce6e8SMatthew G. Knepley @*/ 92036bce6e8SMatthew G. Knepley PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype) 92136bce6e8SMatthew G. Knepley { 92236bce6e8SMatthew G. Knepley PetscFunctionBegin; 92336bce6e8SMatthew G. Knepley if (ptype == PETSC_INT) 92436bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 92536bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_LLONG; 92636bce6e8SMatthew G. Knepley #else 92736bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_INT; 92836bce6e8SMatthew G. Knepley #endif 92936bce6e8SMatthew G. Knepley else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE; 93036bce6e8SMatthew G. Knepley else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG; 93136bce6e8SMatthew G. Knepley else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT; 932c9450970SBarry Smith else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_INT; 933de3ba262SVaclav Hapla else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_INT; 93436bce6e8SMatthew G. Knepley else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT; 93536bce6e8SMatthew G. Knepley else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR; 93636bce6e8SMatthew G. Knepley else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR; 9377e97c476SMatthew G. Knepley else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1); 93836bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype"); 93936bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 94036bce6e8SMatthew G. Knepley } 94136bce6e8SMatthew G. Knepley 94236bce6e8SMatthew G. Knepley /*@C 94336bce6e8SMatthew G. Knepley PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name 94436bce6e8SMatthew G. Knepley 94536bce6e8SMatthew G. Knepley Not collective 94636bce6e8SMatthew G. Knepley 94736bce6e8SMatthew G. Knepley Input Parameter: 94836bce6e8SMatthew G. Knepley . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...) 94936bce6e8SMatthew G. Knepley 95036bce6e8SMatthew G. Knepley Output Parameter: 95136bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 95236bce6e8SMatthew G. Knepley 95336bce6e8SMatthew G. Knepley Level: advanced 95436bce6e8SMatthew G. Knepley 95536bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 95636bce6e8SMatthew G. Knepley @*/ 95736bce6e8SMatthew G. Knepley PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype) 95836bce6e8SMatthew G. Knepley { 95936bce6e8SMatthew G. Knepley PetscFunctionBegin; 96036bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 96136bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG; 96236bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT; 96336bce6e8SMatthew G. Knepley #else 96436bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT; 96536bce6e8SMatthew G. Knepley #endif 96636bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE; 96736bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG; 96836bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT; 96936bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT; 97036bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR; 97136bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR; 9727e97c476SMatthew G. Knepley else if (htype == H5T_C_S1) *ptype = PETSC_STRING; 97336bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype"); 97436bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 97536bce6e8SMatthew G. Knepley } 97636bce6e8SMatthew G. Knepley 977df863907SAlex Fikl /*@C 978b8ee85a1SVaclav Hapla PetscViewerHDF5WriteAttribute - Write an attribute 97936bce6e8SMatthew G. Knepley 980343a20b2SVaclav Hapla Collective 981343a20b2SVaclav Hapla 98236bce6e8SMatthew G. Knepley Input Parameters: 98336bce6e8SMatthew G. Knepley + viewer - The HDF5 viewer 9844302210dSVaclav Hapla . parent - The parent dataset/group name 98536bce6e8SMatthew G. Knepley . name - The attribute name 98636bce6e8SMatthew G. Knepley . datatype - The attribute type 98736bce6e8SMatthew G. Knepley - value - The attribute value 98836bce6e8SMatthew G. Knepley 98936bce6e8SMatthew G. Knepley Level: advanced 99036bce6e8SMatthew G. Knepley 9914302210dSVaclav Hapla Notes: 992343a20b2SVaclav 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. 9934302210dSVaclav Hapla 994577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 99536bce6e8SMatthew G. Knepley @*/ 9964302210dSVaclav Hapla PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value) 99736bce6e8SMatthew G. Knepley { 9984302210dSVaclav Hapla char *parentAbsPath; 99960568592SMatthew G. Knepley hid_t h5, dataspace, obj, attribute, dtype; 1000080f144cSVaclav Hapla PetscBool has; 100136bce6e8SMatthew G. Knepley PetscErrorCode ierr; 100236bce6e8SMatthew G. Knepley 100336bce6e8SMatthew G. Knepley PetscFunctionBegin; 10045cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 10054302210dSVaclav Hapla if (parent) PetscValidCharPointer(parent, 2); 1006c57a1a9aSVaclav Hapla PetscValidCharPointer(name, 3); 10074302210dSVaclav Hapla PetscValidLogicalCollectiveEnum(viewer,datatype,4); 1008b7e81f0fSVaclav Hapla PetscValidPointer(value, 5); 10094302210dSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath);CHKERRQ(ierr); 10104302210dSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_TRUE, NULL, NULL);CHKERRQ(ierr); 10114302210dSVaclav Hapla ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has);CHKERRQ(ierr); 101236bce6e8SMatthew G. Knepley ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 10137e97c476SMatthew G. Knepley if (datatype == PETSC_STRING) { 10147e97c476SMatthew G. Knepley size_t len; 10157e97c476SMatthew G. Knepley ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr); 1016729ab6d9SBarry Smith PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 10177e97c476SMatthew G. Knepley } 101836bce6e8SMatthew G. Knepley ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 1019729ab6d9SBarry Smith PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR)); 10204302210dSVaclav Hapla PetscStackCallHDF5Return(obj,H5Oopen,(h5, parentAbsPath, H5P_DEFAULT)); 1021080f144cSVaclav Hapla if (has) { 1022080f144cSVaclav Hapla PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 1023080f144cSVaclav Hapla } else { 102460568592SMatthew G. Knepley PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT)); 1025080f144cSVaclav Hapla } 1026729ab6d9SBarry Smith PetscStackCallHDF5(H5Awrite,(attribute, dtype, value)); 1027729ab6d9SBarry Smith if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype)); 1028729ab6d9SBarry Smith PetscStackCallHDF5(H5Aclose,(attribute)); 102960568592SMatthew G. Knepley PetscStackCallHDF5(H5Oclose,(obj)); 1030729ab6d9SBarry Smith PetscStackCallHDF5(H5Sclose,(dataspace)); 10314302210dSVaclav Hapla ierr = PetscFree(parentAbsPath);CHKERRQ(ierr); 103236bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 103336bce6e8SMatthew G. Knepley } 103436bce6e8SMatthew G. Knepley 1035df863907SAlex Fikl /*@C 1036577e0e04SVaclav Hapla PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name 1037577e0e04SVaclav Hapla 1038343a20b2SVaclav Hapla Collective 1039343a20b2SVaclav Hapla 1040577e0e04SVaclav Hapla Input Parameters: 1041577e0e04SVaclav Hapla + viewer - The HDF5 viewer 1042577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1043577e0e04SVaclav Hapla . name - The attribute name 1044577e0e04SVaclav Hapla . datatype - The attribute type 1045577e0e04SVaclav Hapla - value - The attribute value 1046577e0e04SVaclav Hapla 1047577e0e04SVaclav Hapla Notes: 1048343a20b2SVaclav 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). 1049577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 1050577e0e04SVaclav Hapla 1051577e0e04SVaclav Hapla Level: advanced 1052577e0e04SVaclav Hapla 1053577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1054577e0e04SVaclav Hapla @*/ 1055577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value) 1056577e0e04SVaclav Hapla { 1057577e0e04SVaclav Hapla PetscErrorCode ierr; 1058577e0e04SVaclav Hapla 1059577e0e04SVaclav Hapla PetscFunctionBegin; 1060577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1061577e0e04SVaclav Hapla PetscValidHeader(obj,2); 1062577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 1063b7e81f0fSVaclav Hapla PetscValidPointer(value,5); 1064577e0e04SVaclav Hapla ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 1065577e0e04SVaclav Hapla ierr = PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr); 1066577e0e04SVaclav Hapla PetscFunctionReturn(0); 1067577e0e04SVaclav Hapla } 1068577e0e04SVaclav Hapla 1069577e0e04SVaclav Hapla /*@C 1070b8ee85a1SVaclav Hapla PetscViewerHDF5ReadAttribute - Read an attribute 1071ad0c61feSMatthew G. Knepley 1072343a20b2SVaclav Hapla Collective 1073343a20b2SVaclav Hapla 1074ad0c61feSMatthew G. Knepley Input Parameters: 1075ad0c61feSMatthew G. Knepley + viewer - The HDF5 viewer 10764302210dSVaclav Hapla . parent - The parent dataset/group name 1077ad0c61feSMatthew G. Knepley . name - The attribute name 1078a2d6be1bSVaclav Hapla . datatype - The attribute type 1079a2d6be1bSVaclav Hapla - defaultValue - The pointer to the default value 1080ad0c61feSMatthew G. Knepley 1081ad0c61feSMatthew G. Knepley Output Parameter: 1082a2d6be1bSVaclav Hapla . value - The pointer to the read HDF5 attribute value 1083ad0c61feSMatthew G. Knepley 1084a2d6be1bSVaclav Hapla Notes: 1085a2d6be1bSVaclav Hapla If defaultValue is NULL and the attribute is not found, an error occurs. 1086a2d6be1bSVaclav Hapla If defaultValue is not NULL and the attribute is not found, defaultValue is copied to value. 1087a2d6be1bSVaclav Hapla The pointers defaultValue and value can be the same; for instance 1088a2d6be1bSVaclav Hapla $ flg = PETSC_FALSE; 1089a2d6be1bSVaclav Hapla $ ierr = PetscViewerHDF5ReadAttribute(viewer,name,"attr",PETSC_BOOL,&flg,&flg);CHKERRQ(ierr); 1090a2d6be1bSVaclav Hapla is valid, but make sure the default value is initialized. 1091a2d6be1bSVaclav Hapla 1092a2d6be1bSVaclav Hapla If the datatype is PETSC_STRING, the output string is newly allocated so one must PetscFree() it when no longer needed. 1093708d4cb3SBarry Smith 1094343a20b2SVaclav 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. 1095ad0c61feSMatthew G. Knepley 1096343a20b2SVaclav Hapla Level: advanced 10974302210dSVaclav Hapla 1098577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1099ad0c61feSMatthew G. Knepley @*/ 1100d70ec8faSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *defaultValue, void *value) 1101ad0c61feSMatthew G. Knepley { 11024302210dSVaclav Hapla char *parentAbsPath; 1103a2d6be1bSVaclav Hapla hid_t h5, obj, attribute, dtype; 1104969748fdSVaclav Hapla PetscBool has; 1105ad0c61feSMatthew G. Knepley PetscErrorCode ierr; 1106ad0c61feSMatthew G. Knepley 1107ad0c61feSMatthew G. Knepley PetscFunctionBegin; 11085cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 11094302210dSVaclav Hapla if (parent) PetscValidCharPointer(parent, 2); 1110c57a1a9aSVaclav Hapla PetscValidCharPointer(name, 3); 1111a2d6be1bSVaclav Hapla if (defaultValue) PetscValidPointer(defaultValue, 5); 1112a2d6be1bSVaclav Hapla PetscValidPointer(value, 6); 1113a2d6be1bSVaclav Hapla ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 11144302210dSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath);CHKERRQ(ierr); 11154302210dSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, &has, NULL);CHKERRQ(ierr); 11164302210dSVaclav Hapla if (has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has);CHKERRQ(ierr);} 1117a2d6be1bSVaclav Hapla if (!has) { 1118a2d6be1bSVaclav Hapla if (defaultValue) { 1119a2d6be1bSVaclav Hapla if (defaultValue != value) { 1120a2d6be1bSVaclav Hapla if (datatype == PETSC_STRING) { 1121a2d6be1bSVaclav Hapla ierr = PetscStrallocpy(*(char**)defaultValue, (char**)value);CHKERRQ(ierr); 1122a2d6be1bSVaclav Hapla } else { 1123a2d6be1bSVaclav Hapla size_t len; 1124a2d6be1bSVaclav Hapla PetscStackCallHDF5Return(len,H5Tget_size,(dtype)); 1125a2d6be1bSVaclav Hapla ierr = PetscMemcpy(value, defaultValue, len);CHKERRQ(ierr); 1126a2d6be1bSVaclav Hapla } 1127a2d6be1bSVaclav Hapla } 1128a2d6be1bSVaclav Hapla ierr = PetscFree(parentAbsPath);CHKERRQ(ierr); 1129a2d6be1bSVaclav Hapla PetscFunctionReturn(0); 1130a2d6be1bSVaclav Hapla } else SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist and default value not provided", parentAbsPath, name); 1131a2d6be1bSVaclav Hapla } 1132ad0c61feSMatthew G. Knepley ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 11334302210dSVaclav Hapla PetscStackCallHDF5Return(obj,H5Oopen,(h5, parentAbsPath, H5P_DEFAULT)); 113460568592SMatthew G. Knepley PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 1135f0b58479SMatthew G. Knepley if (datatype == PETSC_STRING) { 1136f0b58479SMatthew G. Knepley size_t len; 1137a2d6be1bSVaclav Hapla hid_t atype; 1138e90c5075SMarek Pecha PetscStackCallHDF5Return(atype,H5Aget_type,(attribute)); 1139a2d6be1bSVaclav Hapla PetscStackCallHDF5Return(len,H5Tget_size,(atype)); 1140708d4cb3SBarry Smith ierr = PetscMalloc((len+1) * sizeof(char), value);CHKERRQ(ierr); 1141708d4cb3SBarry Smith PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 1142708d4cb3SBarry Smith PetscStackCallHDF5(H5Aread,(attribute, dtype, *(char**)value)); 1143708d4cb3SBarry Smith } else { 114470efba33SBarry Smith PetscStackCallHDF5(H5Aread,(attribute, dtype, value)); 1145708d4cb3SBarry Smith } 1146729ab6d9SBarry Smith PetscStackCallHDF5(H5Aclose,(attribute)); 1147e90c5075SMarek Pecha /* H5Oclose can be used to close groups, datasets, or committed datatypes */ 1148e90c5075SMarek Pecha PetscStackCallHDF5(H5Oclose,(obj)); 11494302210dSVaclav Hapla ierr = PetscFree(parentAbsPath);CHKERRQ(ierr); 1150ad0c61feSMatthew G. Knepley PetscFunctionReturn(0); 1151ad0c61feSMatthew G. Knepley } 1152ad0c61feSMatthew G. Knepley 1153577e0e04SVaclav Hapla /*@C 1154577e0e04SVaclav Hapla PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name 1155577e0e04SVaclav Hapla 1156343a20b2SVaclav Hapla Collective 1157343a20b2SVaclav Hapla 1158577e0e04SVaclav Hapla Input Parameters: 1159577e0e04SVaclav Hapla + viewer - The HDF5 viewer 1160577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1161577e0e04SVaclav Hapla . name - The attribute name 1162577e0e04SVaclav Hapla - datatype - The attribute type 1163577e0e04SVaclav Hapla 1164577e0e04SVaclav Hapla Output Parameter: 1165577e0e04SVaclav Hapla . value - The attribute value 1166577e0e04SVaclav Hapla 1167577e0e04SVaclav Hapla Notes: 1168577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1169577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 1170577e0e04SVaclav Hapla 1171577e0e04SVaclav Hapla Level: advanced 1172577e0e04SVaclav Hapla 1173577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute() PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1174577e0e04SVaclav Hapla @*/ 1175a2d6be1bSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *defaultValue, void *value) 1176577e0e04SVaclav Hapla { 1177577e0e04SVaclav Hapla PetscErrorCode ierr; 1178577e0e04SVaclav Hapla 1179577e0e04SVaclav Hapla PetscFunctionBegin; 1180577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1181577e0e04SVaclav Hapla PetscValidHeader(obj,2); 1182577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 1183064a246eSJacob Faibussowitsch PetscValidPointer(value, 6); 1184577e0e04SVaclav Hapla ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 1185a2d6be1bSVaclav Hapla ierr = PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, defaultValue, value);CHKERRQ(ierr); 1186577e0e04SVaclav Hapla PetscFunctionReturn(0); 1187577e0e04SVaclav Hapla } 1188577e0e04SVaclav Hapla 1189a75c3fd4SVaclav Hapla PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_) 1190a75c3fd4SVaclav Hapla { 1191a75c3fd4SVaclav Hapla htri_t exists; 1192a75c3fd4SVaclav Hapla hid_t group; 1193a75c3fd4SVaclav Hapla 1194a75c3fd4SVaclav Hapla PetscFunctionBegin; 1195a75c3fd4SVaclav Hapla PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT)); 1196a75c3fd4SVaclav Hapla if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT)); 1197a75c3fd4SVaclav Hapla if (!exists && createGroup) { 1198a75c3fd4SVaclav Hapla PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT)); 1199a75c3fd4SVaclav Hapla PetscStackCallHDF5(H5Gclose,(group)); 1200a75c3fd4SVaclav Hapla exists = PETSC_TRUE; 1201a75c3fd4SVaclav Hapla } 1202a75c3fd4SVaclav Hapla *exists_ = (PetscBool) exists; 1203a75c3fd4SVaclav Hapla PetscFunctionReturn(0); 1204a75c3fd4SVaclav Hapla } 1205a75c3fd4SVaclav Hapla 1206bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype) 12075cdeb108SMatthew G. Knepley { 120890e03537SVaclav Hapla const char rootGroupName[] = "/"; 12095cdeb108SMatthew G. Knepley hid_t h5; 1210e5bf9ebcSVaclav Hapla PetscBool exists=PETSC_FALSE; 121115b861d2SVaclav Hapla PetscInt i; 121215b861d2SVaclav Hapla int n; 121385688ae2SVaclav Hapla char **hierarchy; 121485688ae2SVaclav Hapla char buf[PETSC_MAX_PATH_LEN]=""; 12155cdeb108SMatthew G. Knepley PetscErrorCode ierr; 12165cdeb108SMatthew G. Knepley 12175cdeb108SMatthew G. Knepley PetscFunctionBegin; 12185cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 121990e03537SVaclav Hapla if (name) PetscValidCharPointer(name, 2); 122090e03537SVaclav Hapla else name = rootGroupName; 1221ccd66a83SVaclav Hapla if (has) { 1222064a246eSJacob Faibussowitsch PetscValidBoolPointer(has, 4); 12235cdeb108SMatthew G. Knepley *has = PETSC_FALSE; 1224ccd66a83SVaclav Hapla } 1225ccd66a83SVaclav Hapla if (otype) { 1226064a246eSJacob Faibussowitsch PetscValidIntPointer(otype, 5); 122756cc0592SVaclav Hapla *otype = H5O_TYPE_UNKNOWN; 1228ccd66a83SVaclav Hapla } 1229c0cd0301SJed Brown ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 123085688ae2SVaclav Hapla 123185688ae2SVaclav Hapla /* 123285688ae2SVaclav Hapla Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing. 123385688ae2SVaclav Hapla Hence, each of them needs to be tested separately: 123485688ae2SVaclav Hapla 1) whether it's a valid link 123585688ae2SVaclav Hapla 2) whether this link resolves to an object 123685688ae2SVaclav Hapla See H5Oexists_by_name() documentation. 123785688ae2SVaclav Hapla */ 123885688ae2SVaclav Hapla ierr = PetscStrToArray(name,'/',&n,&hierarchy);CHKERRQ(ierr); 123985688ae2SVaclav Hapla if (!n) { 124085688ae2SVaclav Hapla /* Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */ 1241ccd66a83SVaclav Hapla if (has) *has = PETSC_TRUE; 1242ccd66a83SVaclav Hapla if (otype) *otype = H5O_TYPE_GROUP; 124385688ae2SVaclav Hapla ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 124485688ae2SVaclav Hapla PetscFunctionReturn(0); 124585688ae2SVaclav Hapla } 124685688ae2SVaclav Hapla for (i=0; i<n; i++) { 124785688ae2SVaclav Hapla ierr = PetscStrcat(buf,"/");CHKERRQ(ierr); 124885688ae2SVaclav Hapla ierr = PetscStrcat(buf,hierarchy[i]);CHKERRQ(ierr); 1249a75c3fd4SVaclav Hapla ierr = PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);CHKERRQ(ierr); 1250a75c3fd4SVaclav Hapla if (!exists) break; 125185688ae2SVaclav Hapla } 125285688ae2SVaclav Hapla ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 125385688ae2SVaclav Hapla 125485688ae2SVaclav Hapla /* If the object exists, get its type */ 1255ccd66a83SVaclav Hapla if (exists && otype) { 12565cdeb108SMatthew G. Knepley H5O_info_t info; 12575cdeb108SMatthew G. Knepley 125876276748SVaclav Hapla /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */ 125904633f7fSVaclav Hapla PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT)); 126056cc0592SVaclav Hapla *otype = info.type; 12615cdeb108SMatthew G. Knepley } 1262ccd66a83SVaclav Hapla if (has) *has = exists; 12635cdeb108SMatthew G. Knepley PetscFunctionReturn(0); 12645cdeb108SMatthew G. Knepley } 12655cdeb108SMatthew G. Knepley 12664302210dSVaclav Hapla /*@C 12670a9f38efSVaclav Hapla PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file 12680a9f38efSVaclav Hapla 1269343a20b2SVaclav Hapla Collective 1270343a20b2SVaclav Hapla 12710a9f38efSVaclav Hapla Input Parameters: 12724302210dSVaclav Hapla + viewer - The HDF5 viewer 12734302210dSVaclav Hapla - path - The group path 12740a9f38efSVaclav Hapla 12750a9f38efSVaclav Hapla Output Parameter: 12760a9f38efSVaclav Hapla . has - Flag for group existence 12770a9f38efSVaclav Hapla 12780a9f38efSVaclav Hapla Level: advanced 12790a9f38efSVaclav Hapla 12804302210dSVaclav Hapla Notes: 1281785c443eSVaclav 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. 1282785c443eSVaclav Hapla So NULL or empty path means the current pushed group. 1283343a20b2SVaclav Hapla If path exists but is not a group, PETSC_FALSE is returned. 12844302210dSVaclav Hapla 128589e0ef10SVaclav Hapla .seealso: PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasDataset(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5GetGroup(), PetscViewerHDF5OpenGroup() 12860a9f38efSVaclav Hapla @*/ 12874302210dSVaclav Hapla PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, const char path[], PetscBool *has) 12880a9f38efSVaclav Hapla { 12890a9f38efSVaclav Hapla H5O_type_t type; 12904302210dSVaclav Hapla char *abspath; 12910a9f38efSVaclav Hapla PetscErrorCode ierr; 12920a9f38efSVaclav Hapla 12930a9f38efSVaclav Hapla PetscFunctionBegin; 12940a9f38efSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 12954302210dSVaclav Hapla if (path) PetscValidCharPointer(path,2); 12964302210dSVaclav Hapla PetscValidBoolPointer(has,3); 12974302210dSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, path, &abspath);CHKERRQ(ierr); 12984302210dSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type);CHKERRQ(ierr); 12994302210dSVaclav Hapla *has = (PetscBool)(type == H5O_TYPE_GROUP); 13004302210dSVaclav Hapla ierr = PetscFree(abspath);CHKERRQ(ierr); 13010a9f38efSVaclav Hapla PetscFunctionReturn(0); 13020a9f38efSVaclav Hapla } 13030a9f38efSVaclav Hapla 130489e0ef10SVaclav Hapla /*@C 130589e0ef10SVaclav Hapla PetscViewerHDF5HasDataset - Check whether a given dataset exists in the HDF5 file 130689e0ef10SVaclav Hapla 1307343a20b2SVaclav Hapla Collective 1308343a20b2SVaclav Hapla 130989e0ef10SVaclav Hapla Input Parameters: 131089e0ef10SVaclav Hapla + viewer - The HDF5 viewer 131189e0ef10SVaclav Hapla - path - The dataset path 131289e0ef10SVaclav Hapla 131389e0ef10SVaclav Hapla Output Parameter: 131489e0ef10SVaclav Hapla . has - Flag whether dataset exists 131589e0ef10SVaclav Hapla 131689e0ef10SVaclav Hapla Level: advanced 131789e0ef10SVaclav Hapla 131889e0ef10SVaclav Hapla Notes: 1319343a20b2SVaclav 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. 132089e0ef10SVaclav Hapla If path is NULL or empty, has is set to PETSC_FALSE. 1321343a20b2SVaclav Hapla If path exists but is not a dataset, has is set to PETSC_FALSE as well. 132289e0ef10SVaclav Hapla 132389e0ef10SVaclav Hapla .seealso: PetscViewerHDF5HasObject(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasGroup(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5GetGroup() 132489e0ef10SVaclav Hapla @*/ 132589e0ef10SVaclav Hapla PetscErrorCode PetscViewerHDF5HasDataset(PetscViewer viewer, const char path[], PetscBool *has) 132689e0ef10SVaclav Hapla { 132789e0ef10SVaclav Hapla H5O_type_t type; 132889e0ef10SVaclav Hapla char *abspath; 132989e0ef10SVaclav Hapla PetscErrorCode ierr; 133089e0ef10SVaclav Hapla 133189e0ef10SVaclav Hapla PetscFunctionBegin; 133289e0ef10SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 133389e0ef10SVaclav Hapla if (path) PetscValidCharPointer(path,2); 133489e0ef10SVaclav Hapla PetscValidBoolPointer(has,3); 133589e0ef10SVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, path, &abspath);CHKERRQ(ierr); 133689e0ef10SVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type);CHKERRQ(ierr); 133789e0ef10SVaclav Hapla *has = (PetscBool)(type == H5O_TYPE_DATASET); 133889e0ef10SVaclav Hapla ierr = PetscFree(abspath);CHKERRQ(ierr); 133989e0ef10SVaclav Hapla PetscFunctionReturn(0); 134089e0ef10SVaclav Hapla } 134189e0ef10SVaclav Hapla 13420a9f38efSVaclav Hapla /*@ 1343e3f143f7SVaclav Hapla PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group 1344ecce1506SVaclav Hapla 1345343a20b2SVaclav Hapla Collective 1346343a20b2SVaclav Hapla 1347ecce1506SVaclav Hapla Input Parameters: 1348ecce1506SVaclav Hapla + viewer - The HDF5 viewer 1349ecce1506SVaclav Hapla - obj - The named object 1350ecce1506SVaclav Hapla 1351ecce1506SVaclav Hapla Output Parameter: 135289e0ef10SVaclav Hapla . has - Flag for dataset existence 1353ecce1506SVaclav Hapla 1354e3f143f7SVaclav Hapla Notes: 135589e0ef10SVaclav Hapla If the object is unnamed, an error occurs. 1356343a20b2SVaclav Hapla If the path current_group/object_name exists but is not a dataset, has is set to PETSC_FALSE as well. 1357e3f143f7SVaclav Hapla 1358ecce1506SVaclav Hapla Level: advanced 1359ecce1506SVaclav Hapla 136089e0ef10SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasDataset(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5GetGroup() 1361ecce1506SVaclav Hapla @*/ 1362ecce1506SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has) 1363ecce1506SVaclav Hapla { 136489e0ef10SVaclav Hapla size_t len; 1365ecce1506SVaclav Hapla PetscErrorCode ierr; 1366ecce1506SVaclav Hapla 1367ecce1506SVaclav Hapla PetscFunctionBegin; 1368c57a1a9aSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1369c57a1a9aSVaclav Hapla PetscValidHeader(obj,2); 13704302210dSVaclav Hapla PetscValidBoolPointer(has,3); 137189e0ef10SVaclav Hapla ierr = PetscStrlen(obj->name, &len);CHKERRQ(ierr); 137289e0ef10SVaclav Hapla if (!len) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named"); 137389e0ef10SVaclav Hapla ierr = PetscViewerHDF5HasDataset(viewer, obj->name, has);CHKERRQ(ierr); 1374ecce1506SVaclav Hapla PetscFunctionReturn(0); 1375ecce1506SVaclav Hapla } 1376ecce1506SVaclav Hapla 1377df863907SAlex Fikl /*@C 1378ece83b6cSVaclav Hapla PetscViewerHDF5HasAttribute - Check whether an attribute exists 1379e7bdbf8eSMatthew G. Knepley 1380343a20b2SVaclav Hapla Collective 1381343a20b2SVaclav Hapla 1382e7bdbf8eSMatthew G. Knepley Input Parameters: 1383e7bdbf8eSMatthew G. Knepley + viewer - The HDF5 viewer 13844302210dSVaclav Hapla . parent - The parent dataset/group name 1385e7bdbf8eSMatthew G. Knepley - name - The attribute name 1386e7bdbf8eSMatthew G. Knepley 1387e7bdbf8eSMatthew G. Knepley Output Parameter: 1388e7bdbf8eSMatthew G. Knepley . has - Flag for attribute existence 1389e7bdbf8eSMatthew G. Knepley 1390e7bdbf8eSMatthew G. Knepley Level: advanced 1391e7bdbf8eSMatthew G. Knepley 13924302210dSVaclav Hapla Notes: 1393343a20b2SVaclav 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. 13944302210dSVaclav Hapla 1395577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1396e7bdbf8eSMatthew G. Knepley @*/ 13974302210dSVaclav Hapla PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 1398e7bdbf8eSMatthew G. Knepley { 13994302210dSVaclav Hapla char *parentAbsPath; 1400e7bdbf8eSMatthew G. Knepley PetscErrorCode ierr; 1401e7bdbf8eSMatthew G. Knepley 1402e7bdbf8eSMatthew G. Knepley PetscFunctionBegin; 14035cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 14044302210dSVaclav Hapla if (parent) PetscValidCharPointer(parent,2); 1405c57a1a9aSVaclav Hapla PetscValidCharPointer(name,3); 14064302210dSVaclav Hapla PetscValidBoolPointer(has,4); 14074302210dSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath);CHKERRQ(ierr); 14084302210dSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, has, NULL);CHKERRQ(ierr); 14094302210dSVaclav Hapla if (*has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, has);CHKERRQ(ierr);} 14104302210dSVaclav Hapla ierr = PetscFree(parentAbsPath);CHKERRQ(ierr); 141106db490cSVaclav Hapla PetscFunctionReturn(0); 141206db490cSVaclav Hapla } 141306db490cSVaclav Hapla 1414577e0e04SVaclav Hapla /*@C 1415577e0e04SVaclav Hapla PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name 1416577e0e04SVaclav Hapla 1417343a20b2SVaclav Hapla Collective 1418343a20b2SVaclav Hapla 1419577e0e04SVaclav Hapla Input Parameters: 1420577e0e04SVaclav Hapla + viewer - The HDF5 viewer 1421577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1422577e0e04SVaclav Hapla - name - The attribute name 1423577e0e04SVaclav Hapla 1424577e0e04SVaclav Hapla Output Parameter: 1425577e0e04SVaclav Hapla . has - Flag for attribute existence 1426577e0e04SVaclav Hapla 1427577e0e04SVaclav Hapla Notes: 1428577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1429577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 1430577e0e04SVaclav Hapla 1431577e0e04SVaclav Hapla Level: advanced 1432577e0e04SVaclav Hapla 1433577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1434577e0e04SVaclav Hapla @*/ 1435577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has) 1436577e0e04SVaclav Hapla { 1437577e0e04SVaclav Hapla PetscErrorCode ierr; 1438577e0e04SVaclav Hapla 1439577e0e04SVaclav Hapla PetscFunctionBegin; 1440577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1441577e0e04SVaclav Hapla PetscValidHeader(obj,2); 1442577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 14434302210dSVaclav Hapla PetscValidBoolPointer(has,4); 1444577e0e04SVaclav Hapla ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 1445577e0e04SVaclav Hapla ierr = PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);CHKERRQ(ierr); 1446577e0e04SVaclav Hapla PetscFunctionReturn(0); 1447577e0e04SVaclav Hapla } 1448577e0e04SVaclav Hapla 144906db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 145006db490cSVaclav Hapla { 145106db490cSVaclav Hapla hid_t h5; 145206db490cSVaclav Hapla htri_t hhas; 145306db490cSVaclav Hapla PetscErrorCode ierr; 145406db490cSVaclav Hapla 145506db490cSVaclav Hapla PetscFunctionBegin; 145606db490cSVaclav Hapla ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 14572f936e54SVaclav Hapla PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT)); 14585cdeb108SMatthew G. Knepley *has = hhas ? PETSC_TRUE : PETSC_FALSE; 1459e7bdbf8eSMatthew G. Knepley PetscFunctionReturn(0); 1460e7bdbf8eSMatthew G. Knepley } 1461e7bdbf8eSMatthew G. Knepley 1462a75e6a4aSMatthew G. Knepley /* 1463a75e6a4aSMatthew G. Knepley The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that 1464a75e6a4aSMatthew G. Knepley is attached to a communicator, in this case the attribute is a PetscViewer. 1465a75e6a4aSMatthew G. Knepley */ 1466d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID; 1467a75e6a4aSMatthew G. Knepley 1468a75e6a4aSMatthew G. Knepley /*@C 1469a75e6a4aSMatthew G. Knepley PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator. 1470a75e6a4aSMatthew G. Knepley 1471d083f849SBarry Smith Collective 1472a75e6a4aSMatthew G. Knepley 1473a75e6a4aSMatthew G. Knepley Input Parameter: 1474a75e6a4aSMatthew G. Knepley . comm - the MPI communicator to share the HDF5 PetscViewer 1475a75e6a4aSMatthew G. Knepley 1476a75e6a4aSMatthew G. Knepley Level: intermediate 1477a75e6a4aSMatthew G. Knepley 1478a75e6a4aSMatthew G. Knepley Options Database Keys: 147910699b91SBarry Smith . -viewer_hdf5_filename <name> - name of the HDF5 file 1480a75e6a4aSMatthew G. Knepley 1481a75e6a4aSMatthew G. Knepley Environmental variables: 148210699b91SBarry Smith . PETSC_VIEWER_HDF5_FILENAME - name of the HDF5 file 1483a75e6a4aSMatthew G. Knepley 1484a75e6a4aSMatthew G. Knepley Notes: 1485a75e6a4aSMatthew G. Knepley Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return 1486a75e6a4aSMatthew G. Knepley an error code. The HDF5 PetscViewer is usually used in the form 1487a75e6a4aSMatthew G. Knepley $ XXXView(XXX object, PETSC_VIEWER_HDF5_(comm)); 1488a75e6a4aSMatthew G. Knepley 1489a75e6a4aSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy() 1490a75e6a4aSMatthew G. Knepley @*/ 1491a75e6a4aSMatthew G. Knepley PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm) 1492a75e6a4aSMatthew G. Knepley { 1493a75e6a4aSMatthew G. Knepley PetscErrorCode ierr; 1494a75e6a4aSMatthew G. Knepley PetscBool flg; 1495a75e6a4aSMatthew G. Knepley PetscViewer viewer; 1496a75e6a4aSMatthew G. Knepley char fname[PETSC_MAX_PATH_LEN]; 1497a75e6a4aSMatthew G. Knepley MPI_Comm ncomm; 1498a75e6a4aSMatthew G. Knepley 1499a75e6a4aSMatthew G. Knepley PetscFunctionBegin; 1500908793a3SLisandro 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);} 1501a75e6a4aSMatthew G. Knepley if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) { 1502908793a3SLisandro Dalcin ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,NULL); 1503908793a3SLisandro Dalcin if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 1504a75e6a4aSMatthew G. Knepley } 150547435625SJed Brown ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg); 1506908793a3SLisandro Dalcin if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 1507a75e6a4aSMatthew G. Knepley if (!flg) { /* PetscViewer not yet created */ 1508a75e6a4aSMatthew G. Knepley ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg); 15092cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1510a75e6a4aSMatthew G. Knepley if (!flg) { 1511a75e6a4aSMatthew G. Knepley ierr = PetscStrcpy(fname,"output.h5"); 15122cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1513a75e6a4aSMatthew G. Knepley } 1514a75e6a4aSMatthew G. Knepley ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer); 15152cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1516a75e6a4aSMatthew G. Knepley ierr = PetscObjectRegisterDestroy((PetscObject)viewer); 15172cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 151847435625SJed Brown ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer); 1519908793a3SLisandro Dalcin if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 1520a75e6a4aSMatthew G. Knepley } 1521a75e6a4aSMatthew G. Knepley ierr = PetscCommDestroy(&ncomm); 15222cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1523a75e6a4aSMatthew G. Knepley PetscFunctionReturn(viewer); 1524a75e6a4aSMatthew G. Knepley } 1525