1bc307625SVaclav Hapla #include <petsc/private/viewerimpl.h> 27d964842SVaclav Hapla #include <petsc/private/viewerhdf5impl.h> 3d70abbfaSBarry Smith #include <petscviewerhdf5.h> /*I "petscviewerhdf5.h" I*/ 45c6c1daeSBarry Smith 5bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer, const char[], PetscBool, PetscBool*, H5O_type_t*); 606db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer, const char[], const char[], PetscBool*); 706db490cSVaclav Hapla 8*4302210dSVaclav Hapla static PetscErrorCode PetscViewerHDF5GetAbsolutePath_Internal(PetscViewer viewer, const char path[], char *abspath[]) 96c132bc1SVaclav Hapla { 10*4302210dSVaclav Hapla PetscBool relative = PETSC_FALSE; 116c132bc1SVaclav Hapla const char *group; 126c132bc1SVaclav Hapla char buf[PETSC_MAX_PATH_LEN] = ""; 136c132bc1SVaclav Hapla PetscErrorCode ierr; 146c132bc1SVaclav Hapla 156c132bc1SVaclav Hapla PetscFunctionBegin; 166c132bc1SVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr); 17*4302210dSVaclav Hapla ierr = PetscViewerHDF5PathIsRelative(path, PETSC_TRUE, &relative);CHKERRQ(ierr); 18*4302210dSVaclav Hapla if (relative) { 19*4302210dSVaclav Hapla ierr = PetscStrcpy(buf, group);CHKERRQ(ierr); 206c132bc1SVaclav Hapla ierr = PetscStrcat(buf, "/");CHKERRQ(ierr); 21*4302210dSVaclav Hapla ierr = PetscStrcat(buf, path);CHKERRQ(ierr); 22*4302210dSVaclav Hapla ierr = PetscStrallocpy(buf, abspath);CHKERRQ(ierr); 23*4302210dSVaclav Hapla } else { 24*4302210dSVaclav Hapla ierr = PetscStrallocpy(path, abspath);CHKERRQ(ierr); 25*4302210dSVaclav Hapla } 266c132bc1SVaclav Hapla PetscFunctionReturn(0); 276c132bc1SVaclav Hapla } 286c132bc1SVaclav Hapla 29577e0e04SVaclav Hapla static PetscErrorCode PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer, PetscObject obj) 30577e0e04SVaclav Hapla { 31577e0e04SVaclav Hapla PetscBool has; 32577e0e04SVaclav Hapla const char *group; 33577e0e04SVaclav Hapla PetscErrorCode ierr; 34577e0e04SVaclav Hapla 35577e0e04SVaclav Hapla PetscFunctionBegin; 36577e0e04SVaclav Hapla if (!obj->name) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named"); 37577e0e04SVaclav Hapla ierr = PetscViewerHDF5HasObject(viewer, obj, &has);CHKERRQ(ierr); 38577e0e04SVaclav Hapla if (!has) { 39577e0e04SVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr); 40577e0e04SVaclav Hapla SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) %s not stored in group %s", obj->name, group); 41577e0e04SVaclav Hapla } 42577e0e04SVaclav Hapla PetscFunctionReturn(0); 43577e0e04SVaclav Hapla } 44577e0e04SVaclav Hapla 454416b707SBarry Smith static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptionItems *PetscOptionsObject,PetscViewer v) 4682ea9b62SBarry Smith { 4782ea9b62SBarry Smith PetscErrorCode ierr; 48ee8b9732SVaclav Hapla PetscBool flg = PETSC_FALSE, set; 4982ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data; 5082ea9b62SBarry Smith 5182ea9b62SBarry Smith PetscFunctionBegin; 5282ea9b62SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");CHKERRQ(ierr); 5382ea9b62SBarry Smith ierr = PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);CHKERRQ(ierr); 549a0502c6SHåkon Strandenes ierr = PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL);CHKERRQ(ierr); 55ee8b9732SVaclav Hapla ierr = PetscOptionsBool("-viewer_hdf5_collective","Enable collective transfer mode","PetscViewerHDF5SetCollective",flg,&flg,&set);CHKERRQ(ierr); 56ee8b9732SVaclav Hapla if (set) {ierr = PetscViewerHDF5SetCollective(v,flg);CHKERRQ(ierr);} 5782ea9b62SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 5882ea9b62SBarry Smith PetscFunctionReturn(0); 5982ea9b62SBarry Smith } 6082ea9b62SBarry Smith 611b793a25SVaclav Hapla static PetscErrorCode PetscViewerView_HDF5(PetscViewer v,PetscViewer viewer) 621b793a25SVaclav Hapla { 631b793a25SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data; 6403fa8834SVaclav Hapla PetscBool flg; 651b793a25SVaclav Hapla PetscErrorCode ierr; 661b793a25SVaclav Hapla 671b793a25SVaclav Hapla PetscFunctionBegin; 681b793a25SVaclav Hapla if (hdf5->filename) { 691b793a25SVaclav Hapla ierr = PetscViewerASCIIPrintf(viewer,"Filename: %s\n",hdf5->filename);CHKERRQ(ierr); 701b793a25SVaclav Hapla } 711b793a25SVaclav Hapla ierr = PetscViewerASCIIPrintf(viewer,"Vectors with blocksize 1 saved as 2D datasets: %s\n",PetscBools[hdf5->basedimension2]);CHKERRQ(ierr); 721b793a25SVaclav Hapla ierr = PetscViewerASCIIPrintf(viewer,"Enforce single precision storage: %s\n",PetscBools[hdf5->spoutput]);CHKERRQ(ierr); 7303fa8834SVaclav Hapla ierr = PetscViewerHDF5GetCollective(v,&flg);CHKERRQ(ierr); 7403fa8834SVaclav Hapla ierr = PetscViewerASCIIPrintf(viewer,"MPI-IO transfer mode: %s\n",flg ? "collective" : "independent");CHKERRQ(ierr); 751b793a25SVaclav Hapla PetscFunctionReturn(0); 761b793a25SVaclav Hapla } 771b793a25SVaclav Hapla 785c6c1daeSBarry Smith static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer) 795c6c1daeSBarry Smith { 805c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data; 815c6c1daeSBarry Smith PetscErrorCode ierr; 825c6c1daeSBarry Smith 835c6c1daeSBarry Smith PetscFunctionBegin; 845c6c1daeSBarry Smith ierr = PetscFree(hdf5->filename);CHKERRQ(ierr); 85729ab6d9SBarry Smith if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 865c6c1daeSBarry Smith PetscFunctionReturn(0); 875c6c1daeSBarry Smith } 885c6c1daeSBarry Smith 899b1403beSVaclav Hapla static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer) 905c6c1daeSBarry Smith { 915c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 925c6c1daeSBarry Smith PetscErrorCode ierr; 935c6c1daeSBarry Smith 945c6c1daeSBarry Smith PetscFunctionBegin; 959c5072fbSVaclav Hapla PetscStackCallHDF5(H5Pclose,(hdf5->dxpl_id)); 965c6c1daeSBarry Smith ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr); 975c6c1daeSBarry Smith while (hdf5->groups) { 987d964842SVaclav Hapla PetscViewerHDF5GroupList *tmp = hdf5->groups->next; 995c6c1daeSBarry Smith 1005c6c1daeSBarry Smith ierr = PetscFree(hdf5->groups->name);CHKERRQ(ierr); 1015c6c1daeSBarry Smith ierr = PetscFree(hdf5->groups);CHKERRQ(ierr); 1025c6c1daeSBarry Smith hdf5->groups = tmp; 1035c6c1daeSBarry Smith } 1045c6c1daeSBarry Smith ierr = PetscFree(hdf5);CHKERRQ(ierr); 1050b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr); 106d1232d7fSToby Isaac ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr); 1070b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr); 108058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetBaseDimension2_C",NULL);CHKERRQ(ierr); 109058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetSPOutput_C",NULL);CHKERRQ(ierr); 1105c6c1daeSBarry Smith PetscFunctionReturn(0); 1115c6c1daeSBarry Smith } 1125c6c1daeSBarry Smith 1139b1403beSVaclav Hapla static PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type) 1145c6c1daeSBarry Smith { 1155c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 1165c6c1daeSBarry Smith 1175c6c1daeSBarry Smith PetscFunctionBegin; 1185c6c1daeSBarry Smith hdf5->btype = type; 1195c6c1daeSBarry Smith PetscFunctionReturn(0); 1205c6c1daeSBarry Smith } 1215c6c1daeSBarry Smith 1220b62783dSVaclav Hapla static PetscErrorCode PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type) 1230b62783dSVaclav Hapla { 1240b62783dSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 1250b62783dSVaclav Hapla 1260b62783dSVaclav Hapla PetscFunctionBegin; 1270b62783dSVaclav Hapla *type = hdf5->btype; 1280b62783dSVaclav Hapla PetscFunctionReturn(0); 1290b62783dSVaclav Hapla } 1300b62783dSVaclav Hapla 1319b1403beSVaclav Hapla static PetscErrorCode PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg) 13282ea9b62SBarry Smith { 13382ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 13482ea9b62SBarry Smith 13582ea9b62SBarry Smith PetscFunctionBegin; 13682ea9b62SBarry Smith hdf5->basedimension2 = flg; 13782ea9b62SBarry Smith PetscFunctionReturn(0); 13882ea9b62SBarry Smith } 13982ea9b62SBarry Smith 140df863907SAlex Fikl /*@ 14182ea9b62SBarry Smith PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 14282ea9b62SBarry Smith dimension of 2. 14382ea9b62SBarry Smith 14482ea9b62SBarry Smith Logically Collective on PetscViewer 14582ea9b62SBarry Smith 14682ea9b62SBarry Smith Input Parameters: 14782ea9b62SBarry Smith + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 14882ea9b62SBarry 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 14982ea9b62SBarry Smith 15082ea9b62SBarry Smith Options Database: 15182ea9b62SBarry 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 15282ea9b62SBarry Smith 15382ea9b62SBarry Smith 15495452b02SPatrick Sanan Notes: 15595452b02SPatrick 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 15682ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 15782ea9b62SBarry Smith 15882ea9b62SBarry Smith Level: intermediate 15982ea9b62SBarry Smith 16082ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 16182ea9b62SBarry Smith 16282ea9b62SBarry Smith @*/ 16382ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg) 16482ea9b62SBarry Smith { 16582ea9b62SBarry Smith PetscErrorCode ierr; 16682ea9b62SBarry Smith 16782ea9b62SBarry Smith PetscFunctionBegin; 16882ea9b62SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 16982ea9b62SBarry Smith ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 17082ea9b62SBarry Smith PetscFunctionReturn(0); 17182ea9b62SBarry Smith } 17282ea9b62SBarry Smith 173df863907SAlex Fikl /*@ 17482ea9b62SBarry Smith PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 17582ea9b62SBarry Smith dimension of 2. 17682ea9b62SBarry Smith 17782ea9b62SBarry Smith Logically Collective on PetscViewer 17882ea9b62SBarry Smith 17982ea9b62SBarry Smith Input Parameter: 18082ea9b62SBarry Smith . viewer - the PetscViewer, must be of type HDF5 18182ea9b62SBarry Smith 18282ea9b62SBarry Smith Output Parameter: 18382ea9b62SBarry 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 18482ea9b62SBarry Smith 18595452b02SPatrick Sanan Notes: 18695452b02SPatrick 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 18782ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 18882ea9b62SBarry Smith 18982ea9b62SBarry Smith Level: intermediate 19082ea9b62SBarry Smith 19182ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 19282ea9b62SBarry Smith 19382ea9b62SBarry Smith @*/ 19482ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg) 19582ea9b62SBarry Smith { 19682ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 19782ea9b62SBarry Smith 19882ea9b62SBarry Smith PetscFunctionBegin; 19982ea9b62SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 20082ea9b62SBarry Smith *flg = hdf5->basedimension2; 20182ea9b62SBarry Smith PetscFunctionReturn(0); 20282ea9b62SBarry Smith } 20382ea9b62SBarry Smith 2049b1403beSVaclav Hapla static PetscErrorCode PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg) 2059a0502c6SHåkon Strandenes { 2069a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 2079a0502c6SHåkon Strandenes 2089a0502c6SHåkon Strandenes PetscFunctionBegin; 2099a0502c6SHåkon Strandenes hdf5->spoutput = flg; 2109a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2119a0502c6SHåkon Strandenes } 2129a0502c6SHåkon Strandenes 213df863907SAlex Fikl /*@ 2149a0502c6SHåkon Strandenes PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is 2159a0502c6SHåkon Strandenes compiled with double precision PetscReal. 2169a0502c6SHåkon Strandenes 2179a0502c6SHåkon Strandenes Logically Collective on PetscViewer 2189a0502c6SHåkon Strandenes 2199a0502c6SHåkon Strandenes Input Parameters: 2209a0502c6SHåkon Strandenes + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 2219a0502c6SHåkon Strandenes - flg - if PETSC_TRUE the data will be written to disk with single precision 2229a0502c6SHåkon Strandenes 2239a0502c6SHåkon Strandenes Options Database: 2249a0502c6SHåkon Strandenes . -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision 2259a0502c6SHåkon Strandenes 2269a0502c6SHåkon Strandenes 22795452b02SPatrick Sanan Notes: 22895452b02SPatrick Sanan Setting this option does not make any difference if PETSc is compiled with single precision 2299a0502c6SHåkon Strandenes in the first place. It does not affect reading datasets (HDF5 handle this internally). 2309a0502c6SHåkon Strandenes 2319a0502c6SHåkon Strandenes Level: intermediate 2329a0502c6SHåkon Strandenes 2339a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 2349a0502c6SHåkon Strandenes PetscReal 2359a0502c6SHåkon Strandenes 2369a0502c6SHåkon Strandenes @*/ 2379a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg) 2389a0502c6SHåkon Strandenes { 2399a0502c6SHåkon Strandenes PetscErrorCode ierr; 2409a0502c6SHåkon Strandenes 2419a0502c6SHåkon Strandenes PetscFunctionBegin; 2429a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 2439a0502c6SHåkon Strandenes ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 2449a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2459a0502c6SHåkon Strandenes } 2469a0502c6SHåkon Strandenes 247df863907SAlex Fikl /*@ 2489a0502c6SHåkon Strandenes PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is 2499a0502c6SHåkon Strandenes compiled with double precision PetscReal. 2509a0502c6SHåkon Strandenes 2519a0502c6SHåkon Strandenes Logically Collective on PetscViewer 2529a0502c6SHåkon Strandenes 2539a0502c6SHåkon Strandenes Input Parameter: 2549a0502c6SHåkon Strandenes . viewer - the PetscViewer, must be of type HDF5 2559a0502c6SHåkon Strandenes 2569a0502c6SHåkon Strandenes Output Parameter: 2579a0502c6SHåkon Strandenes . flg - if PETSC_TRUE the data will be written to disk with single precision 2589a0502c6SHåkon Strandenes 25995452b02SPatrick Sanan Notes: 26095452b02SPatrick Sanan Setting this option does not make any difference if PETSc is compiled with single precision 2619a0502c6SHåkon Strandenes in the first place. It does not affect reading datasets (HDF5 handle this internally). 2629a0502c6SHåkon Strandenes 2639a0502c6SHåkon Strandenes Level: intermediate 2649a0502c6SHåkon Strandenes 2659a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 2669a0502c6SHåkon Strandenes PetscReal 2679a0502c6SHåkon Strandenes 2689a0502c6SHåkon Strandenes @*/ 2699a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg) 2709a0502c6SHåkon Strandenes { 2719a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 2729a0502c6SHåkon Strandenes 2739a0502c6SHåkon Strandenes PetscFunctionBegin; 2749a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 2759a0502c6SHåkon Strandenes *flg = hdf5->spoutput; 2769a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2779a0502c6SHåkon Strandenes } 2789a0502c6SHåkon Strandenes 279ee8b9732SVaclav Hapla static PetscErrorCode PetscViewerHDF5SetCollective_HDF5(PetscViewer viewer, PetscBool flg) 280ee8b9732SVaclav Hapla { 281ee8b9732SVaclav Hapla PetscFunctionBegin; 282ee8b9732SVaclav Hapla /* H5FD_MPIO_COLLECTIVE is wrong in hdf5 1.10.2, and is the same as H5FD_MPIO_INDEPENDENT in earlier versions 283ee8b9732SVaclav Hapla - see e.g. https://gitlab.cosma.dur.ac.uk/swift/swiftsim/issues/431 */ 284c796909eSBarry Smith #if H5_VERSION_GE(1,10,3) && defined(H5_HAVE_PARALLEL) 285ee8b9732SVaclav Hapla { 286ee8b9732SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 287ee8b9732SVaclav Hapla PetscStackCallHDF5(H5Pset_dxpl_mpio,(hdf5->dxpl_id, flg ? H5FD_MPIO_COLLECTIVE : H5FD_MPIO_INDEPENDENT)); 288ee8b9732SVaclav Hapla } 289ee8b9732SVaclav Hapla #else 290ee8b9732SVaclav Hapla if (flg) { 291ee8b9732SVaclav Hapla PetscErrorCode ierr; 292c796909eSBarry 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); 293ee8b9732SVaclav Hapla } 294ee8b9732SVaclav Hapla #endif 295ee8b9732SVaclav Hapla PetscFunctionReturn(0); 296ee8b9732SVaclav Hapla } 297ee8b9732SVaclav Hapla 298ee8b9732SVaclav Hapla /*@ 299ee8b9732SVaclav Hapla PetscViewerHDF5SetCollective - Use collective MPI-IO transfer mode for HDF5 reads and writes. 300ee8b9732SVaclav Hapla 301ee8b9732SVaclav Hapla Logically Collective; flg must contain common value 302ee8b9732SVaclav Hapla 303ee8b9732SVaclav Hapla Input Parameters: 304ee8b9732SVaclav Hapla + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 305ee8b9732SVaclav Hapla - flg - PETSC_TRUE for collective mode; PETSC_FALSE for independent mode (default) 306ee8b9732SVaclav Hapla 307ee8b9732SVaclav Hapla Options Database: 308ee8b9732SVaclav Hapla . -viewer_hdf5_collective - turns on (true) or off (false) collective transfers 309ee8b9732SVaclav Hapla 310ee8b9732SVaclav Hapla Notes: 311ee8b9732SVaclav Hapla Collective mode gives the MPI-IO layer underneath HDF5 a chance to do some additional collective optimizations and hence can perform better. 31253effed7SBarry 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. 313ee8b9732SVaclav Hapla 314ee8b9732SVaclav Hapla Developer notes: 315ee8b9732SVaclav Hapla In the HDF5 layer, PETSC_TRUE / PETSC_FALSE means H5Pset_dxpl_mpio() is called with H5FD_MPIO_COLLECTIVE / H5FD_MPIO_INDEPENDENT, respectively. 316ee8b9732SVaclav Hapla This in turn means use of MPI_File_{read,write}_all / MPI_File_{read,write} in the MPI-IO layer, respectively. 317ee8b9732SVaclav Hapla See HDF5 documentation and MPI-IO documentation for details. 318ee8b9732SVaclav Hapla 319ee8b9732SVaclav Hapla Level: intermediate 320ee8b9732SVaclav Hapla 321ee8b9732SVaclav Hapla .seealso: PetscViewerHDF5GetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open() 322ee8b9732SVaclav Hapla 323ee8b9732SVaclav Hapla @*/ 324ee8b9732SVaclav Hapla PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer viewer,PetscBool flg) 325ee8b9732SVaclav Hapla { 326ee8b9732SVaclav Hapla PetscErrorCode ierr; 327ee8b9732SVaclav Hapla 328ee8b9732SVaclav Hapla PetscFunctionBegin; 329ee8b9732SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 330ee8b9732SVaclav Hapla PetscValidLogicalCollectiveBool(viewer,flg,2); 331ee8b9732SVaclav Hapla ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetCollective_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 332ee8b9732SVaclav Hapla PetscFunctionReturn(0); 333ee8b9732SVaclav Hapla } 334ee8b9732SVaclav Hapla 335ee8b9732SVaclav Hapla static PetscErrorCode PetscViewerHDF5GetCollective_HDF5(PetscViewer viewer, PetscBool *flg) 336ee8b9732SVaclav Hapla { 337c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL) 338ee8b9732SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 339ee8b9732SVaclav Hapla H5FD_mpio_xfer_t mode; 340c796909eSBarry Smith #endif 341ee8b9732SVaclav Hapla 342ee8b9732SVaclav Hapla PetscFunctionBegin; 343c796909eSBarry Smith #if !defined(H5_HAVE_PARALLEL) 344c796909eSBarry Smith *flg = PETSC_FALSE; 345c796909eSBarry Smith #else 346ee8b9732SVaclav Hapla PetscStackCallHDF5(H5Pget_dxpl_mpio,(hdf5->dxpl_id, &mode)); 347ee8b9732SVaclav Hapla *flg = (mode == H5FD_MPIO_COLLECTIVE) ? PETSC_TRUE : PETSC_FALSE; 348c796909eSBarry Smith #endif 349ee8b9732SVaclav Hapla PetscFunctionReturn(0); 350ee8b9732SVaclav Hapla } 351ee8b9732SVaclav Hapla 352ee8b9732SVaclav Hapla /*@ 353ee8b9732SVaclav Hapla PetscViewerHDF5GetCollective - Return flag whether collective MPI-IO transfer mode is used for HDF5 reads and writes. 354ee8b9732SVaclav Hapla 355ee8b9732SVaclav Hapla Not Collective 356ee8b9732SVaclav Hapla 357ee8b9732SVaclav Hapla Input Parameters: 358ee8b9732SVaclav Hapla . viewer - the HDF5 PetscViewer 359ee8b9732SVaclav Hapla 360ee8b9732SVaclav Hapla Output Parameters: 361ee8b9732SVaclav Hapla . flg - the flag 362ee8b9732SVaclav Hapla 363ee8b9732SVaclav Hapla Level: intermediate 364ee8b9732SVaclav Hapla 365ee8b9732SVaclav Hapla Notes: 366c796909eSBarry 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. 367ee8b9732SVaclav Hapla For more details, see PetscViewerHDF5SetCollective(). 368ee8b9732SVaclav Hapla 369ee8b9732SVaclav Hapla .seealso: PetscViewerHDF5SetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open() 370ee8b9732SVaclav Hapla 371ee8b9732SVaclav Hapla @*/ 372ee8b9732SVaclav Hapla PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer viewer,PetscBool *flg) 373ee8b9732SVaclav Hapla { 374ee8b9732SVaclav Hapla PetscErrorCode ierr; 375ee8b9732SVaclav Hapla 376ee8b9732SVaclav Hapla PetscFunctionBegin; 377ee8b9732SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 378534a8f05SLisandro Dalcin PetscValidBoolPointer(flg,2); 379ee8b9732SVaclav Hapla 380ee8b9732SVaclav Hapla ierr = PetscUseMethod(viewer,"PetscViewerHDF5GetCollective_C",(PetscViewer,PetscBool*),(viewer,flg));CHKERRQ(ierr); 381ee8b9732SVaclav Hapla PetscFunctionReturn(0); 382ee8b9732SVaclav Hapla } 383ee8b9732SVaclav Hapla 3849b1403beSVaclav Hapla static PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[]) 3855c6c1daeSBarry Smith { 3865c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 3875c6c1daeSBarry Smith hid_t plist_id; 3885c6c1daeSBarry Smith PetscErrorCode ierr; 3895c6c1daeSBarry Smith 3905c6c1daeSBarry Smith PetscFunctionBegin; 391f683cc58SToby Isaac if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 392f683cc58SToby Isaac if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);} 3935c6c1daeSBarry Smith ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr); 3945c6c1daeSBarry Smith /* Set up file access property list with parallel I/O access */ 395729ab6d9SBarry Smith PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS)); 396c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL) 397d28bfa55SVaclav Hapla PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), MPI_INFO_NULL)); 398c796909eSBarry Smith #endif 3995c6c1daeSBarry Smith /* Create or open the file collectively */ 4005c6c1daeSBarry Smith switch (hdf5->btype) { 4015c6c1daeSBarry Smith case FILE_MODE_READ: 402729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id)); 4035c6c1daeSBarry Smith break; 4045c6c1daeSBarry Smith case FILE_MODE_APPEND: 4057e4fd573SVaclav Hapla case FILE_MODE_UPDATE: 4067e4fd573SVaclav Hapla { 4077e4fd573SVaclav Hapla PetscBool flg; 4087e4fd573SVaclav Hapla ierr = PetscTestFile(hdf5->filename, 'r', &flg);CHKERRQ(ierr); 4097e4fd573SVaclav Hapla if (flg) PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id)); 4107e4fd573SVaclav Hapla else PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_EXCL, H5P_DEFAULT, plist_id)); 4115c6c1daeSBarry Smith break; 4127e4fd573SVaclav Hapla } 4135c6c1daeSBarry Smith case FILE_MODE_WRITE: 414729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id)); 4155c6c1daeSBarry Smith break; 4167e4fd573SVaclav Hapla case FILE_MODE_UNDEFINED: 4177e4fd573SVaclav Hapla SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 4185c6c1daeSBarry Smith default: 4197e4fd573SVaclav Hapla SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP, "Unsupported file mode %s",PetscFileModes[hdf5->btype]); 4205c6c1daeSBarry Smith } 4215c6c1daeSBarry Smith if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name); 422729ab6d9SBarry Smith PetscStackCallHDF5(H5Pclose,(plist_id)); 4235c6c1daeSBarry Smith PetscFunctionReturn(0); 4245c6c1daeSBarry Smith } 4255c6c1daeSBarry Smith 426d1232d7fSToby Isaac static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name) 427d1232d7fSToby Isaac { 428d1232d7fSToby Isaac PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data; 429d1232d7fSToby Isaac 430d1232d7fSToby Isaac PetscFunctionBegin; 431d1232d7fSToby Isaac *name = vhdf5->filename; 432d1232d7fSToby Isaac PetscFunctionReturn(0); 433d1232d7fSToby Isaac } 434d1232d7fSToby Isaac 435b723ab35SVaclav Hapla static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer) 436b723ab35SVaclav Hapla { 4375dc64a97SVaclav Hapla /* 438b723ab35SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 439b723ab35SVaclav Hapla PetscErrorCode ierr; 4405dc64a97SVaclav Hapla */ 441b723ab35SVaclav Hapla 442b723ab35SVaclav Hapla PetscFunctionBegin; 443b723ab35SVaclav Hapla PetscFunctionReturn(0); 444b723ab35SVaclav Hapla } 445b723ab35SVaclav Hapla 4468556b5ebSBarry Smith /*MC 4478556b5ebSBarry Smith PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file 4488556b5ebSBarry Smith 4498556b5ebSBarry Smith 4508556b5ebSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET, 4518556b5ebSBarry Smith PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING, 4528556b5ebSBarry Smith PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB, 4538556b5ebSBarry Smith PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType() 4548556b5ebSBarry Smith 4551b266c99SBarry Smith Level: beginner 4568556b5ebSBarry Smith M*/ 457d1232d7fSToby Isaac 4588cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v) 4595c6c1daeSBarry Smith { 4605c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5; 4615c6c1daeSBarry Smith PetscErrorCode ierr; 4625c6c1daeSBarry Smith 4635c6c1daeSBarry Smith PetscFunctionBegin; 464b00a9115SJed Brown ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr); 4655c6c1daeSBarry Smith 4665c6c1daeSBarry Smith v->data = (void*) hdf5; 4675c6c1daeSBarry Smith v->ops->destroy = PetscViewerDestroy_HDF5; 46882ea9b62SBarry Smith v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5; 469b723ab35SVaclav Hapla v->ops->setup = PetscViewerSetUp_HDF5; 4701b793a25SVaclav Hapla v->ops->view = PetscViewerView_HDF5; 471908793a3SLisandro Dalcin v->ops->flush = NULL; 4727e4fd573SVaclav Hapla hdf5->btype = FILE_MODE_UNDEFINED; 473908793a3SLisandro Dalcin hdf5->filename = NULL; 4745c6c1daeSBarry Smith hdf5->timestep = -1; 4750298fd71SBarry Smith hdf5->groups = NULL; 4765c6c1daeSBarry Smith 4779c5072fbSVaclav Hapla PetscStackCallHDF5Return(hdf5->dxpl_id,H5Pcreate,(H5P_DATASET_XFER)); 4789c5072fbSVaclav Hapla 4790b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr); 480d1232d7fSToby Isaac ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr); 4810b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr); 4820b62783dSVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_HDF5);CHKERRQ(ierr); 48382ea9b62SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr); 4849a0502c6SHåkon Strandenes ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr); 485ee8b9732SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetCollective_C",PetscViewerHDF5SetCollective_HDF5);CHKERRQ(ierr); 486ee8b9732SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetCollective_C",PetscViewerHDF5GetCollective_HDF5);CHKERRQ(ierr); 4875c6c1daeSBarry Smith PetscFunctionReturn(0); 4885c6c1daeSBarry Smith } 4895c6c1daeSBarry Smith 4905c6c1daeSBarry Smith /*@C 4915c6c1daeSBarry Smith PetscViewerHDF5Open - Opens a file for HDF5 input/output. 4925c6c1daeSBarry Smith 493d083f849SBarry Smith Collective 4945c6c1daeSBarry Smith 4955c6c1daeSBarry Smith Input Parameters: 4965c6c1daeSBarry Smith + comm - MPI communicator 4975c6c1daeSBarry Smith . name - name of file 4985c6c1daeSBarry Smith - type - type of file 4995c6c1daeSBarry Smith 5005c6c1daeSBarry Smith Output Parameter: 5015c6c1daeSBarry Smith . hdf5v - PetscViewer for HDF5 input/output to use with the specified file 5025c6c1daeSBarry Smith 50382ea9b62SBarry Smith Options Database: 504a2b725a8SWilliam 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 505a2b725a8SWilliam Gropp - -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal 50682ea9b62SBarry Smith 5075c6c1daeSBarry Smith Level: beginner 5085c6c1daeSBarry Smith 5097e4fd573SVaclav Hapla Notes: 5107e4fd573SVaclav Hapla Reading is always available, regardless of the mode. Available modes are 5117e4fd573SVaclav Hapla + FILE_MODE_READ - open existing HDF5 file for read only access, fail if file does not exist [H5Fopen() with H5F_ACC_RDONLY] 5127e4fd573SVaclav Hapla . FILE_MODE_WRITE - if file exists, fully overwrite it, else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_TRUNC] 5137e4fd573SVaclav 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] 5147e4fd573SVaclav Hapla - FILE_MODE_UPDATE - same as FILE_MODE_APPEND 5157e4fd573SVaclav Hapla 5167e4fd573SVaclav 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. 5177e4fd573SVaclav Hapla 5185c6c1daeSBarry Smith This PetscViewer should be destroyed with PetscViewerDestroy(). 5195c6c1daeSBarry Smith 5205c6c1daeSBarry Smith 5216a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(), 5229a0502c6SHåkon Strandenes PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(), 523a56f64adSBarry Smith MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName() 5245c6c1daeSBarry Smith @*/ 5255c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v) 5265c6c1daeSBarry Smith { 5275c6c1daeSBarry Smith PetscErrorCode ierr; 5285c6c1daeSBarry Smith 5295c6c1daeSBarry Smith PetscFunctionBegin; 5305c6c1daeSBarry Smith ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr); 5315c6c1daeSBarry Smith ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr); 5325c6c1daeSBarry Smith ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr); 5335c6c1daeSBarry Smith ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr); 534b663d2d2SVaclav Hapla ierr = PetscViewerSetFromOptions(*hdf5v);CHKERRQ(ierr); 5355c6c1daeSBarry Smith PetscFunctionReturn(0); 5365c6c1daeSBarry Smith } 5375c6c1daeSBarry Smith 5385c6c1daeSBarry Smith /*@C 5395c6c1daeSBarry Smith PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls 5405c6c1daeSBarry Smith 5415c6c1daeSBarry Smith Not collective 5425c6c1daeSBarry Smith 5435c6c1daeSBarry Smith Input Parameter: 5445c6c1daeSBarry Smith . viewer - the PetscViewer 5455c6c1daeSBarry Smith 5465c6c1daeSBarry Smith Output Parameter: 5475c6c1daeSBarry Smith . file_id - The file id 5485c6c1daeSBarry Smith 5495c6c1daeSBarry Smith Level: intermediate 5505c6c1daeSBarry Smith 5515c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open() 5525c6c1daeSBarry Smith @*/ 5535c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id) 5545c6c1daeSBarry Smith { 5555c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 5565c6c1daeSBarry Smith 5575c6c1daeSBarry Smith PetscFunctionBegin; 5585c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 5595c6c1daeSBarry Smith if (file_id) *file_id = hdf5->file_id; 5605c6c1daeSBarry Smith PetscFunctionReturn(0); 5615c6c1daeSBarry Smith } 5625c6c1daeSBarry Smith 5635c6c1daeSBarry Smith /*@C 5645c6c1daeSBarry Smith PetscViewerHDF5PushGroup - Set the current HDF5 group for output 5655c6c1daeSBarry Smith 5665c6c1daeSBarry Smith Not collective 5675c6c1daeSBarry Smith 5685c6c1daeSBarry Smith Input Parameters: 5695c6c1daeSBarry Smith + viewer - the PetscViewer 5705c6c1daeSBarry Smith - name - The group name 5715c6c1daeSBarry Smith 5725c6c1daeSBarry Smith Level: intermediate 5735c6c1daeSBarry Smith 5742e361344SVaclav Hapla Notes: 5752e361344SVaclav Hapla This is designed to mnemonically resemble the Unix cd command. 5762e361344SVaclav 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. 5772e361344SVaclav Hapla . NULL, empty string, or any sequence of all slashes (e.g. "///") is interpreted as the root group "/". 5782e361344SVaclav Hapla - "." means the current group is pushed again. 5792e361344SVaclav Hapla 5802e361344SVaclav Hapla Example: 5812e361344SVaclav Hapla Suppose the current group is "/a". 5822e361344SVaclav Hapla + If name is NULL, empty string, or a sequence of all slashes (e.g. "///"), then the new group will be "/". 5832e361344SVaclav Hapla . If name is ".", then the new group will be "/a". 5842e361344SVaclav Hapla . If name is "b", then the new group will be "/a/b". 5852e361344SVaclav Hapla - If name is "/b", then the new group will be "/b". 5862e361344SVaclav Hapla 5872e361344SVaclav Hapla Developer Notes: 5882e361344SVaclav Hapla The root group "/" is internally stored as NULL. 589820fc2d1SVaclav Hapla 590874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 5915c6c1daeSBarry Smith @*/ 592be7aa430SVaclav Hapla PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char name[]) 5935c6c1daeSBarry Smith { 5945c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 5957d964842SVaclav Hapla PetscViewerHDF5GroupList *groupNode; 5962e361344SVaclav Hapla size_t i,len; 5972e361344SVaclav Hapla char buf[PETSC_MAX_PATH_LEN]; 5982e361344SVaclav Hapla const char *gname; 5995c6c1daeSBarry Smith PetscErrorCode ierr; 6005c6c1daeSBarry Smith 6015c6c1daeSBarry Smith PetscFunctionBegin; 6025c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 603820fc2d1SVaclav Hapla if (name) PetscValidCharPointer(name,2); 604820fc2d1SVaclav Hapla ierr = PetscStrlen(name, &len);CHKERRQ(ierr); 6052e361344SVaclav Hapla gname = NULL; 6062e361344SVaclav Hapla if (len) { 6072e361344SVaclav Hapla if (len == 1 && name[0] == '.') { 6082e361344SVaclav Hapla /* use current name */ 6092e361344SVaclav Hapla gname = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : NULL; 6102e361344SVaclav Hapla } else if (name[0] == '/') { 6112e361344SVaclav Hapla /* absolute */ 6122e361344SVaclav Hapla for (i=1; i<len; i++) { 6132e361344SVaclav Hapla if (name[i] != '/') { 6142e361344SVaclav Hapla gname = name; 6152e361344SVaclav Hapla break; 6162e361344SVaclav Hapla } 6172e361344SVaclav Hapla } 6182e361344SVaclav Hapla } else { 6192e361344SVaclav Hapla /* relative */ 6202e361344SVaclav Hapla const char *parent = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : ""; 6212e361344SVaclav Hapla ierr = PetscSNPrintf(buf, sizeof(buf), "%s/%s", parent, name);CHKERRQ(ierr); 6222e361344SVaclav Hapla gname = buf; 6232e361344SVaclav Hapla } 6242e361344SVaclav Hapla } 62595dccacaSBarry Smith ierr = PetscNew(&groupNode);CHKERRQ(ierr); 6262e361344SVaclav Hapla ierr = PetscStrallocpy(gname, (char**) &groupNode->name);CHKERRQ(ierr); 6275c6c1daeSBarry Smith groupNode->next = hdf5->groups; 6285c6c1daeSBarry Smith hdf5->groups = groupNode; 6295c6c1daeSBarry Smith PetscFunctionReturn(0); 6305c6c1daeSBarry Smith } 6315c6c1daeSBarry Smith 6323ef9c667SSatish Balay /*@ 6335c6c1daeSBarry Smith PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value 6345c6c1daeSBarry Smith 6355c6c1daeSBarry Smith Not collective 6365c6c1daeSBarry Smith 6375c6c1daeSBarry Smith Input Parameter: 6385c6c1daeSBarry Smith . viewer - the PetscViewer 6395c6c1daeSBarry Smith 6405c6c1daeSBarry Smith Level: intermediate 6415c6c1daeSBarry Smith 642874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 6435c6c1daeSBarry Smith @*/ 6445c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer) 6455c6c1daeSBarry Smith { 6465c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 6477d964842SVaclav Hapla PetscViewerHDF5GroupList *groupNode; 6485c6c1daeSBarry Smith PetscErrorCode ierr; 6495c6c1daeSBarry Smith 6505c6c1daeSBarry Smith PetscFunctionBegin; 6515c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 65282f516ccSBarry Smith if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop"); 6535c6c1daeSBarry Smith groupNode = hdf5->groups; 6545c6c1daeSBarry Smith hdf5->groups = hdf5->groups->next; 6555c6c1daeSBarry Smith ierr = PetscFree(groupNode->name);CHKERRQ(ierr); 6565c6c1daeSBarry Smith ierr = PetscFree(groupNode);CHKERRQ(ierr); 6575c6c1daeSBarry Smith PetscFunctionReturn(0); 6585c6c1daeSBarry Smith } 6595c6c1daeSBarry Smith 6605c6c1daeSBarry Smith /*@C 661874270d9SVaclav Hapla PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup(). 662874270d9SVaclav Hapla If none has been assigned, returns NULL. 6635c6c1daeSBarry Smith 6645c6c1daeSBarry Smith Not collective 6655c6c1daeSBarry Smith 6665c6c1daeSBarry Smith Input Parameter: 6675c6c1daeSBarry Smith . viewer - the PetscViewer 6685c6c1daeSBarry Smith 6695c6c1daeSBarry Smith Output Parameter: 6705c6c1daeSBarry Smith . name - The group name 6715c6c1daeSBarry Smith 6725c6c1daeSBarry Smith Level: intermediate 6735c6c1daeSBarry Smith 674874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup() 6755c6c1daeSBarry Smith @*/ 676be7aa430SVaclav Hapla PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char *name[]) 6775c6c1daeSBarry Smith { 6785c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 6795c6c1daeSBarry Smith 6805c6c1daeSBarry Smith PetscFunctionBegin; 6815c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 682c959eef4SJed Brown PetscValidPointer(name,2); 683a297a907SKarl Rupp if (hdf5->groups) *name = hdf5->groups->name; 6840298fd71SBarry Smith else *name = NULL; 6855c6c1daeSBarry Smith PetscFunctionReturn(0); 6865c6c1daeSBarry Smith } 6875c6c1daeSBarry Smith 6888c557b5aSVaclav Hapla /*@ 689874270d9SVaclav Hapla PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(), 690874270d9SVaclav Hapla and return this group's ID and file ID. 691874270d9SVaclav Hapla If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID. 692874270d9SVaclav Hapla 693874270d9SVaclav Hapla Not collective 694874270d9SVaclav Hapla 695874270d9SVaclav Hapla Input Parameter: 696874270d9SVaclav Hapla . viewer - the PetscViewer 697874270d9SVaclav Hapla 698874270d9SVaclav Hapla Output Parameter: 699874270d9SVaclav Hapla + fileId - The HDF5 file ID 700874270d9SVaclav Hapla - groupId - The HDF5 group ID 701874270d9SVaclav Hapla 702e74616beSVaclav Hapla Notes: 703e74616beSVaclav Hapla If the viewer is writable, the group is created if it doesn't exist yet. 704e74616beSVaclav Hapla 705874270d9SVaclav Hapla Level: intermediate 706874270d9SVaclav Hapla 707874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 708874270d9SVaclav Hapla @*/ 70954dbf706SVaclav Hapla PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId) 71054dbf706SVaclav Hapla { 71190e03537SVaclav Hapla hid_t file_id; 71290e03537SVaclav Hapla H5O_type_t type; 71354dbf706SVaclav Hapla const char *groupName = NULL; 714e74616beSVaclav Hapla PetscBool create; 71554dbf706SVaclav Hapla PetscErrorCode ierr; 71654dbf706SVaclav Hapla 71754dbf706SVaclav Hapla PetscFunctionBegin; 718e74616beSVaclav Hapla ierr = PetscViewerWritable(viewer, &create);CHKERRQ(ierr); 71954dbf706SVaclav Hapla ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr); 72054dbf706SVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr); 721e74616beSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, groupName, create, NULL, &type);CHKERRQ(ierr); 72290e03537SVaclav 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); 72390e03537SVaclav Hapla PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT)); 72454dbf706SVaclav Hapla *fileId = file_id; 72554dbf706SVaclav Hapla PetscFunctionReturn(0); 72654dbf706SVaclav Hapla } 72754dbf706SVaclav Hapla 7283ef9c667SSatish Balay /*@ 7295c6c1daeSBarry Smith PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time. 7305c6c1daeSBarry Smith 7315c6c1daeSBarry Smith Not collective 7325c6c1daeSBarry Smith 7335c6c1daeSBarry Smith Input Parameter: 7345c6c1daeSBarry Smith . viewer - the PetscViewer 7355c6c1daeSBarry Smith 7365c6c1daeSBarry Smith Level: intermediate 7375c6c1daeSBarry Smith 7385c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep() 7395c6c1daeSBarry Smith @*/ 7405c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer) 7415c6c1daeSBarry Smith { 7425c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 7435c6c1daeSBarry Smith 7445c6c1daeSBarry Smith PetscFunctionBegin; 7455c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 7465c6c1daeSBarry Smith ++hdf5->timestep; 7475c6c1daeSBarry Smith PetscFunctionReturn(0); 7485c6c1daeSBarry Smith } 7495c6c1daeSBarry Smith 7503ef9c667SSatish Balay /*@ 7515c6c1daeSBarry Smith PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep 7525c6c1daeSBarry Smith of -1 disables blocking with timesteps. 7535c6c1daeSBarry Smith 7545c6c1daeSBarry Smith Not collective 7555c6c1daeSBarry Smith 7565c6c1daeSBarry Smith Input Parameters: 7575c6c1daeSBarry Smith + viewer - the PetscViewer 7585c6c1daeSBarry Smith - timestep - The timestep number 7595c6c1daeSBarry Smith 7605c6c1daeSBarry Smith Level: intermediate 7615c6c1daeSBarry Smith 7625c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 7635c6c1daeSBarry Smith @*/ 7645c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep) 7655c6c1daeSBarry Smith { 7665c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 7675c6c1daeSBarry Smith 7685c6c1daeSBarry Smith PetscFunctionBegin; 7695c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 7705c6c1daeSBarry Smith hdf5->timestep = timestep; 7715c6c1daeSBarry Smith PetscFunctionReturn(0); 7725c6c1daeSBarry Smith } 7735c6c1daeSBarry Smith 7743ef9c667SSatish Balay /*@ 7755c6c1daeSBarry Smith PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time. 7765c6c1daeSBarry Smith 7775c6c1daeSBarry Smith Not collective 7785c6c1daeSBarry Smith 7795c6c1daeSBarry Smith Input Parameter: 7805c6c1daeSBarry Smith . viewer - the PetscViewer 7815c6c1daeSBarry Smith 7825c6c1daeSBarry Smith Output Parameter: 7835c6c1daeSBarry Smith . timestep - The timestep number 7845c6c1daeSBarry Smith 7855c6c1daeSBarry Smith Level: intermediate 7865c6c1daeSBarry Smith 7875c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep() 7885c6c1daeSBarry Smith @*/ 7895c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep) 7905c6c1daeSBarry Smith { 7915c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 7925c6c1daeSBarry Smith 7935c6c1daeSBarry Smith PetscFunctionBegin; 7945c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 7955c6c1daeSBarry Smith PetscValidPointer(timestep,2); 7965c6c1daeSBarry Smith *timestep = hdf5->timestep; 7975c6c1daeSBarry Smith PetscFunctionReturn(0); 7985c6c1daeSBarry Smith } 7995c6c1daeSBarry Smith 80036bce6e8SMatthew G. Knepley /*@C 80136bce6e8SMatthew G. Knepley PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name. 80236bce6e8SMatthew G. Knepley 80336bce6e8SMatthew G. Knepley Not collective 80436bce6e8SMatthew G. Knepley 80536bce6e8SMatthew G. Knepley Input Parameter: 80636bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 80736bce6e8SMatthew G. Knepley 80836bce6e8SMatthew G. Knepley Output Parameter: 80936bce6e8SMatthew G. Knepley . mtype - the MPI datatype (for example MPI_DOUBLE, ...) 81036bce6e8SMatthew G. Knepley 81136bce6e8SMatthew G. Knepley Level: advanced 81236bce6e8SMatthew G. Knepley 81336bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 81436bce6e8SMatthew G. Knepley @*/ 81536bce6e8SMatthew G. Knepley PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype) 81636bce6e8SMatthew G. Knepley { 81736bce6e8SMatthew G. Knepley PetscFunctionBegin; 81836bce6e8SMatthew G. Knepley if (ptype == PETSC_INT) 81936bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 82036bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_LLONG; 82136bce6e8SMatthew G. Knepley #else 82236bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_INT; 82336bce6e8SMatthew G. Knepley #endif 82436bce6e8SMatthew G. Knepley else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE; 82536bce6e8SMatthew G. Knepley else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG; 82636bce6e8SMatthew G. Knepley else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT; 827c9450970SBarry Smith else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_INT; 828de3ba262SVaclav Hapla else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_INT; 82936bce6e8SMatthew G. Knepley else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT; 83036bce6e8SMatthew G. Knepley else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR; 83136bce6e8SMatthew G. Knepley else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR; 8327e97c476SMatthew G. Knepley else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1); 83336bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype"); 83436bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 83536bce6e8SMatthew G. Knepley } 83636bce6e8SMatthew G. Knepley 83736bce6e8SMatthew G. Knepley /*@C 83836bce6e8SMatthew G. Knepley PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name 83936bce6e8SMatthew G. Knepley 84036bce6e8SMatthew G. Knepley Not collective 84136bce6e8SMatthew G. Knepley 84236bce6e8SMatthew G. Knepley Input Parameter: 84336bce6e8SMatthew G. Knepley . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...) 84436bce6e8SMatthew G. Knepley 84536bce6e8SMatthew G. Knepley Output Parameter: 84636bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 84736bce6e8SMatthew G. Knepley 84836bce6e8SMatthew G. Knepley Level: advanced 84936bce6e8SMatthew G. Knepley 85036bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 85136bce6e8SMatthew G. Knepley @*/ 85236bce6e8SMatthew G. Knepley PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype) 85336bce6e8SMatthew G. Knepley { 85436bce6e8SMatthew G. Knepley PetscFunctionBegin; 85536bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 85636bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG; 85736bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT; 85836bce6e8SMatthew G. Knepley #else 85936bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT; 86036bce6e8SMatthew G. Knepley #endif 86136bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE; 86236bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG; 86336bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT; 86436bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT; 86536bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR; 86636bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR; 8677e97c476SMatthew G. Knepley else if (htype == H5T_C_S1) *ptype = PETSC_STRING; 86836bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype"); 86936bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 87036bce6e8SMatthew G. Knepley } 87136bce6e8SMatthew G. Knepley 872df863907SAlex Fikl /*@C 873b8ee85a1SVaclav Hapla PetscViewerHDF5WriteAttribute - Write an attribute 87436bce6e8SMatthew G. Knepley 87536bce6e8SMatthew G. Knepley Input Parameters: 87636bce6e8SMatthew G. Knepley + viewer - The HDF5 viewer 877*4302210dSVaclav Hapla . parent - The parent dataset/group name 87836bce6e8SMatthew G. Knepley . name - The attribute name 87936bce6e8SMatthew G. Knepley . datatype - The attribute type 88036bce6e8SMatthew G. Knepley - value - The attribute value 88136bce6e8SMatthew G. Knepley 88236bce6e8SMatthew G. Knepley Level: advanced 88336bce6e8SMatthew G. Knepley 884*4302210dSVaclav Hapla Notes: 885*4302210dSVaclav Hapla If parent starts with '/', it is taken as an absolute path overriding currently pushed group, else the dataset/group is relative to the current pushed group. NULL means the current pushed group. 886*4302210dSVaclav Hapla 887577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 88836bce6e8SMatthew G. Knepley @*/ 889*4302210dSVaclav Hapla PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value) 89036bce6e8SMatthew G. Knepley { 891*4302210dSVaclav Hapla char *parentAbsPath; 89260568592SMatthew G. Knepley hid_t h5, dataspace, obj, attribute, dtype; 893080f144cSVaclav Hapla PetscBool has; 89436bce6e8SMatthew G. Knepley PetscErrorCode ierr; 89536bce6e8SMatthew G. Knepley 89636bce6e8SMatthew G. Knepley PetscFunctionBegin; 8975cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 898*4302210dSVaclav Hapla if (parent) PetscValidCharPointer(parent, 2); 899c57a1a9aSVaclav Hapla PetscValidCharPointer(name, 3); 900*4302210dSVaclav Hapla PetscValidLogicalCollectiveEnum(viewer,datatype,4); 901b7e81f0fSVaclav Hapla PetscValidPointer(value, 5); 902*4302210dSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath);CHKERRQ(ierr); 903*4302210dSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_TRUE, NULL, NULL);CHKERRQ(ierr); 904*4302210dSVaclav Hapla ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has);CHKERRQ(ierr); 90536bce6e8SMatthew G. Knepley ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 9067e97c476SMatthew G. Knepley if (datatype == PETSC_STRING) { 9077e97c476SMatthew G. Knepley size_t len; 9087e97c476SMatthew G. Knepley ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr); 909729ab6d9SBarry Smith PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 9107e97c476SMatthew G. Knepley } 91136bce6e8SMatthew G. Knepley ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 912729ab6d9SBarry Smith PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR)); 913*4302210dSVaclav Hapla PetscStackCallHDF5Return(obj,H5Oopen,(h5, parentAbsPath, H5P_DEFAULT)); 914080f144cSVaclav Hapla if (has) { 915080f144cSVaclav Hapla PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 916080f144cSVaclav Hapla } else { 91760568592SMatthew G. Knepley PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT)); 918080f144cSVaclav Hapla } 919729ab6d9SBarry Smith PetscStackCallHDF5(H5Awrite,(attribute, dtype, value)); 920729ab6d9SBarry Smith if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype)); 921729ab6d9SBarry Smith PetscStackCallHDF5(H5Aclose,(attribute)); 92260568592SMatthew G. Knepley PetscStackCallHDF5(H5Oclose,(obj)); 923729ab6d9SBarry Smith PetscStackCallHDF5(H5Sclose,(dataspace)); 924*4302210dSVaclav Hapla ierr = PetscFree(parentAbsPath);CHKERRQ(ierr); 92536bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 92636bce6e8SMatthew G. Knepley } 92736bce6e8SMatthew G. Knepley 928df863907SAlex Fikl /*@C 929577e0e04SVaclav Hapla PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name 930577e0e04SVaclav Hapla 931577e0e04SVaclav Hapla Input Parameters: 932577e0e04SVaclav Hapla + viewer - The HDF5 viewer 933577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 934577e0e04SVaclav Hapla . name - The attribute name 935577e0e04SVaclav Hapla . datatype - The attribute type 936577e0e04SVaclav Hapla - value - The attribute value 937577e0e04SVaclav Hapla 938577e0e04SVaclav Hapla Notes: 939*4302210dSVaclav Hapla This fails if currentGroup/objectName doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 940577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 941577e0e04SVaclav Hapla 942577e0e04SVaclav Hapla Level: advanced 943577e0e04SVaclav Hapla 944577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 945577e0e04SVaclav Hapla @*/ 946577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value) 947577e0e04SVaclav Hapla { 948577e0e04SVaclav Hapla PetscErrorCode ierr; 949577e0e04SVaclav Hapla 950577e0e04SVaclav Hapla PetscFunctionBegin; 951577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 952577e0e04SVaclav Hapla PetscValidHeader(obj,2); 953577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 954b7e81f0fSVaclav Hapla PetscValidPointer(value,5); 955577e0e04SVaclav Hapla ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 956577e0e04SVaclav Hapla ierr = PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr); 957577e0e04SVaclav Hapla PetscFunctionReturn(0); 958577e0e04SVaclav Hapla } 959577e0e04SVaclav Hapla 960577e0e04SVaclav Hapla /*@C 961b8ee85a1SVaclav Hapla PetscViewerHDF5ReadAttribute - Read an attribute 962ad0c61feSMatthew G. Knepley 963ad0c61feSMatthew G. Knepley Input Parameters: 964ad0c61feSMatthew G. Knepley + viewer - The HDF5 viewer 965*4302210dSVaclav Hapla . parent - The parent dataset/group name 966ad0c61feSMatthew G. Knepley . name - The attribute name 967ad0c61feSMatthew G. Knepley - datatype - The attribute type 968ad0c61feSMatthew G. Knepley 969ad0c61feSMatthew G. Knepley Output Parameter: 970ad0c61feSMatthew G. Knepley . value - The attribute value 971ad0c61feSMatthew G. Knepley 972708d4cb3SBarry Smith Notes: If the datatype is PETSC_STRING one must PetscFree() the obtained value when it is no longer needed. 973708d4cb3SBarry Smith 974ad0c61feSMatthew G. Knepley Level: advanced 975ad0c61feSMatthew G. Knepley 976*4302210dSVaclav Hapla Notes: 977*4302210dSVaclav Hapla If parent starts with '/', it is taken as an absolute path overriding currently pushed group, else the dataset/group is relative to the current pushed group. NULL means the current pushed group. 978*4302210dSVaclav Hapla 979577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 980ad0c61feSMatthew G. Knepley @*/ 981*4302210dSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, void *value) 982ad0c61feSMatthew G. Knepley { 983*4302210dSVaclav Hapla char *parentAbsPath; 984f0b58479SMatthew G. Knepley hid_t h5, obj, attribute, atype, dtype; 985969748fdSVaclav Hapla PetscBool has; 986ad0c61feSMatthew G. Knepley PetscErrorCode ierr; 987ad0c61feSMatthew G. Knepley 988ad0c61feSMatthew G. Knepley PetscFunctionBegin; 9895cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 990*4302210dSVaclav Hapla if (parent) PetscValidCharPointer(parent, 2); 991c57a1a9aSVaclav Hapla PetscValidCharPointer(name, 3); 992b7e81f0fSVaclav Hapla PetscValidPointer(value, 5); 993*4302210dSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath);CHKERRQ(ierr); 994*4302210dSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, &has, NULL);CHKERRQ(ierr); 995*4302210dSVaclav Hapla if (has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has);CHKERRQ(ierr);} 996*4302210dSVaclav Hapla if (!has) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist", parentAbsPath, name); 997ad0c61feSMatthew G. Knepley ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 998ad0c61feSMatthew G. Knepley ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 999*4302210dSVaclav Hapla PetscStackCallHDF5Return(obj,H5Oopen,(h5, parentAbsPath, H5P_DEFAULT)); 100060568592SMatthew G. Knepley PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 1001f0b58479SMatthew G. Knepley if (datatype == PETSC_STRING) { 1002f0b58479SMatthew G. Knepley size_t len; 1003e90c5075SMarek Pecha PetscStackCallHDF5Return(atype,H5Aget_type,(attribute)); 100415b861d2SVaclav Hapla PetscStackCall("H5Tget_size",len = H5Tget_size(atype)); 1005708d4cb3SBarry Smith ierr = PetscMalloc((len+1) * sizeof(char), value);CHKERRQ(ierr); 1006708d4cb3SBarry Smith PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 1007708d4cb3SBarry Smith PetscStackCallHDF5(H5Aread,(attribute, dtype, *(char**)value)); 1008708d4cb3SBarry Smith } else { 100970efba33SBarry Smith PetscStackCallHDF5(H5Aread,(attribute, dtype, value)); 1010708d4cb3SBarry Smith } 1011729ab6d9SBarry Smith PetscStackCallHDF5(H5Aclose,(attribute)); 1012e90c5075SMarek Pecha /* H5Oclose can be used to close groups, datasets, or committed datatypes */ 1013e90c5075SMarek Pecha PetscStackCallHDF5(H5Oclose,(obj)); 1014*4302210dSVaclav Hapla ierr = PetscFree(parentAbsPath);CHKERRQ(ierr); 1015ad0c61feSMatthew G. Knepley PetscFunctionReturn(0); 1016ad0c61feSMatthew G. Knepley } 1017ad0c61feSMatthew G. Knepley 1018577e0e04SVaclav Hapla /*@C 1019577e0e04SVaclav Hapla PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name 1020577e0e04SVaclav Hapla 1021577e0e04SVaclav Hapla Input Parameters: 1022577e0e04SVaclav Hapla + viewer - The HDF5 viewer 1023577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1024577e0e04SVaclav Hapla . name - The attribute name 1025577e0e04SVaclav Hapla - datatype - The attribute type 1026577e0e04SVaclav Hapla 1027577e0e04SVaclav Hapla Output Parameter: 1028577e0e04SVaclav Hapla . value - The attribute value 1029577e0e04SVaclav Hapla 1030577e0e04SVaclav Hapla Notes: 1031577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1032577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 1033577e0e04SVaclav Hapla 1034577e0e04SVaclav Hapla Level: advanced 1035577e0e04SVaclav Hapla 1036577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute() PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1037577e0e04SVaclav Hapla @*/ 1038577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *value) 1039577e0e04SVaclav Hapla { 1040577e0e04SVaclav Hapla PetscErrorCode ierr; 1041577e0e04SVaclav Hapla 1042577e0e04SVaclav Hapla PetscFunctionBegin; 1043577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1044577e0e04SVaclav Hapla PetscValidHeader(obj,2); 1045577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 1046b7e81f0fSVaclav Hapla PetscValidPointer(value, 5); 1047577e0e04SVaclav Hapla ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 1048577e0e04SVaclav Hapla ierr = PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr); 1049577e0e04SVaclav Hapla PetscFunctionReturn(0); 1050577e0e04SVaclav Hapla } 1051577e0e04SVaclav Hapla 1052a75c3fd4SVaclav Hapla PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_) 1053a75c3fd4SVaclav Hapla { 1054a75c3fd4SVaclav Hapla htri_t exists; 1055a75c3fd4SVaclav Hapla hid_t group; 1056a75c3fd4SVaclav Hapla 1057a75c3fd4SVaclav Hapla PetscFunctionBegin; 1058a75c3fd4SVaclav Hapla PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT)); 1059a75c3fd4SVaclav Hapla if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT)); 1060a75c3fd4SVaclav Hapla if (!exists && createGroup) { 1061a75c3fd4SVaclav Hapla PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT)); 1062a75c3fd4SVaclav Hapla PetscStackCallHDF5(H5Gclose,(group)); 1063a75c3fd4SVaclav Hapla exists = PETSC_TRUE; 1064a75c3fd4SVaclav Hapla } 1065a75c3fd4SVaclav Hapla *exists_ = (PetscBool) exists; 1066a75c3fd4SVaclav Hapla PetscFunctionReturn(0); 1067a75c3fd4SVaclav Hapla } 1068a75c3fd4SVaclav Hapla 1069bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype) 10705cdeb108SMatthew G. Knepley { 107190e03537SVaclav Hapla const char rootGroupName[] = "/"; 10725cdeb108SMatthew G. Knepley hid_t h5; 1073e5bf9ebcSVaclav Hapla PetscBool exists=PETSC_FALSE; 107415b861d2SVaclav Hapla PetscInt i; 107515b861d2SVaclav Hapla int n; 107685688ae2SVaclav Hapla char **hierarchy; 107785688ae2SVaclav Hapla char buf[PETSC_MAX_PATH_LEN]=""; 10785cdeb108SMatthew G. Knepley PetscErrorCode ierr; 10795cdeb108SMatthew G. Knepley 10805cdeb108SMatthew G. Knepley PetscFunctionBegin; 10815cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 108290e03537SVaclav Hapla if (name) PetscValidCharPointer(name, 2); 108390e03537SVaclav Hapla else name = rootGroupName; 1084ccd66a83SVaclav Hapla if (has) { 108556cc0592SVaclav Hapla PetscValidIntPointer(has, 3); 10865cdeb108SMatthew G. Knepley *has = PETSC_FALSE; 1087ccd66a83SVaclav Hapla } 1088ccd66a83SVaclav Hapla if (otype) { 1089ccd66a83SVaclav Hapla PetscValidIntPointer(otype, 4); 109056cc0592SVaclav Hapla *otype = H5O_TYPE_UNKNOWN; 1091ccd66a83SVaclav Hapla } 1092c0cd0301SJed Brown ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 109385688ae2SVaclav Hapla 109485688ae2SVaclav Hapla /* 109585688ae2SVaclav Hapla Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing. 109685688ae2SVaclav Hapla Hence, each of them needs to be tested separately: 109785688ae2SVaclav Hapla 1) whether it's a valid link 109885688ae2SVaclav Hapla 2) whether this link resolves to an object 109985688ae2SVaclav Hapla See H5Oexists_by_name() documentation. 110085688ae2SVaclav Hapla */ 110185688ae2SVaclav Hapla ierr = PetscStrToArray(name,'/',&n,&hierarchy);CHKERRQ(ierr); 110285688ae2SVaclav Hapla if (!n) { 110385688ae2SVaclav Hapla /* Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */ 1104ccd66a83SVaclav Hapla if (has) *has = PETSC_TRUE; 1105ccd66a83SVaclav Hapla if (otype) *otype = H5O_TYPE_GROUP; 110685688ae2SVaclav Hapla ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 110785688ae2SVaclav Hapla PetscFunctionReturn(0); 110885688ae2SVaclav Hapla } 110985688ae2SVaclav Hapla for (i=0; i<n; i++) { 111085688ae2SVaclav Hapla ierr = PetscStrcat(buf,"/");CHKERRQ(ierr); 111185688ae2SVaclav Hapla ierr = PetscStrcat(buf,hierarchy[i]);CHKERRQ(ierr); 1112a75c3fd4SVaclav Hapla ierr = PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);CHKERRQ(ierr); 1113a75c3fd4SVaclav Hapla if (!exists) break; 111485688ae2SVaclav Hapla } 111585688ae2SVaclav Hapla ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 111685688ae2SVaclav Hapla 111785688ae2SVaclav Hapla /* If the object exists, get its type */ 1118ccd66a83SVaclav Hapla if (exists && otype) { 11195cdeb108SMatthew G. Knepley H5O_info_t info; 11205cdeb108SMatthew G. Knepley 112176276748SVaclav Hapla /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */ 112204633f7fSVaclav Hapla PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT)); 112356cc0592SVaclav Hapla *otype = info.type; 11245cdeb108SMatthew G. Knepley } 1125ccd66a83SVaclav Hapla if (has) *has = exists; 11265cdeb108SMatthew G. Knepley PetscFunctionReturn(0); 11275cdeb108SMatthew G. Knepley } 11285cdeb108SMatthew G. Knepley 1129*4302210dSVaclav Hapla /*@C 11300a9f38efSVaclav Hapla PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file 11310a9f38efSVaclav Hapla 11320a9f38efSVaclav Hapla Input Parameters: 1133*4302210dSVaclav Hapla + viewer - The HDF5 viewer 1134*4302210dSVaclav Hapla - path - The group path 11350a9f38efSVaclav Hapla 11360a9f38efSVaclav Hapla Output Parameter: 11370a9f38efSVaclav Hapla . has - Flag for group existence 11380a9f38efSVaclav Hapla 11390a9f38efSVaclav Hapla Level: advanced 11400a9f38efSVaclav Hapla 1141*4302210dSVaclav Hapla Notes: 1142*4302210dSVaclav Hapla If path starts with '/', it is taken as an absolute path overriding currently pushed group, else the dataset/group is relative to the current pushed group. NULL means the current pushed group. 1143*4302210dSVaclav Hapla If the path exists but is not a group, PETSC_FALSE is returned. 1144*4302210dSVaclav Hapla 1145*4302210dSVaclav Hapla .seealso: PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5GetGroup(), PetscViewerHDF5OpenGroup() 11460a9f38efSVaclav Hapla @*/ 1147*4302210dSVaclav Hapla PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, const char path[], PetscBool *has) 11480a9f38efSVaclav Hapla { 11490a9f38efSVaclav Hapla H5O_type_t type; 1150*4302210dSVaclav Hapla char *abspath; 11510a9f38efSVaclav Hapla PetscErrorCode ierr; 11520a9f38efSVaclav Hapla 11530a9f38efSVaclav Hapla PetscFunctionBegin; 11540a9f38efSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1155*4302210dSVaclav Hapla if (path) PetscValidCharPointer(path,2); 1156*4302210dSVaclav Hapla PetscValidBoolPointer(has,3); 1157*4302210dSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, path, &abspath);CHKERRQ(ierr); 1158*4302210dSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type);CHKERRQ(ierr); 1159*4302210dSVaclav Hapla *has = (PetscBool)(type == H5O_TYPE_GROUP); 1160*4302210dSVaclav Hapla ierr = PetscFree(abspath);CHKERRQ(ierr); 11610a9f38efSVaclav Hapla PetscFunctionReturn(0); 11620a9f38efSVaclav Hapla } 11630a9f38efSVaclav Hapla 11640a9f38efSVaclav Hapla /*@ 1165e3f143f7SVaclav Hapla PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group 1166ecce1506SVaclav Hapla 1167ecce1506SVaclav Hapla Input Parameters: 1168ecce1506SVaclav Hapla + viewer - The HDF5 viewer 1169ecce1506SVaclav Hapla - obj - The named object 1170ecce1506SVaclav Hapla 1171ecce1506SVaclav Hapla Output Parameter: 1172ecce1506SVaclav Hapla . has - Flag for dataset existence; PETSC_FALSE for unnamed object 1173ecce1506SVaclav Hapla 1174e3f143f7SVaclav Hapla Notes: 1175*4302210dSVaclav Hapla If the object is unnamed, has is set to PETSC_FALSE. 1176*4302210dSVaclav Hapla If the path currentGroup/objectName exists but is not a dataset, has is set to PETSC_FALSE as well. 1177e3f143f7SVaclav Hapla 1178ecce1506SVaclav Hapla Level: advanced 1179ecce1506SVaclav Hapla 1180e3f143f7SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5GetGroup() 1181ecce1506SVaclav Hapla @*/ 1182ecce1506SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has) 1183ecce1506SVaclav Hapla { 118456cc0592SVaclav Hapla H5O_type_t type; 11856c132bc1SVaclav Hapla char *path; 1186ecce1506SVaclav Hapla PetscErrorCode ierr; 1187ecce1506SVaclav Hapla 1188ecce1506SVaclav Hapla PetscFunctionBegin; 1189c57a1a9aSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1190c57a1a9aSVaclav Hapla PetscValidHeader(obj,2); 1191*4302210dSVaclav Hapla PetscValidBoolPointer(has,3); 1192ecce1506SVaclav Hapla *has = PETSC_FALSE; 1193e3f143f7SVaclav Hapla if (!obj->name) PetscFunctionReturn(0); 11946c132bc1SVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, obj->name, &path);CHKERRQ(ierr); 11956c132bc1SVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, path, PETSC_FALSE, has, &type);CHKERRQ(ierr); 119656cc0592SVaclav Hapla *has = (type == H5O_TYPE_DATASET) ? PETSC_TRUE : PETSC_FALSE; 11976c132bc1SVaclav Hapla ierr = PetscFree(path);CHKERRQ(ierr); 1198ecce1506SVaclav Hapla PetscFunctionReturn(0); 1199ecce1506SVaclav Hapla } 1200ecce1506SVaclav Hapla 1201df863907SAlex Fikl /*@C 1202ece83b6cSVaclav Hapla PetscViewerHDF5HasAttribute - Check whether an attribute exists 1203e7bdbf8eSMatthew G. Knepley 1204e7bdbf8eSMatthew G. Knepley Input Parameters: 1205e7bdbf8eSMatthew G. Knepley + viewer - The HDF5 viewer 1206*4302210dSVaclav Hapla . parent - The parent dataset/group name 1207e7bdbf8eSMatthew G. Knepley - name - The attribute name 1208e7bdbf8eSMatthew G. Knepley 1209e7bdbf8eSMatthew G. Knepley Output Parameter: 1210e7bdbf8eSMatthew G. Knepley . has - Flag for attribute existence 1211e7bdbf8eSMatthew G. Knepley 1212e7bdbf8eSMatthew G. Knepley Level: advanced 1213e7bdbf8eSMatthew G. Knepley 1214*4302210dSVaclav Hapla Notes: 1215*4302210dSVaclav Hapla If parent starts with '/', it is taken as an absolute path overriding currently pushed group, else the dataset/group is relative to the current pushed group. NULL means the current pushed group. 1216*4302210dSVaclav Hapla 1217577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1218e7bdbf8eSMatthew G. Knepley @*/ 1219*4302210dSVaclav Hapla PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 1220e7bdbf8eSMatthew G. Knepley { 1221*4302210dSVaclav Hapla char *parentAbsPath; 1222e7bdbf8eSMatthew G. Knepley PetscErrorCode ierr; 1223e7bdbf8eSMatthew G. Knepley 1224e7bdbf8eSMatthew G. Knepley PetscFunctionBegin; 12255cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1226*4302210dSVaclav Hapla if (parent) PetscValidCharPointer(parent,2); 1227c57a1a9aSVaclav Hapla PetscValidCharPointer(name,3); 1228*4302210dSVaclav Hapla PetscValidBoolPointer(has,4); 1229*4302210dSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath);CHKERRQ(ierr); 1230*4302210dSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, has, NULL);CHKERRQ(ierr); 1231*4302210dSVaclav Hapla if (*has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, has);CHKERRQ(ierr);} 1232*4302210dSVaclav Hapla ierr = PetscFree(parentAbsPath);CHKERRQ(ierr); 123306db490cSVaclav Hapla PetscFunctionReturn(0); 123406db490cSVaclav Hapla } 123506db490cSVaclav Hapla 1236577e0e04SVaclav Hapla /*@C 1237577e0e04SVaclav Hapla PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name 1238577e0e04SVaclav Hapla 1239577e0e04SVaclav Hapla Input Parameters: 1240577e0e04SVaclav Hapla + viewer - The HDF5 viewer 1241577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1242577e0e04SVaclav Hapla - name - The attribute name 1243577e0e04SVaclav Hapla 1244577e0e04SVaclav Hapla Output Parameter: 1245577e0e04SVaclav Hapla . has - Flag for attribute existence 1246577e0e04SVaclav Hapla 1247577e0e04SVaclav Hapla Notes: 1248577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1249577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 1250577e0e04SVaclav Hapla 1251577e0e04SVaclav Hapla Level: advanced 1252577e0e04SVaclav Hapla 1253577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1254577e0e04SVaclav Hapla @*/ 1255577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has) 1256577e0e04SVaclav Hapla { 1257577e0e04SVaclav Hapla PetscErrorCode ierr; 1258577e0e04SVaclav Hapla 1259577e0e04SVaclav Hapla PetscFunctionBegin; 1260577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1261577e0e04SVaclav Hapla PetscValidHeader(obj,2); 1262577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 1263*4302210dSVaclav Hapla PetscValidBoolPointer(has,4); 1264577e0e04SVaclav Hapla ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 1265577e0e04SVaclav Hapla ierr = PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);CHKERRQ(ierr); 1266577e0e04SVaclav Hapla PetscFunctionReturn(0); 1267577e0e04SVaclav Hapla } 1268577e0e04SVaclav Hapla 126906db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 127006db490cSVaclav Hapla { 127106db490cSVaclav Hapla hid_t h5; 127206db490cSVaclav Hapla htri_t hhas; 127306db490cSVaclav Hapla PetscErrorCode ierr; 127406db490cSVaclav Hapla 127506db490cSVaclav Hapla PetscFunctionBegin; 127606db490cSVaclav Hapla ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 12772f936e54SVaclav Hapla PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT)); 12785cdeb108SMatthew G. Knepley *has = hhas ? PETSC_TRUE : PETSC_FALSE; 1279e7bdbf8eSMatthew G. Knepley PetscFunctionReturn(0); 1280e7bdbf8eSMatthew G. Knepley } 1281e7bdbf8eSMatthew G. Knepley 1282a75e6a4aSMatthew G. Knepley /* 1283a75e6a4aSMatthew G. Knepley The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that 1284a75e6a4aSMatthew G. Knepley is attached to a communicator, in this case the attribute is a PetscViewer. 1285a75e6a4aSMatthew G. Knepley */ 1286d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID; 1287a75e6a4aSMatthew G. Knepley 1288a75e6a4aSMatthew G. Knepley /*@C 1289a75e6a4aSMatthew G. Knepley PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator. 1290a75e6a4aSMatthew G. Knepley 1291d083f849SBarry Smith Collective 1292a75e6a4aSMatthew G. Knepley 1293a75e6a4aSMatthew G. Knepley Input Parameter: 1294a75e6a4aSMatthew G. Knepley . comm - the MPI communicator to share the HDF5 PetscViewer 1295a75e6a4aSMatthew G. Knepley 1296a75e6a4aSMatthew G. Knepley Level: intermediate 1297a75e6a4aSMatthew G. Knepley 1298a75e6a4aSMatthew G. Knepley Options Database Keys: 1299a75e6a4aSMatthew G. Knepley . -viewer_hdf5_filename <name> 1300a75e6a4aSMatthew G. Knepley 1301a75e6a4aSMatthew G. Knepley Environmental variables: 1302a75e6a4aSMatthew G. Knepley . PETSC_VIEWER_HDF5_FILENAME 1303a75e6a4aSMatthew G. Knepley 1304a75e6a4aSMatthew G. Knepley Notes: 1305a75e6a4aSMatthew G. Knepley Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return 1306a75e6a4aSMatthew G. Knepley an error code. The HDF5 PetscViewer is usually used in the form 1307a75e6a4aSMatthew G. Knepley $ XXXView(XXX object, PETSC_VIEWER_HDF5_(comm)); 1308a75e6a4aSMatthew G. Knepley 1309a75e6a4aSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy() 1310a75e6a4aSMatthew G. Knepley @*/ 1311a75e6a4aSMatthew G. Knepley PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm) 1312a75e6a4aSMatthew G. Knepley { 1313a75e6a4aSMatthew G. Knepley PetscErrorCode ierr; 1314a75e6a4aSMatthew G. Knepley PetscBool flg; 1315a75e6a4aSMatthew G. Knepley PetscViewer viewer; 1316a75e6a4aSMatthew G. Knepley char fname[PETSC_MAX_PATH_LEN]; 1317a75e6a4aSMatthew G. Knepley MPI_Comm ncomm; 1318a75e6a4aSMatthew G. Knepley 1319a75e6a4aSMatthew G. Knepley PetscFunctionBegin; 1320908793a3SLisandro 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);} 1321a75e6a4aSMatthew G. Knepley if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) { 1322908793a3SLisandro Dalcin ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,NULL); 1323908793a3SLisandro Dalcin if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 1324a75e6a4aSMatthew G. Knepley } 132547435625SJed Brown ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg); 1326908793a3SLisandro Dalcin if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 1327a75e6a4aSMatthew G. Knepley if (!flg) { /* PetscViewer not yet created */ 1328a75e6a4aSMatthew G. Knepley ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg); 13292cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1330a75e6a4aSMatthew G. Knepley if (!flg) { 1331a75e6a4aSMatthew G. Knepley ierr = PetscStrcpy(fname,"output.h5"); 13322cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1333a75e6a4aSMatthew G. Knepley } 1334a75e6a4aSMatthew G. Knepley ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer); 13352cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1336a75e6a4aSMatthew G. Knepley ierr = PetscObjectRegisterDestroy((PetscObject)viewer); 13372cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 133847435625SJed Brown ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer); 1339908793a3SLisandro Dalcin if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 1340a75e6a4aSMatthew G. Knepley } 1341a75e6a4aSMatthew G. Knepley ierr = PetscCommDestroy(&ncomm); 13422cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1343a75e6a4aSMatthew G. Knepley PetscFunctionReturn(viewer); 1344a75e6a4aSMatthew G. Knepley } 1345