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 86c132bc1SVaclav Hapla static PetscErrorCode PetscViewerHDF5GetAbsolutePath_Internal(PetscViewer viewer, const char objname[], char **fullpath) 96c132bc1SVaclav Hapla { 106c132bc1SVaclav Hapla const char *group; 116c132bc1SVaclav Hapla char buf[PETSC_MAX_PATH_LEN]=""; 126c132bc1SVaclav Hapla PetscErrorCode ierr; 136c132bc1SVaclav Hapla 146c132bc1SVaclav Hapla PetscFunctionBegin; 156c132bc1SVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr); 166c132bc1SVaclav Hapla ierr = PetscStrcat(buf, group);CHKERRQ(ierr); 176c132bc1SVaclav Hapla ierr = PetscStrcat(buf, "/");CHKERRQ(ierr); 186c132bc1SVaclav Hapla ierr = PetscStrcat(buf, objname);CHKERRQ(ierr); 196c132bc1SVaclav Hapla ierr = PetscStrallocpy(buf, fullpath);CHKERRQ(ierr); 206c132bc1SVaclav Hapla PetscFunctionReturn(0); 216c132bc1SVaclav Hapla } 226c132bc1SVaclav Hapla 23577e0e04SVaclav Hapla static PetscErrorCode PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer, PetscObject obj) 24577e0e04SVaclav Hapla { 25577e0e04SVaclav Hapla PetscBool has; 26577e0e04SVaclav Hapla const char *group; 27577e0e04SVaclav Hapla PetscErrorCode ierr; 28577e0e04SVaclav Hapla 29577e0e04SVaclav Hapla PetscFunctionBegin; 30577e0e04SVaclav Hapla if (!obj->name) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named"); 31577e0e04SVaclav Hapla ierr = PetscViewerHDF5HasObject(viewer, obj, &has);CHKERRQ(ierr); 32577e0e04SVaclav Hapla if (!has) { 33577e0e04SVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr); 34577e0e04SVaclav Hapla SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) %s not stored in group %s", obj->name, group); 35577e0e04SVaclav Hapla } 36577e0e04SVaclav Hapla PetscFunctionReturn(0); 37577e0e04SVaclav Hapla } 38577e0e04SVaclav Hapla 394416b707SBarry Smith static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptionItems *PetscOptionsObject,PetscViewer v) 4082ea9b62SBarry Smith { 4182ea9b62SBarry Smith PetscErrorCode ierr; 42ee8b9732SVaclav Hapla PetscBool flg = PETSC_FALSE, set; 4382ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data; 4482ea9b62SBarry Smith 4582ea9b62SBarry Smith PetscFunctionBegin; 4682ea9b62SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");CHKERRQ(ierr); 4782ea9b62SBarry Smith ierr = PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);CHKERRQ(ierr); 489a0502c6SHåkon Strandenes ierr = PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL);CHKERRQ(ierr); 49ee8b9732SVaclav Hapla ierr = PetscOptionsBool("-viewer_hdf5_collective","Enable collective transfer mode","PetscViewerHDF5SetCollective",flg,&flg,&set);CHKERRQ(ierr); 50ee8b9732SVaclav Hapla if (set) {ierr = PetscViewerHDF5SetCollective(v,flg);CHKERRQ(ierr);} 5182ea9b62SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 5282ea9b62SBarry Smith PetscFunctionReturn(0); 5382ea9b62SBarry Smith } 5482ea9b62SBarry Smith 551b793a25SVaclav Hapla static PetscErrorCode PetscViewerView_HDF5(PetscViewer v,PetscViewer viewer) 561b793a25SVaclav Hapla { 571b793a25SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data; 5803fa8834SVaclav Hapla PetscBool flg; 591b793a25SVaclav Hapla PetscErrorCode ierr; 601b793a25SVaclav Hapla 611b793a25SVaclav Hapla PetscFunctionBegin; 621b793a25SVaclav Hapla if (hdf5->filename) { 631b793a25SVaclav Hapla ierr = PetscViewerASCIIPrintf(viewer,"Filename: %s\n",hdf5->filename);CHKERRQ(ierr); 641b793a25SVaclav Hapla } 651b793a25SVaclav Hapla ierr = PetscViewerASCIIPrintf(viewer,"Vectors with blocksize 1 saved as 2D datasets: %s\n",PetscBools[hdf5->basedimension2]);CHKERRQ(ierr); 661b793a25SVaclav Hapla ierr = PetscViewerASCIIPrintf(viewer,"Enforce single precision storage: %s\n",PetscBools[hdf5->spoutput]);CHKERRQ(ierr); 6703fa8834SVaclav Hapla ierr = PetscViewerHDF5GetCollective(v,&flg);CHKERRQ(ierr); 6803fa8834SVaclav Hapla ierr = PetscViewerASCIIPrintf(viewer,"MPI-IO transfer mode: %s\n",flg ? "collective" : "independent");CHKERRQ(ierr); 691b793a25SVaclav Hapla PetscFunctionReturn(0); 701b793a25SVaclav Hapla } 711b793a25SVaclav Hapla 725c6c1daeSBarry Smith static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer) 735c6c1daeSBarry Smith { 745c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data; 755c6c1daeSBarry Smith PetscErrorCode ierr; 765c6c1daeSBarry Smith 775c6c1daeSBarry Smith PetscFunctionBegin; 785c6c1daeSBarry Smith ierr = PetscFree(hdf5->filename);CHKERRQ(ierr); 79729ab6d9SBarry Smith if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 805c6c1daeSBarry Smith PetscFunctionReturn(0); 815c6c1daeSBarry Smith } 825c6c1daeSBarry Smith 839b1403beSVaclav Hapla static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer) 845c6c1daeSBarry Smith { 855c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 865c6c1daeSBarry Smith PetscErrorCode ierr; 875c6c1daeSBarry Smith 885c6c1daeSBarry Smith PetscFunctionBegin; 899c5072fbSVaclav Hapla PetscStackCallHDF5(H5Pclose,(hdf5->dxpl_id)); 905c6c1daeSBarry Smith ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr); 915c6c1daeSBarry Smith while (hdf5->groups) { 927d964842SVaclav Hapla PetscViewerHDF5GroupList *tmp = hdf5->groups->next; 935c6c1daeSBarry Smith 945c6c1daeSBarry Smith ierr = PetscFree(hdf5->groups->name);CHKERRQ(ierr); 955c6c1daeSBarry Smith ierr = PetscFree(hdf5->groups);CHKERRQ(ierr); 965c6c1daeSBarry Smith hdf5->groups = tmp; 975c6c1daeSBarry Smith } 985c6c1daeSBarry Smith ierr = PetscFree(hdf5);CHKERRQ(ierr); 990b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr); 100d1232d7fSToby Isaac ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr); 1010b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr); 102058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetBaseDimension2_C",NULL);CHKERRQ(ierr); 103058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetSPOutput_C",NULL);CHKERRQ(ierr); 1045c6c1daeSBarry Smith PetscFunctionReturn(0); 1055c6c1daeSBarry Smith } 1065c6c1daeSBarry Smith 1079b1403beSVaclav Hapla static PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type) 1085c6c1daeSBarry Smith { 1095c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 1105c6c1daeSBarry Smith 1115c6c1daeSBarry Smith PetscFunctionBegin; 1125c6c1daeSBarry Smith hdf5->btype = type; 1135c6c1daeSBarry Smith PetscFunctionReturn(0); 1145c6c1daeSBarry Smith } 1155c6c1daeSBarry Smith 1160b62783dSVaclav Hapla static PetscErrorCode PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type) 1170b62783dSVaclav Hapla { 1180b62783dSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 1190b62783dSVaclav Hapla 1200b62783dSVaclav Hapla PetscFunctionBegin; 1210b62783dSVaclav Hapla *type = hdf5->btype; 1220b62783dSVaclav Hapla PetscFunctionReturn(0); 1230b62783dSVaclav Hapla } 1240b62783dSVaclav Hapla 1259b1403beSVaclav Hapla static PetscErrorCode PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg) 12682ea9b62SBarry Smith { 12782ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 12882ea9b62SBarry Smith 12982ea9b62SBarry Smith PetscFunctionBegin; 13082ea9b62SBarry Smith hdf5->basedimension2 = flg; 13182ea9b62SBarry Smith PetscFunctionReturn(0); 13282ea9b62SBarry Smith } 13382ea9b62SBarry Smith 134df863907SAlex Fikl /*@ 13582ea9b62SBarry Smith PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 13682ea9b62SBarry Smith dimension of 2. 13782ea9b62SBarry Smith 13882ea9b62SBarry Smith Logically Collective on PetscViewer 13982ea9b62SBarry Smith 14082ea9b62SBarry Smith Input Parameters: 14182ea9b62SBarry Smith + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 14282ea9b62SBarry 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 14382ea9b62SBarry Smith 14482ea9b62SBarry Smith Options Database: 14582ea9b62SBarry 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 14682ea9b62SBarry Smith 14782ea9b62SBarry Smith 14895452b02SPatrick Sanan Notes: 14995452b02SPatrick 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 15082ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 15182ea9b62SBarry Smith 15282ea9b62SBarry Smith Level: intermediate 15382ea9b62SBarry Smith 15482ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 15582ea9b62SBarry Smith 15682ea9b62SBarry Smith @*/ 15782ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg) 15882ea9b62SBarry Smith { 15982ea9b62SBarry Smith PetscErrorCode ierr; 16082ea9b62SBarry Smith 16182ea9b62SBarry Smith PetscFunctionBegin; 16282ea9b62SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 16382ea9b62SBarry Smith ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 16482ea9b62SBarry Smith PetscFunctionReturn(0); 16582ea9b62SBarry Smith } 16682ea9b62SBarry Smith 167df863907SAlex Fikl /*@ 16882ea9b62SBarry Smith PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 16982ea9b62SBarry Smith dimension of 2. 17082ea9b62SBarry Smith 17182ea9b62SBarry Smith Logically Collective on PetscViewer 17282ea9b62SBarry Smith 17382ea9b62SBarry Smith Input Parameter: 17482ea9b62SBarry Smith . viewer - the PetscViewer, must be of type HDF5 17582ea9b62SBarry Smith 17682ea9b62SBarry Smith Output Parameter: 17782ea9b62SBarry 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 17882ea9b62SBarry Smith 17995452b02SPatrick Sanan Notes: 18095452b02SPatrick 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 18182ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 18282ea9b62SBarry Smith 18382ea9b62SBarry Smith Level: intermediate 18482ea9b62SBarry Smith 18582ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 18682ea9b62SBarry Smith 18782ea9b62SBarry Smith @*/ 18882ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg) 18982ea9b62SBarry Smith { 19082ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 19182ea9b62SBarry Smith 19282ea9b62SBarry Smith PetscFunctionBegin; 19382ea9b62SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 19482ea9b62SBarry Smith *flg = hdf5->basedimension2; 19582ea9b62SBarry Smith PetscFunctionReturn(0); 19682ea9b62SBarry Smith } 19782ea9b62SBarry Smith 1989b1403beSVaclav Hapla static PetscErrorCode PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg) 1999a0502c6SHåkon Strandenes { 2009a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 2019a0502c6SHåkon Strandenes 2029a0502c6SHåkon Strandenes PetscFunctionBegin; 2039a0502c6SHåkon Strandenes hdf5->spoutput = flg; 2049a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2059a0502c6SHåkon Strandenes } 2069a0502c6SHåkon Strandenes 207df863907SAlex Fikl /*@ 2089a0502c6SHåkon Strandenes PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is 2099a0502c6SHåkon Strandenes compiled with double precision PetscReal. 2109a0502c6SHåkon Strandenes 2119a0502c6SHåkon Strandenes Logically Collective on PetscViewer 2129a0502c6SHåkon Strandenes 2139a0502c6SHåkon Strandenes Input Parameters: 2149a0502c6SHåkon Strandenes + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 2159a0502c6SHåkon Strandenes - flg - if PETSC_TRUE the data will be written to disk with single precision 2169a0502c6SHåkon Strandenes 2179a0502c6SHåkon Strandenes Options Database: 2189a0502c6SHåkon Strandenes . -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision 2199a0502c6SHåkon Strandenes 2209a0502c6SHåkon Strandenes 22195452b02SPatrick Sanan Notes: 22295452b02SPatrick Sanan Setting this option does not make any difference if PETSc is compiled with single precision 2239a0502c6SHåkon Strandenes in the first place. It does not affect reading datasets (HDF5 handle this internally). 2249a0502c6SHåkon Strandenes 2259a0502c6SHåkon Strandenes Level: intermediate 2269a0502c6SHåkon Strandenes 2279a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 2289a0502c6SHåkon Strandenes PetscReal 2299a0502c6SHåkon Strandenes 2309a0502c6SHåkon Strandenes @*/ 2319a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg) 2329a0502c6SHåkon Strandenes { 2339a0502c6SHåkon Strandenes PetscErrorCode ierr; 2349a0502c6SHåkon Strandenes 2359a0502c6SHåkon Strandenes PetscFunctionBegin; 2369a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 2379a0502c6SHåkon Strandenes ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 2389a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2399a0502c6SHåkon Strandenes } 2409a0502c6SHåkon Strandenes 241df863907SAlex Fikl /*@ 2429a0502c6SHåkon Strandenes PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is 2439a0502c6SHåkon Strandenes compiled with double precision PetscReal. 2449a0502c6SHåkon Strandenes 2459a0502c6SHåkon Strandenes Logically Collective on PetscViewer 2469a0502c6SHåkon Strandenes 2479a0502c6SHåkon Strandenes Input Parameter: 2489a0502c6SHåkon Strandenes . viewer - the PetscViewer, must be of type HDF5 2499a0502c6SHåkon Strandenes 2509a0502c6SHåkon Strandenes Output Parameter: 2519a0502c6SHåkon Strandenes . flg - if PETSC_TRUE the data will be written to disk with single precision 2529a0502c6SHåkon Strandenes 25395452b02SPatrick Sanan Notes: 25495452b02SPatrick Sanan Setting this option does not make any difference if PETSc is compiled with single precision 2559a0502c6SHåkon Strandenes in the first place. It does not affect reading datasets (HDF5 handle this internally). 2569a0502c6SHåkon Strandenes 2579a0502c6SHåkon Strandenes Level: intermediate 2589a0502c6SHåkon Strandenes 2599a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 2609a0502c6SHåkon Strandenes PetscReal 2619a0502c6SHåkon Strandenes 2629a0502c6SHåkon Strandenes @*/ 2639a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg) 2649a0502c6SHåkon Strandenes { 2659a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 2669a0502c6SHåkon Strandenes 2679a0502c6SHåkon Strandenes PetscFunctionBegin; 2689a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 2699a0502c6SHåkon Strandenes *flg = hdf5->spoutput; 2709a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2719a0502c6SHåkon Strandenes } 2729a0502c6SHåkon Strandenes 273ee8b9732SVaclav Hapla static PetscErrorCode PetscViewerHDF5SetCollective_HDF5(PetscViewer viewer, PetscBool flg) 274ee8b9732SVaclav Hapla { 275ee8b9732SVaclav Hapla PetscFunctionBegin; 276ee8b9732SVaclav Hapla /* H5FD_MPIO_COLLECTIVE is wrong in hdf5 1.10.2, and is the same as H5FD_MPIO_INDEPENDENT in earlier versions 277ee8b9732SVaclav Hapla - see e.g. https://gitlab.cosma.dur.ac.uk/swift/swiftsim/issues/431 */ 278c796909eSBarry Smith #if H5_VERSION_GE(1,10,3) && defined(H5_HAVE_PARALLEL) 279ee8b9732SVaclav Hapla { 280ee8b9732SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 281ee8b9732SVaclav Hapla PetscStackCallHDF5(H5Pset_dxpl_mpio,(hdf5->dxpl_id, flg ? H5FD_MPIO_COLLECTIVE : H5FD_MPIO_INDEPENDENT)); 282ee8b9732SVaclav Hapla } 283ee8b9732SVaclav Hapla #else 284ee8b9732SVaclav Hapla if (flg) { 285ee8b9732SVaclav Hapla PetscErrorCode ierr; 286c796909eSBarry 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); 287ee8b9732SVaclav Hapla } 288ee8b9732SVaclav Hapla #endif 289ee8b9732SVaclav Hapla PetscFunctionReturn(0); 290ee8b9732SVaclav Hapla } 291ee8b9732SVaclav Hapla 292ee8b9732SVaclav Hapla /*@ 293ee8b9732SVaclav Hapla PetscViewerHDF5SetCollective - Use collective MPI-IO transfer mode for HDF5 reads and writes. 294ee8b9732SVaclav Hapla 295ee8b9732SVaclav Hapla Logically Collective; flg must contain common value 296ee8b9732SVaclav Hapla 297ee8b9732SVaclav Hapla Input Parameters: 298ee8b9732SVaclav Hapla + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 299ee8b9732SVaclav Hapla - flg - PETSC_TRUE for collective mode; PETSC_FALSE for independent mode (default) 300ee8b9732SVaclav Hapla 301ee8b9732SVaclav Hapla Options Database: 302ee8b9732SVaclav Hapla . -viewer_hdf5_collective - turns on (true) or off (false) collective transfers 303ee8b9732SVaclav Hapla 304ee8b9732SVaclav Hapla Notes: 305ee8b9732SVaclav Hapla Collective mode gives the MPI-IO layer underneath HDF5 a chance to do some additional collective optimizations and hence can perform better. 30653effed7SBarry 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. 307ee8b9732SVaclav Hapla 308ee8b9732SVaclav Hapla Developer notes: 309ee8b9732SVaclav Hapla In the HDF5 layer, PETSC_TRUE / PETSC_FALSE means H5Pset_dxpl_mpio() is called with H5FD_MPIO_COLLECTIVE / H5FD_MPIO_INDEPENDENT, respectively. 310ee8b9732SVaclav Hapla This in turn means use of MPI_File_{read,write}_all / MPI_File_{read,write} in the MPI-IO layer, respectively. 311ee8b9732SVaclav Hapla See HDF5 documentation and MPI-IO documentation for details. 312ee8b9732SVaclav Hapla 313ee8b9732SVaclav Hapla Level: intermediate 314ee8b9732SVaclav Hapla 315ee8b9732SVaclav Hapla .seealso: PetscViewerHDF5GetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open() 316ee8b9732SVaclav Hapla 317ee8b9732SVaclav Hapla @*/ 318ee8b9732SVaclav Hapla PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer viewer,PetscBool flg) 319ee8b9732SVaclav Hapla { 320ee8b9732SVaclav Hapla PetscErrorCode ierr; 321ee8b9732SVaclav Hapla 322ee8b9732SVaclav Hapla PetscFunctionBegin; 323ee8b9732SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 324ee8b9732SVaclav Hapla PetscValidLogicalCollectiveBool(viewer,flg,2); 325ee8b9732SVaclav Hapla ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetCollective_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 326ee8b9732SVaclav Hapla PetscFunctionReturn(0); 327ee8b9732SVaclav Hapla } 328ee8b9732SVaclav Hapla 329ee8b9732SVaclav Hapla static PetscErrorCode PetscViewerHDF5GetCollective_HDF5(PetscViewer viewer, PetscBool *flg) 330ee8b9732SVaclav Hapla { 331c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL) 332ee8b9732SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 333ee8b9732SVaclav Hapla H5FD_mpio_xfer_t mode; 334c796909eSBarry Smith #endif 335ee8b9732SVaclav Hapla 336ee8b9732SVaclav Hapla PetscFunctionBegin; 337c796909eSBarry Smith #if !defined(H5_HAVE_PARALLEL) 338c796909eSBarry Smith *flg = PETSC_FALSE; 339c796909eSBarry Smith #else 340ee8b9732SVaclav Hapla PetscStackCallHDF5(H5Pget_dxpl_mpio,(hdf5->dxpl_id, &mode)); 341ee8b9732SVaclav Hapla *flg = (mode == H5FD_MPIO_COLLECTIVE) ? PETSC_TRUE : PETSC_FALSE; 342c796909eSBarry Smith #endif 343ee8b9732SVaclav Hapla PetscFunctionReturn(0); 344ee8b9732SVaclav Hapla } 345ee8b9732SVaclav Hapla 346ee8b9732SVaclav Hapla /*@ 347ee8b9732SVaclav Hapla PetscViewerHDF5GetCollective - Return flag whether collective MPI-IO transfer mode is used for HDF5 reads and writes. 348ee8b9732SVaclav Hapla 349ee8b9732SVaclav Hapla Not Collective 350ee8b9732SVaclav Hapla 351ee8b9732SVaclav Hapla Input Parameters: 352ee8b9732SVaclav Hapla . viewer - the HDF5 PetscViewer 353ee8b9732SVaclav Hapla 354ee8b9732SVaclav Hapla Output Parameters: 355ee8b9732SVaclav Hapla . flg - the flag 356ee8b9732SVaclav Hapla 357ee8b9732SVaclav Hapla Level: intermediate 358ee8b9732SVaclav Hapla 359ee8b9732SVaclav Hapla Notes: 360c796909eSBarry 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. 361ee8b9732SVaclav Hapla For more details, see PetscViewerHDF5SetCollective(). 362ee8b9732SVaclav Hapla 363ee8b9732SVaclav Hapla .seealso: PetscViewerHDF5SetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open() 364ee8b9732SVaclav Hapla 365ee8b9732SVaclav Hapla @*/ 366ee8b9732SVaclav Hapla PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer viewer,PetscBool *flg) 367ee8b9732SVaclav Hapla { 368ee8b9732SVaclav Hapla PetscErrorCode ierr; 369ee8b9732SVaclav Hapla 370ee8b9732SVaclav Hapla PetscFunctionBegin; 371ee8b9732SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 372534a8f05SLisandro Dalcin PetscValidBoolPointer(flg,2); 373ee8b9732SVaclav Hapla 374ee8b9732SVaclav Hapla ierr = PetscUseMethod(viewer,"PetscViewerHDF5GetCollective_C",(PetscViewer,PetscBool*),(viewer,flg));CHKERRQ(ierr); 375ee8b9732SVaclav Hapla PetscFunctionReturn(0); 376ee8b9732SVaclav Hapla } 377ee8b9732SVaclav Hapla 3789b1403beSVaclav Hapla static PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[]) 3795c6c1daeSBarry Smith { 3805c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 3815c6c1daeSBarry Smith hid_t plist_id; 3825c6c1daeSBarry Smith PetscErrorCode ierr; 3835c6c1daeSBarry Smith 3845c6c1daeSBarry Smith PetscFunctionBegin; 385f683cc58SToby Isaac if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 386f683cc58SToby Isaac if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);} 3875c6c1daeSBarry Smith ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr); 3885c6c1daeSBarry Smith /* Set up file access property list with parallel I/O access */ 389729ab6d9SBarry Smith PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS)); 390c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL) 391d28bfa55SVaclav Hapla PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), MPI_INFO_NULL)); 392c796909eSBarry Smith #endif 3935c6c1daeSBarry Smith /* Create or open the file collectively */ 3945c6c1daeSBarry Smith switch (hdf5->btype) { 3955c6c1daeSBarry Smith case FILE_MODE_READ: 396729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id)); 3975c6c1daeSBarry Smith break; 3985c6c1daeSBarry Smith case FILE_MODE_APPEND: 3997e4fd573SVaclav Hapla case FILE_MODE_UPDATE: 4007e4fd573SVaclav Hapla { 4017e4fd573SVaclav Hapla PetscBool flg; 4027e4fd573SVaclav Hapla ierr = PetscTestFile(hdf5->filename, 'r', &flg);CHKERRQ(ierr); 4037e4fd573SVaclav Hapla if (flg) PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id)); 4047e4fd573SVaclav Hapla else PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_EXCL, H5P_DEFAULT, plist_id)); 4055c6c1daeSBarry Smith break; 4067e4fd573SVaclav Hapla } 4075c6c1daeSBarry Smith case FILE_MODE_WRITE: 408729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id)); 4095c6c1daeSBarry Smith break; 4107e4fd573SVaclav Hapla case FILE_MODE_UNDEFINED: 4117e4fd573SVaclav Hapla SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 4125c6c1daeSBarry Smith default: 4137e4fd573SVaclav Hapla SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP, "Unsupported file mode %s",PetscFileModes[hdf5->btype]); 4145c6c1daeSBarry Smith } 4155c6c1daeSBarry Smith if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name); 416729ab6d9SBarry Smith PetscStackCallHDF5(H5Pclose,(plist_id)); 4175c6c1daeSBarry Smith PetscFunctionReturn(0); 4185c6c1daeSBarry Smith } 4195c6c1daeSBarry Smith 420d1232d7fSToby Isaac static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name) 421d1232d7fSToby Isaac { 422d1232d7fSToby Isaac PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data; 423d1232d7fSToby Isaac 424d1232d7fSToby Isaac PetscFunctionBegin; 425d1232d7fSToby Isaac *name = vhdf5->filename; 426d1232d7fSToby Isaac PetscFunctionReturn(0); 427d1232d7fSToby Isaac } 428d1232d7fSToby Isaac 429b723ab35SVaclav Hapla static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer) 430b723ab35SVaclav Hapla { 4315dc64a97SVaclav Hapla /* 432b723ab35SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 433b723ab35SVaclav Hapla PetscErrorCode ierr; 4345dc64a97SVaclav Hapla */ 435b723ab35SVaclav Hapla 436b723ab35SVaclav Hapla PetscFunctionBegin; 437b723ab35SVaclav Hapla PetscFunctionReturn(0); 438b723ab35SVaclav Hapla } 439b723ab35SVaclav Hapla 4408556b5ebSBarry Smith /*MC 4418556b5ebSBarry Smith PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file 4428556b5ebSBarry Smith 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; 458b00a9115SJed Brown ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr); 4595c6c1daeSBarry Smith 4605c6c1daeSBarry Smith v->data = (void*) hdf5; 4615c6c1daeSBarry Smith v->ops->destroy = PetscViewerDestroy_HDF5; 46282ea9b62SBarry Smith v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5; 463b723ab35SVaclav Hapla v->ops->setup = PetscViewerSetUp_HDF5; 4641b793a25SVaclav Hapla v->ops->view = PetscViewerView_HDF5; 465908793a3SLisandro Dalcin v->ops->flush = NULL; 4667e4fd573SVaclav Hapla hdf5->btype = FILE_MODE_UNDEFINED; 467908793a3SLisandro Dalcin hdf5->filename = NULL; 4685c6c1daeSBarry Smith hdf5->timestep = -1; 4690298fd71SBarry Smith hdf5->groups = NULL; 4705c6c1daeSBarry Smith 4719c5072fbSVaclav Hapla PetscStackCallHDF5Return(hdf5->dxpl_id,H5Pcreate,(H5P_DATASET_XFER)); 4729c5072fbSVaclav Hapla 4730b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr); 474d1232d7fSToby Isaac ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr); 4750b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr); 4760b62783dSVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_HDF5);CHKERRQ(ierr); 47782ea9b62SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr); 4789a0502c6SHåkon Strandenes ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr); 479ee8b9732SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetCollective_C",PetscViewerHDF5SetCollective_HDF5);CHKERRQ(ierr); 480ee8b9732SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetCollective_C",PetscViewerHDF5GetCollective_HDF5);CHKERRQ(ierr); 4815c6c1daeSBarry Smith PetscFunctionReturn(0); 4825c6c1daeSBarry Smith } 4835c6c1daeSBarry Smith 4845c6c1daeSBarry Smith /*@C 4855c6c1daeSBarry Smith PetscViewerHDF5Open - Opens a file for HDF5 input/output. 4865c6c1daeSBarry Smith 487d083f849SBarry Smith Collective 4885c6c1daeSBarry Smith 4895c6c1daeSBarry Smith Input Parameters: 4905c6c1daeSBarry Smith + comm - MPI communicator 4915c6c1daeSBarry Smith . name - name of file 4925c6c1daeSBarry Smith - type - type of file 4935c6c1daeSBarry Smith 4945c6c1daeSBarry Smith Output Parameter: 4955c6c1daeSBarry Smith . hdf5v - PetscViewer for HDF5 input/output to use with the specified file 4965c6c1daeSBarry Smith 49782ea9b62SBarry Smith Options Database: 498a2b725a8SWilliam 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 499a2b725a8SWilliam Gropp - -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal 50082ea9b62SBarry Smith 5015c6c1daeSBarry Smith Level: beginner 5025c6c1daeSBarry Smith 5037e4fd573SVaclav Hapla Notes: 5047e4fd573SVaclav Hapla Reading is always available, regardless of the mode. Available modes are 5057e4fd573SVaclav Hapla + FILE_MODE_READ - open existing HDF5 file for read only access, fail if file does not exist [H5Fopen() with H5F_ACC_RDONLY] 5067e4fd573SVaclav Hapla . FILE_MODE_WRITE - if file exists, fully overwrite it, else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_TRUNC] 5077e4fd573SVaclav 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] 5087e4fd573SVaclav Hapla - FILE_MODE_UPDATE - same as FILE_MODE_APPEND 5097e4fd573SVaclav Hapla 5107e4fd573SVaclav 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. 5117e4fd573SVaclav Hapla 5125c6c1daeSBarry Smith This PetscViewer should be destroyed with PetscViewerDestroy(). 5135c6c1daeSBarry Smith 5145c6c1daeSBarry Smith 5156a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(), 5169a0502c6SHåkon Strandenes PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(), 517a56f64adSBarry Smith MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName() 5185c6c1daeSBarry Smith @*/ 5195c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v) 5205c6c1daeSBarry Smith { 5215c6c1daeSBarry Smith PetscErrorCode ierr; 5225c6c1daeSBarry Smith 5235c6c1daeSBarry Smith PetscFunctionBegin; 5245c6c1daeSBarry Smith ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr); 5255c6c1daeSBarry Smith ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr); 5265c6c1daeSBarry Smith ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr); 5275c6c1daeSBarry Smith ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr); 528b663d2d2SVaclav Hapla ierr = PetscViewerSetFromOptions(*hdf5v);CHKERRQ(ierr); 5295c6c1daeSBarry Smith PetscFunctionReturn(0); 5305c6c1daeSBarry Smith } 5315c6c1daeSBarry Smith 5325c6c1daeSBarry Smith /*@C 5335c6c1daeSBarry Smith PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls 5345c6c1daeSBarry Smith 5355c6c1daeSBarry Smith Not collective 5365c6c1daeSBarry Smith 5375c6c1daeSBarry Smith Input Parameter: 5385c6c1daeSBarry Smith . viewer - the PetscViewer 5395c6c1daeSBarry Smith 5405c6c1daeSBarry Smith Output Parameter: 5415c6c1daeSBarry Smith . file_id - The file id 5425c6c1daeSBarry Smith 5435c6c1daeSBarry Smith Level: intermediate 5445c6c1daeSBarry Smith 5455c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open() 5465c6c1daeSBarry Smith @*/ 5475c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id) 5485c6c1daeSBarry Smith { 5495c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 5505c6c1daeSBarry Smith 5515c6c1daeSBarry Smith PetscFunctionBegin; 5525c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 5535c6c1daeSBarry Smith if (file_id) *file_id = hdf5->file_id; 5545c6c1daeSBarry Smith PetscFunctionReturn(0); 5555c6c1daeSBarry Smith } 5565c6c1daeSBarry Smith 5575c6c1daeSBarry Smith /*@C 5585c6c1daeSBarry Smith PetscViewerHDF5PushGroup - Set the current HDF5 group for output 5595c6c1daeSBarry Smith 5605c6c1daeSBarry Smith Not collective 5615c6c1daeSBarry Smith 5625c6c1daeSBarry Smith Input Parameters: 5635c6c1daeSBarry Smith + viewer - the PetscViewer 5645c6c1daeSBarry Smith - name - The group name 5655c6c1daeSBarry Smith 5665c6c1daeSBarry Smith Level: intermediate 5675c6c1daeSBarry Smith 568*2e361344SVaclav Hapla Notes: 569*2e361344SVaclav Hapla This is designed to mnemonically resemble the Unix cd command. 570*2e361344SVaclav 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. 571*2e361344SVaclav Hapla . NULL, empty string, or any sequence of all slashes (e.g. "///") is interpreted as the root group "/". 572*2e361344SVaclav Hapla - "." means the current group is pushed again. 573*2e361344SVaclav Hapla 574*2e361344SVaclav Hapla Example: 575*2e361344SVaclav Hapla Suppose the current group is "/a". 576*2e361344SVaclav Hapla + If name is NULL, empty string, or a sequence of all slashes (e.g. "///"), then the new group will be "/". 577*2e361344SVaclav Hapla . If name is ".", then the new group will be "/a". 578*2e361344SVaclav Hapla . If name is "b", then the new group will be "/a/b". 579*2e361344SVaclav Hapla - If name is "/b", then the new group will be "/b". 580*2e361344SVaclav Hapla 581*2e361344SVaclav Hapla Developer Notes: 582*2e361344SVaclav Hapla The root group "/" is internally stored as NULL. 583820fc2d1SVaclav Hapla 584874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 5855c6c1daeSBarry Smith @*/ 586be7aa430SVaclav Hapla PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char name[]) 5875c6c1daeSBarry Smith { 5885c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 5897d964842SVaclav Hapla PetscViewerHDF5GroupList *groupNode; 590*2e361344SVaclav Hapla size_t i,len; 591*2e361344SVaclav Hapla char buf[PETSC_MAX_PATH_LEN]; 592*2e361344SVaclav Hapla const char *gname; 5935c6c1daeSBarry Smith PetscErrorCode ierr; 5945c6c1daeSBarry Smith 5955c6c1daeSBarry Smith PetscFunctionBegin; 5965c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 597820fc2d1SVaclav Hapla if (name) PetscValidCharPointer(name,2); 598820fc2d1SVaclav Hapla ierr = PetscStrlen(name, &len);CHKERRQ(ierr); 599*2e361344SVaclav Hapla gname = NULL; 600*2e361344SVaclav Hapla if (len) { 601*2e361344SVaclav Hapla if (len == 1 && name[0] == '.') { 602*2e361344SVaclav Hapla /* use current name */ 603*2e361344SVaclav Hapla gname = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : NULL; 604*2e361344SVaclav Hapla } else if (name[0] == '/') { 605*2e361344SVaclav Hapla /* absolute */ 606*2e361344SVaclav Hapla for (i=1; i<len; i++) { 607*2e361344SVaclav Hapla if (name[i] != '/') { 608*2e361344SVaclav Hapla gname = name; 609*2e361344SVaclav Hapla break; 610*2e361344SVaclav Hapla } 611*2e361344SVaclav Hapla } 612*2e361344SVaclav Hapla } else { 613*2e361344SVaclav Hapla /* relative */ 614*2e361344SVaclav Hapla const char *parent = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : ""; 615*2e361344SVaclav Hapla ierr = PetscSNPrintf(buf, sizeof(buf), "%s/%s", parent, name);CHKERRQ(ierr); 616*2e361344SVaclav Hapla gname = buf; 617*2e361344SVaclav Hapla } 618*2e361344SVaclav Hapla } 61995dccacaSBarry Smith ierr = PetscNew(&groupNode);CHKERRQ(ierr); 620*2e361344SVaclav Hapla ierr = PetscStrallocpy(gname, (char**) &groupNode->name);CHKERRQ(ierr); 6215c6c1daeSBarry Smith groupNode->next = hdf5->groups; 6225c6c1daeSBarry Smith hdf5->groups = groupNode; 6235c6c1daeSBarry Smith PetscFunctionReturn(0); 6245c6c1daeSBarry Smith } 6255c6c1daeSBarry Smith 6263ef9c667SSatish Balay /*@ 6275c6c1daeSBarry Smith PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value 6285c6c1daeSBarry Smith 6295c6c1daeSBarry Smith Not collective 6305c6c1daeSBarry Smith 6315c6c1daeSBarry Smith Input Parameter: 6325c6c1daeSBarry Smith . viewer - the PetscViewer 6335c6c1daeSBarry Smith 6345c6c1daeSBarry Smith Level: intermediate 6355c6c1daeSBarry Smith 636874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 6375c6c1daeSBarry Smith @*/ 6385c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer) 6395c6c1daeSBarry Smith { 6405c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 6417d964842SVaclav Hapla PetscViewerHDF5GroupList *groupNode; 6425c6c1daeSBarry Smith PetscErrorCode ierr; 6435c6c1daeSBarry Smith 6445c6c1daeSBarry Smith PetscFunctionBegin; 6455c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 64682f516ccSBarry Smith if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop"); 6475c6c1daeSBarry Smith groupNode = hdf5->groups; 6485c6c1daeSBarry Smith hdf5->groups = hdf5->groups->next; 6495c6c1daeSBarry Smith ierr = PetscFree(groupNode->name);CHKERRQ(ierr); 6505c6c1daeSBarry Smith ierr = PetscFree(groupNode);CHKERRQ(ierr); 6515c6c1daeSBarry Smith PetscFunctionReturn(0); 6525c6c1daeSBarry Smith } 6535c6c1daeSBarry Smith 6545c6c1daeSBarry Smith /*@C 655874270d9SVaclav Hapla PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup(). 656874270d9SVaclav Hapla If none has been assigned, returns NULL. 6575c6c1daeSBarry Smith 6585c6c1daeSBarry Smith Not collective 6595c6c1daeSBarry Smith 6605c6c1daeSBarry Smith Input Parameter: 6615c6c1daeSBarry Smith . viewer - the PetscViewer 6625c6c1daeSBarry Smith 6635c6c1daeSBarry Smith Output Parameter: 6645c6c1daeSBarry Smith . name - The group name 6655c6c1daeSBarry Smith 6665c6c1daeSBarry Smith Level: intermediate 6675c6c1daeSBarry Smith 668874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup() 6695c6c1daeSBarry Smith @*/ 670be7aa430SVaclav Hapla PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char *name[]) 6715c6c1daeSBarry Smith { 6725c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 6735c6c1daeSBarry Smith 6745c6c1daeSBarry Smith PetscFunctionBegin; 6755c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 676c959eef4SJed Brown PetscValidPointer(name,2); 677a297a907SKarl Rupp if (hdf5->groups) *name = hdf5->groups->name; 6780298fd71SBarry Smith else *name = NULL; 6795c6c1daeSBarry Smith PetscFunctionReturn(0); 6805c6c1daeSBarry Smith } 6815c6c1daeSBarry Smith 6828c557b5aSVaclav Hapla /*@ 683874270d9SVaclav Hapla PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(), 684874270d9SVaclav Hapla and return this group's ID and file ID. 685874270d9SVaclav Hapla If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID. 686874270d9SVaclav Hapla 687874270d9SVaclav Hapla Not collective 688874270d9SVaclav Hapla 689874270d9SVaclav Hapla Input Parameter: 690874270d9SVaclav Hapla . viewer - the PetscViewer 691874270d9SVaclav Hapla 692874270d9SVaclav Hapla Output Parameter: 693874270d9SVaclav Hapla + fileId - The HDF5 file ID 694874270d9SVaclav Hapla - groupId - The HDF5 group ID 695874270d9SVaclav Hapla 696e74616beSVaclav Hapla Notes: 697e74616beSVaclav Hapla If the viewer is writable, the group is created if it doesn't exist yet. 698e74616beSVaclav Hapla 699874270d9SVaclav Hapla Level: intermediate 700874270d9SVaclav Hapla 701874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 702874270d9SVaclav Hapla @*/ 70354dbf706SVaclav Hapla PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId) 70454dbf706SVaclav Hapla { 70590e03537SVaclav Hapla hid_t file_id; 70690e03537SVaclav Hapla H5O_type_t type; 70754dbf706SVaclav Hapla const char *groupName = NULL; 708e74616beSVaclav Hapla PetscBool create; 70954dbf706SVaclav Hapla PetscErrorCode ierr; 71054dbf706SVaclav Hapla 71154dbf706SVaclav Hapla PetscFunctionBegin; 712e74616beSVaclav Hapla ierr = PetscViewerWritable(viewer, &create);CHKERRQ(ierr); 71354dbf706SVaclav Hapla ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr); 71454dbf706SVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr); 715e74616beSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, groupName, create, NULL, &type);CHKERRQ(ierr); 71690e03537SVaclav 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); 71790e03537SVaclav Hapla PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT)); 71854dbf706SVaclav Hapla *fileId = file_id; 71954dbf706SVaclav Hapla PetscFunctionReturn(0); 72054dbf706SVaclav Hapla } 72154dbf706SVaclav Hapla 7223ef9c667SSatish Balay /*@ 7235c6c1daeSBarry Smith PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time. 7245c6c1daeSBarry Smith 7255c6c1daeSBarry Smith Not collective 7265c6c1daeSBarry Smith 7275c6c1daeSBarry Smith Input Parameter: 7285c6c1daeSBarry Smith . viewer - the PetscViewer 7295c6c1daeSBarry Smith 7305c6c1daeSBarry Smith Level: intermediate 7315c6c1daeSBarry Smith 7325c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep() 7335c6c1daeSBarry Smith @*/ 7345c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer) 7355c6c1daeSBarry Smith { 7365c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 7375c6c1daeSBarry Smith 7385c6c1daeSBarry Smith PetscFunctionBegin; 7395c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 7405c6c1daeSBarry Smith ++hdf5->timestep; 7415c6c1daeSBarry Smith PetscFunctionReturn(0); 7425c6c1daeSBarry Smith } 7435c6c1daeSBarry Smith 7443ef9c667SSatish Balay /*@ 7455c6c1daeSBarry Smith PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep 7465c6c1daeSBarry Smith of -1 disables blocking with timesteps. 7475c6c1daeSBarry Smith 7485c6c1daeSBarry Smith Not collective 7495c6c1daeSBarry Smith 7505c6c1daeSBarry Smith Input Parameters: 7515c6c1daeSBarry Smith + viewer - the PetscViewer 7525c6c1daeSBarry Smith - timestep - The timestep number 7535c6c1daeSBarry Smith 7545c6c1daeSBarry Smith Level: intermediate 7555c6c1daeSBarry Smith 7565c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 7575c6c1daeSBarry Smith @*/ 7585c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep) 7595c6c1daeSBarry Smith { 7605c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 7615c6c1daeSBarry Smith 7625c6c1daeSBarry Smith PetscFunctionBegin; 7635c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 7645c6c1daeSBarry Smith hdf5->timestep = timestep; 7655c6c1daeSBarry Smith PetscFunctionReturn(0); 7665c6c1daeSBarry Smith } 7675c6c1daeSBarry Smith 7683ef9c667SSatish Balay /*@ 7695c6c1daeSBarry Smith PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time. 7705c6c1daeSBarry Smith 7715c6c1daeSBarry Smith Not collective 7725c6c1daeSBarry Smith 7735c6c1daeSBarry Smith Input Parameter: 7745c6c1daeSBarry Smith . viewer - the PetscViewer 7755c6c1daeSBarry Smith 7765c6c1daeSBarry Smith Output Parameter: 7775c6c1daeSBarry Smith . timestep - The timestep number 7785c6c1daeSBarry Smith 7795c6c1daeSBarry Smith Level: intermediate 7805c6c1daeSBarry Smith 7815c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep() 7825c6c1daeSBarry Smith @*/ 7835c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep) 7845c6c1daeSBarry Smith { 7855c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 7865c6c1daeSBarry Smith 7875c6c1daeSBarry Smith PetscFunctionBegin; 7885c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 7895c6c1daeSBarry Smith PetscValidPointer(timestep,2); 7905c6c1daeSBarry Smith *timestep = hdf5->timestep; 7915c6c1daeSBarry Smith PetscFunctionReturn(0); 7925c6c1daeSBarry Smith } 7935c6c1daeSBarry Smith 79436bce6e8SMatthew G. Knepley /*@C 79536bce6e8SMatthew G. Knepley PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name. 79636bce6e8SMatthew G. Knepley 79736bce6e8SMatthew G. Knepley Not collective 79836bce6e8SMatthew G. Knepley 79936bce6e8SMatthew G. Knepley Input Parameter: 80036bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 80136bce6e8SMatthew G. Knepley 80236bce6e8SMatthew G. Knepley Output Parameter: 80336bce6e8SMatthew G. Knepley . mtype - the MPI datatype (for example MPI_DOUBLE, ...) 80436bce6e8SMatthew G. Knepley 80536bce6e8SMatthew G. Knepley Level: advanced 80636bce6e8SMatthew G. Knepley 80736bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 80836bce6e8SMatthew G. Knepley @*/ 80936bce6e8SMatthew G. Knepley PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype) 81036bce6e8SMatthew G. Knepley { 81136bce6e8SMatthew G. Knepley PetscFunctionBegin; 81236bce6e8SMatthew G. Knepley if (ptype == PETSC_INT) 81336bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 81436bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_LLONG; 81536bce6e8SMatthew G. Knepley #else 81636bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_INT; 81736bce6e8SMatthew G. Knepley #endif 81836bce6e8SMatthew G. Knepley else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE; 81936bce6e8SMatthew G. Knepley else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG; 82036bce6e8SMatthew G. Knepley else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT; 821c9450970SBarry Smith else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_INT; 822de3ba262SVaclav Hapla else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_INT; 82336bce6e8SMatthew G. Knepley else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT; 82436bce6e8SMatthew G. Knepley else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR; 82536bce6e8SMatthew G. Knepley else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR; 8267e97c476SMatthew G. Knepley else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1); 82736bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype"); 82836bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 82936bce6e8SMatthew G. Knepley } 83036bce6e8SMatthew G. Knepley 83136bce6e8SMatthew G. Knepley /*@C 83236bce6e8SMatthew G. Knepley PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name 83336bce6e8SMatthew G. Knepley 83436bce6e8SMatthew G. Knepley Not collective 83536bce6e8SMatthew G. Knepley 83636bce6e8SMatthew G. Knepley Input Parameter: 83736bce6e8SMatthew G. Knepley . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...) 83836bce6e8SMatthew G. Knepley 83936bce6e8SMatthew G. Knepley Output Parameter: 84036bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 84136bce6e8SMatthew G. Knepley 84236bce6e8SMatthew G. Knepley Level: advanced 84336bce6e8SMatthew G. Knepley 84436bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 84536bce6e8SMatthew G. Knepley @*/ 84636bce6e8SMatthew G. Knepley PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype) 84736bce6e8SMatthew G. Knepley { 84836bce6e8SMatthew G. Knepley PetscFunctionBegin; 84936bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 85036bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG; 85136bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT; 85236bce6e8SMatthew G. Knepley #else 85336bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT; 85436bce6e8SMatthew G. Knepley #endif 85536bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE; 85636bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG; 85736bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT; 85836bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT; 85936bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR; 86036bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR; 8617e97c476SMatthew G. Knepley else if (htype == H5T_C_S1) *ptype = PETSC_STRING; 86236bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype"); 86336bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 86436bce6e8SMatthew G. Knepley } 86536bce6e8SMatthew G. Knepley 866df863907SAlex Fikl /*@C 867b8ee85a1SVaclav Hapla PetscViewerHDF5WriteAttribute - Write an attribute 86836bce6e8SMatthew G. Knepley 86936bce6e8SMatthew G. Knepley Input Parameters: 87036bce6e8SMatthew G. Knepley + viewer - The HDF5 viewer 87157d22f7dSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 87236bce6e8SMatthew G. Knepley . name - The attribute name 87336bce6e8SMatthew G. Knepley . datatype - The attribute type 87436bce6e8SMatthew G. Knepley - value - The attribute value 87536bce6e8SMatthew G. Knepley 87636bce6e8SMatthew G. Knepley Level: advanced 87736bce6e8SMatthew G. Knepley 878577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 87936bce6e8SMatthew G. Knepley @*/ 88057d22f7dSVaclav Hapla PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, const void *value) 88136bce6e8SMatthew G. Knepley { 88257d22f7dSVaclav Hapla char *parent; 88360568592SMatthew G. Knepley hid_t h5, dataspace, obj, attribute, dtype; 884080f144cSVaclav Hapla PetscBool has; 88536bce6e8SMatthew G. Knepley PetscErrorCode ierr; 88636bce6e8SMatthew G. Knepley 88736bce6e8SMatthew G. Knepley PetscFunctionBegin; 8885cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 88957d22f7dSVaclav Hapla if (dataset) PetscValidCharPointer(dataset, 2); 890c57a1a9aSVaclav Hapla PetscValidCharPointer(name, 3); 891b7e81f0fSVaclav Hapla PetscValidPointer(value, 5); 89257d22f7dSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 893bb286ee1SVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_TRUE, NULL, NULL);CHKERRQ(ierr); 894b67695ecSVaclav Hapla ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr); 89536bce6e8SMatthew G. Knepley ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 8967e97c476SMatthew G. Knepley if (datatype == PETSC_STRING) { 8977e97c476SMatthew G. Knepley size_t len; 8987e97c476SMatthew G. Knepley ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr); 899729ab6d9SBarry Smith PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 9007e97c476SMatthew G. Knepley } 90136bce6e8SMatthew G. Knepley ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 902729ab6d9SBarry Smith PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR)); 90360568592SMatthew G. Knepley PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT)); 904080f144cSVaclav Hapla if (has) { 905080f144cSVaclav Hapla PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 906080f144cSVaclav Hapla } else { 90760568592SMatthew G. Knepley PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT)); 908080f144cSVaclav Hapla } 909729ab6d9SBarry Smith PetscStackCallHDF5(H5Awrite,(attribute, dtype, value)); 910729ab6d9SBarry Smith if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype)); 911729ab6d9SBarry Smith PetscStackCallHDF5(H5Aclose,(attribute)); 91260568592SMatthew G. Knepley PetscStackCallHDF5(H5Oclose,(obj)); 913729ab6d9SBarry Smith PetscStackCallHDF5(H5Sclose,(dataspace)); 91457d22f7dSVaclav Hapla ierr = PetscFree(parent);CHKERRQ(ierr); 91536bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 91636bce6e8SMatthew G. Knepley } 91736bce6e8SMatthew G. Knepley 918df863907SAlex Fikl /*@C 919577e0e04SVaclav Hapla PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name 920577e0e04SVaclav Hapla 921577e0e04SVaclav Hapla Input Parameters: 922577e0e04SVaclav Hapla + viewer - The HDF5 viewer 923577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 924577e0e04SVaclav Hapla . name - The attribute name 925577e0e04SVaclav Hapla . datatype - The attribute type 926577e0e04SVaclav Hapla - value - The attribute value 927577e0e04SVaclav Hapla 928577e0e04SVaclav Hapla Notes: 929577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 930577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 931577e0e04SVaclav Hapla 932577e0e04SVaclav Hapla Level: advanced 933577e0e04SVaclav Hapla 934577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 935577e0e04SVaclav Hapla @*/ 936577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value) 937577e0e04SVaclav Hapla { 938577e0e04SVaclav Hapla PetscErrorCode ierr; 939577e0e04SVaclav Hapla 940577e0e04SVaclav Hapla PetscFunctionBegin; 941577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 942577e0e04SVaclav Hapla PetscValidHeader(obj,2); 943577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 944b7e81f0fSVaclav Hapla PetscValidPointer(value,5); 945577e0e04SVaclav Hapla ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 946577e0e04SVaclav Hapla ierr = PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr); 947577e0e04SVaclav Hapla PetscFunctionReturn(0); 948577e0e04SVaclav Hapla } 949577e0e04SVaclav Hapla 950577e0e04SVaclav Hapla /*@C 951b8ee85a1SVaclav Hapla PetscViewerHDF5ReadAttribute - Read an attribute 952ad0c61feSMatthew G. Knepley 953ad0c61feSMatthew G. Knepley Input Parameters: 954ad0c61feSMatthew G. Knepley + viewer - The HDF5 viewer 95557d22f7dSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 956ad0c61feSMatthew G. Knepley . name - The attribute name 957ad0c61feSMatthew G. Knepley - datatype - The attribute type 958ad0c61feSMatthew G. Knepley 959ad0c61feSMatthew G. Knepley Output Parameter: 960ad0c61feSMatthew G. Knepley . value - The attribute value 961ad0c61feSMatthew G. Knepley 962708d4cb3SBarry Smith Notes: If the datatype is PETSC_STRING one must PetscFree() the obtained value when it is no longer needed. 963708d4cb3SBarry Smith 964ad0c61feSMatthew G. Knepley Level: advanced 965ad0c61feSMatthew G. Knepley 966577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 967ad0c61feSMatthew G. Knepley @*/ 96857d22f7dSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, void *value) 969ad0c61feSMatthew G. Knepley { 97057d22f7dSVaclav Hapla char *parent; 971f0b58479SMatthew G. Knepley hid_t h5, obj, attribute, atype, dtype; 972969748fdSVaclav Hapla PetscBool has; 973ad0c61feSMatthew G. Knepley PetscErrorCode ierr; 974ad0c61feSMatthew G. Knepley 975ad0c61feSMatthew G. Knepley PetscFunctionBegin; 9765cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 97757d22f7dSVaclav Hapla if (dataset) PetscValidCharPointer(dataset, 2); 978c57a1a9aSVaclav Hapla PetscValidCharPointer(name, 3); 979b7e81f0fSVaclav Hapla PetscValidPointer(value, 5); 98057d22f7dSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 981969748fdSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, &has, NULL);CHKERRQ(ierr); 982969748fdSVaclav Hapla if (has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr);} 983969748fdSVaclav Hapla if (!has) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist", parent, name); 984ad0c61feSMatthew G. Knepley ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 985ad0c61feSMatthew G. Knepley ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 98660568592SMatthew G. Knepley PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT)); 98760568592SMatthew G. Knepley PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 988f0b58479SMatthew G. Knepley if (datatype == PETSC_STRING) { 989f0b58479SMatthew G. Knepley size_t len; 990e90c5075SMarek Pecha PetscStackCallHDF5Return(atype,H5Aget_type,(attribute)); 99115b861d2SVaclav Hapla PetscStackCall("H5Tget_size",len = H5Tget_size(atype)); 992708d4cb3SBarry Smith ierr = PetscMalloc((len+1) * sizeof(char), value);CHKERRQ(ierr); 993708d4cb3SBarry Smith PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 994708d4cb3SBarry Smith PetscStackCallHDF5(H5Aread,(attribute, dtype, *(char**)value)); 995708d4cb3SBarry Smith } else { 99670efba33SBarry Smith PetscStackCallHDF5(H5Aread,(attribute, dtype, value)); 997708d4cb3SBarry Smith } 998729ab6d9SBarry Smith PetscStackCallHDF5(H5Aclose,(attribute)); 999e90c5075SMarek Pecha /* H5Oclose can be used to close groups, datasets, or committed datatypes */ 1000e90c5075SMarek Pecha PetscStackCallHDF5(H5Oclose,(obj)); 100157d22f7dSVaclav Hapla ierr = PetscFree(parent);CHKERRQ(ierr); 1002ad0c61feSMatthew G. Knepley PetscFunctionReturn(0); 1003ad0c61feSMatthew G. Knepley } 1004ad0c61feSMatthew G. Knepley 1005577e0e04SVaclav Hapla /*@C 1006577e0e04SVaclav Hapla PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name 1007577e0e04SVaclav Hapla 1008577e0e04SVaclav Hapla Input Parameters: 1009577e0e04SVaclav Hapla + viewer - The HDF5 viewer 1010577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1011577e0e04SVaclav Hapla . name - The attribute name 1012577e0e04SVaclav Hapla - datatype - The attribute type 1013577e0e04SVaclav Hapla 1014577e0e04SVaclav Hapla Output Parameter: 1015577e0e04SVaclav Hapla . value - The attribute value 1016577e0e04SVaclav Hapla 1017577e0e04SVaclav Hapla Notes: 1018577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1019577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 1020577e0e04SVaclav Hapla 1021577e0e04SVaclav Hapla Level: advanced 1022577e0e04SVaclav Hapla 1023577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute() PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1024577e0e04SVaclav Hapla @*/ 1025577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *value) 1026577e0e04SVaclav Hapla { 1027577e0e04SVaclav Hapla PetscErrorCode ierr; 1028577e0e04SVaclav Hapla 1029577e0e04SVaclav Hapla PetscFunctionBegin; 1030577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1031577e0e04SVaclav Hapla PetscValidHeader(obj,2); 1032577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 1033b7e81f0fSVaclav Hapla PetscValidPointer(value, 5); 1034577e0e04SVaclav Hapla ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 1035577e0e04SVaclav Hapla ierr = PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr); 1036577e0e04SVaclav Hapla PetscFunctionReturn(0); 1037577e0e04SVaclav Hapla } 1038577e0e04SVaclav Hapla 1039a75c3fd4SVaclav Hapla PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_) 1040a75c3fd4SVaclav Hapla { 1041a75c3fd4SVaclav Hapla htri_t exists; 1042a75c3fd4SVaclav Hapla hid_t group; 1043a75c3fd4SVaclav Hapla 1044a75c3fd4SVaclav Hapla PetscFunctionBegin; 1045a75c3fd4SVaclav Hapla PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT)); 1046a75c3fd4SVaclav Hapla if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT)); 1047a75c3fd4SVaclav Hapla if (!exists && createGroup) { 1048a75c3fd4SVaclav Hapla PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT)); 1049a75c3fd4SVaclav Hapla PetscStackCallHDF5(H5Gclose,(group)); 1050a75c3fd4SVaclav Hapla exists = PETSC_TRUE; 1051a75c3fd4SVaclav Hapla } 1052a75c3fd4SVaclav Hapla *exists_ = (PetscBool) exists; 1053a75c3fd4SVaclav Hapla PetscFunctionReturn(0); 1054a75c3fd4SVaclav Hapla } 1055a75c3fd4SVaclav Hapla 1056bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype) 10575cdeb108SMatthew G. Knepley { 105890e03537SVaclav Hapla const char rootGroupName[] = "/"; 10595cdeb108SMatthew G. Knepley hid_t h5; 1060e5bf9ebcSVaclav Hapla PetscBool exists=PETSC_FALSE; 106115b861d2SVaclav Hapla PetscInt i; 106215b861d2SVaclav Hapla int n; 106385688ae2SVaclav Hapla char **hierarchy; 106485688ae2SVaclav Hapla char buf[PETSC_MAX_PATH_LEN]=""; 10655cdeb108SMatthew G. Knepley PetscErrorCode ierr; 10665cdeb108SMatthew G. Knepley 10675cdeb108SMatthew G. Knepley PetscFunctionBegin; 10685cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 106990e03537SVaclav Hapla if (name) PetscValidCharPointer(name, 2); 107090e03537SVaclav Hapla else name = rootGroupName; 1071ccd66a83SVaclav Hapla if (has) { 107256cc0592SVaclav Hapla PetscValidIntPointer(has, 3); 10735cdeb108SMatthew G. Knepley *has = PETSC_FALSE; 1074ccd66a83SVaclav Hapla } 1075ccd66a83SVaclav Hapla if (otype) { 1076ccd66a83SVaclav Hapla PetscValidIntPointer(otype, 4); 107756cc0592SVaclav Hapla *otype = H5O_TYPE_UNKNOWN; 1078ccd66a83SVaclav Hapla } 1079c0cd0301SJed Brown ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 108085688ae2SVaclav Hapla 108185688ae2SVaclav Hapla /* 108285688ae2SVaclav Hapla Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing. 108385688ae2SVaclav Hapla Hence, each of them needs to be tested separately: 108485688ae2SVaclav Hapla 1) whether it's a valid link 108585688ae2SVaclav Hapla 2) whether this link resolves to an object 108685688ae2SVaclav Hapla See H5Oexists_by_name() documentation. 108785688ae2SVaclav Hapla */ 108885688ae2SVaclav Hapla ierr = PetscStrToArray(name,'/',&n,&hierarchy);CHKERRQ(ierr); 108985688ae2SVaclav Hapla if (!n) { 109085688ae2SVaclav Hapla /* Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */ 1091ccd66a83SVaclav Hapla if (has) *has = PETSC_TRUE; 1092ccd66a83SVaclav Hapla if (otype) *otype = H5O_TYPE_GROUP; 109385688ae2SVaclav Hapla ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 109485688ae2SVaclav Hapla PetscFunctionReturn(0); 109585688ae2SVaclav Hapla } 109685688ae2SVaclav Hapla for (i=0; i<n; i++) { 109785688ae2SVaclav Hapla ierr = PetscStrcat(buf,"/");CHKERRQ(ierr); 109885688ae2SVaclav Hapla ierr = PetscStrcat(buf,hierarchy[i]);CHKERRQ(ierr); 1099a75c3fd4SVaclav Hapla ierr = PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);CHKERRQ(ierr); 1100a75c3fd4SVaclav Hapla if (!exists) break; 110185688ae2SVaclav Hapla } 110285688ae2SVaclav Hapla ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 110385688ae2SVaclav Hapla 110485688ae2SVaclav Hapla /* If the object exists, get its type */ 1105ccd66a83SVaclav Hapla if (exists && otype) { 11065cdeb108SMatthew G. Knepley H5O_info_t info; 11075cdeb108SMatthew G. Knepley 110876276748SVaclav Hapla /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */ 110904633f7fSVaclav Hapla PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT)); 111056cc0592SVaclav Hapla *otype = info.type; 11115cdeb108SMatthew G. Knepley } 1112ccd66a83SVaclav Hapla if (has) *has = exists; 11135cdeb108SMatthew G. Knepley PetscFunctionReturn(0); 11145cdeb108SMatthew G. Knepley } 11155cdeb108SMatthew G. Knepley 1116ecce1506SVaclav Hapla /*@ 11170a9f38efSVaclav Hapla PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file 11180a9f38efSVaclav Hapla 11190a9f38efSVaclav Hapla Input Parameters: 11200a9f38efSVaclav Hapla . viewer - The HDF5 viewer 11210a9f38efSVaclav Hapla 11220a9f38efSVaclav Hapla Output Parameter: 11230a9f38efSVaclav Hapla . has - Flag for group existence 11240a9f38efSVaclav Hapla 11250a9f38efSVaclav Hapla Notes: 11260a9f38efSVaclav Hapla If the path exists but is not a group, this returns PETSC_FALSE as well. 11270a9f38efSVaclav Hapla 11280a9f38efSVaclav Hapla Level: advanced 11290a9f38efSVaclav Hapla 11300a9f38efSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5OpenGroup() 11310a9f38efSVaclav Hapla @*/ 11320a9f38efSVaclav Hapla PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, PetscBool *has) 11330a9f38efSVaclav Hapla { 11340a9f38efSVaclav Hapla H5O_type_t type; 11350a9f38efSVaclav Hapla const char *name; 11360a9f38efSVaclav Hapla PetscErrorCode ierr; 11370a9f38efSVaclav Hapla 11380a9f38efSVaclav Hapla PetscFunctionBegin; 11390a9f38efSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 11400a9f38efSVaclav Hapla PetscValidIntPointer(has,2); 11410a9f38efSVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &name);CHKERRQ(ierr); 11420a9f38efSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, name, PETSC_FALSE, has, &type);CHKERRQ(ierr); 11430a9f38efSVaclav Hapla *has = (type == H5O_TYPE_GROUP) ? PETSC_TRUE : PETSC_FALSE; 11440a9f38efSVaclav Hapla PetscFunctionReturn(0); 11450a9f38efSVaclav Hapla } 11460a9f38efSVaclav Hapla 11470a9f38efSVaclav Hapla /*@ 1148e3f143f7SVaclav Hapla PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group 1149ecce1506SVaclav Hapla 1150ecce1506SVaclav Hapla Input Parameters: 1151ecce1506SVaclav Hapla + viewer - The HDF5 viewer 1152ecce1506SVaclav Hapla - obj - The named object 1153ecce1506SVaclav Hapla 1154ecce1506SVaclav Hapla Output Parameter: 1155ecce1506SVaclav Hapla . has - Flag for dataset existence; PETSC_FALSE for unnamed object 1156ecce1506SVaclav Hapla 1157e3f143f7SVaclav Hapla Notes: 1158e3f143f7SVaclav Hapla If the path exists but is not a dataset, this returns PETSC_FALSE as well. 1159e3f143f7SVaclav Hapla 1160ecce1506SVaclav Hapla Level: advanced 1161ecce1506SVaclav Hapla 1162e3f143f7SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1163ecce1506SVaclav Hapla @*/ 1164ecce1506SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has) 1165ecce1506SVaclav Hapla { 116656cc0592SVaclav Hapla H5O_type_t type; 11676c132bc1SVaclav Hapla char *path; 1168ecce1506SVaclav Hapla PetscErrorCode ierr; 1169ecce1506SVaclav Hapla 1170ecce1506SVaclav Hapla PetscFunctionBegin; 1171c57a1a9aSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1172c57a1a9aSVaclav Hapla PetscValidHeader(obj,2); 1173c57a1a9aSVaclav Hapla PetscValidIntPointer(has,3); 1174ecce1506SVaclav Hapla *has = PETSC_FALSE; 1175e3f143f7SVaclav Hapla if (!obj->name) PetscFunctionReturn(0); 11766c132bc1SVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, obj->name, &path);CHKERRQ(ierr); 11776c132bc1SVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, path, PETSC_FALSE, has, &type);CHKERRQ(ierr); 117856cc0592SVaclav Hapla *has = (type == H5O_TYPE_DATASET) ? PETSC_TRUE : PETSC_FALSE; 11796c132bc1SVaclav Hapla ierr = PetscFree(path);CHKERRQ(ierr); 1180ecce1506SVaclav Hapla PetscFunctionReturn(0); 1181ecce1506SVaclav Hapla } 1182ecce1506SVaclav Hapla 1183df863907SAlex Fikl /*@C 1184ece83b6cSVaclav Hapla PetscViewerHDF5HasAttribute - Check whether an attribute exists 1185e7bdbf8eSMatthew G. Knepley 1186e7bdbf8eSMatthew G. Knepley Input Parameters: 1187e7bdbf8eSMatthew G. Knepley + viewer - The HDF5 viewer 118810e69e8fSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 1189e7bdbf8eSMatthew G. Knepley - name - The attribute name 1190e7bdbf8eSMatthew G. Knepley 1191e7bdbf8eSMatthew G. Knepley Output Parameter: 1192e7bdbf8eSMatthew G. Knepley . has - Flag for attribute existence 1193e7bdbf8eSMatthew G. Knepley 1194e7bdbf8eSMatthew G. Knepley Level: advanced 1195e7bdbf8eSMatthew G. Knepley 1196577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1197e7bdbf8eSMatthew G. Knepley @*/ 119810e69e8fSVaclav Hapla PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscBool *has) 1199e7bdbf8eSMatthew G. Knepley { 120010e69e8fSVaclav Hapla char *parent; 1201e7bdbf8eSMatthew G. Knepley PetscErrorCode ierr; 1202e7bdbf8eSMatthew G. Knepley 1203e7bdbf8eSMatthew G. Knepley PetscFunctionBegin; 12045cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 120510e69e8fSVaclav Hapla if (dataset) PetscValidCharPointer(dataset,2); 1206c57a1a9aSVaclav Hapla PetscValidCharPointer(name,3); 1207c57a1a9aSVaclav Hapla PetscValidIntPointer(has,4); 120810e69e8fSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 1209bb286ee1SVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, has, NULL);CHKERRQ(ierr); 121010e69e8fSVaclav Hapla if (*has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, has);CHKERRQ(ierr);} 121110e69e8fSVaclav Hapla ierr = PetscFree(parent);CHKERRQ(ierr); 121206db490cSVaclav Hapla PetscFunctionReturn(0); 121306db490cSVaclav Hapla } 121406db490cSVaclav Hapla 1215577e0e04SVaclav Hapla /*@C 1216577e0e04SVaclav Hapla PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name 1217577e0e04SVaclav Hapla 1218577e0e04SVaclav Hapla Input Parameters: 1219577e0e04SVaclav Hapla + viewer - The HDF5 viewer 1220577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1221577e0e04SVaclav Hapla - name - The attribute name 1222577e0e04SVaclav Hapla 1223577e0e04SVaclav Hapla Output Parameter: 1224577e0e04SVaclav Hapla . has - Flag for attribute existence 1225577e0e04SVaclav Hapla 1226577e0e04SVaclav Hapla Notes: 1227577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1228577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 1229577e0e04SVaclav Hapla 1230577e0e04SVaclav Hapla Level: advanced 1231577e0e04SVaclav Hapla 1232577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1233577e0e04SVaclav Hapla @*/ 1234577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has) 1235577e0e04SVaclav Hapla { 1236577e0e04SVaclav Hapla PetscErrorCode ierr; 1237577e0e04SVaclav Hapla 1238577e0e04SVaclav Hapla PetscFunctionBegin; 1239577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1240577e0e04SVaclav Hapla PetscValidHeader(obj,2); 1241577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 1242577e0e04SVaclav Hapla PetscValidIntPointer(has,4); 1243577e0e04SVaclav Hapla ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 1244577e0e04SVaclav Hapla ierr = PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);CHKERRQ(ierr); 1245577e0e04SVaclav Hapla PetscFunctionReturn(0); 1246577e0e04SVaclav Hapla } 1247577e0e04SVaclav Hapla 124806db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 124906db490cSVaclav Hapla { 125006db490cSVaclav Hapla hid_t h5; 125106db490cSVaclav Hapla htri_t hhas; 125206db490cSVaclav Hapla PetscErrorCode ierr; 125306db490cSVaclav Hapla 125406db490cSVaclav Hapla PetscFunctionBegin; 125506db490cSVaclav Hapla ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 12562f936e54SVaclav Hapla PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT)); 12575cdeb108SMatthew G. Knepley *has = hhas ? PETSC_TRUE : PETSC_FALSE; 1258e7bdbf8eSMatthew G. Knepley PetscFunctionReturn(0); 1259e7bdbf8eSMatthew G. Knepley } 1260e7bdbf8eSMatthew G. Knepley 1261a75e6a4aSMatthew G. Knepley /* 1262a75e6a4aSMatthew G. Knepley The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that 1263a75e6a4aSMatthew G. Knepley is attached to a communicator, in this case the attribute is a PetscViewer. 1264a75e6a4aSMatthew G. Knepley */ 1265d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID; 1266a75e6a4aSMatthew G. Knepley 1267a75e6a4aSMatthew G. Knepley /*@C 1268a75e6a4aSMatthew G. Knepley PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator. 1269a75e6a4aSMatthew G. Knepley 1270d083f849SBarry Smith Collective 1271a75e6a4aSMatthew G. Knepley 1272a75e6a4aSMatthew G. Knepley Input Parameter: 1273a75e6a4aSMatthew G. Knepley . comm - the MPI communicator to share the HDF5 PetscViewer 1274a75e6a4aSMatthew G. Knepley 1275a75e6a4aSMatthew G. Knepley Level: intermediate 1276a75e6a4aSMatthew G. Knepley 1277a75e6a4aSMatthew G. Knepley Options Database Keys: 1278a75e6a4aSMatthew G. Knepley . -viewer_hdf5_filename <name> 1279a75e6a4aSMatthew G. Knepley 1280a75e6a4aSMatthew G. Knepley Environmental variables: 1281a75e6a4aSMatthew G. Knepley . PETSC_VIEWER_HDF5_FILENAME 1282a75e6a4aSMatthew G. Knepley 1283a75e6a4aSMatthew G. Knepley Notes: 1284a75e6a4aSMatthew G. Knepley Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return 1285a75e6a4aSMatthew G. Knepley an error code. The HDF5 PetscViewer is usually used in the form 1286a75e6a4aSMatthew G. Knepley $ XXXView(XXX object, PETSC_VIEWER_HDF5_(comm)); 1287a75e6a4aSMatthew G. Knepley 1288a75e6a4aSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy() 1289a75e6a4aSMatthew G. Knepley @*/ 1290a75e6a4aSMatthew G. Knepley PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm) 1291a75e6a4aSMatthew G. Knepley { 1292a75e6a4aSMatthew G. Knepley PetscErrorCode ierr; 1293a75e6a4aSMatthew G. Knepley PetscBool flg; 1294a75e6a4aSMatthew G. Knepley PetscViewer viewer; 1295a75e6a4aSMatthew G. Knepley char fname[PETSC_MAX_PATH_LEN]; 1296a75e6a4aSMatthew G. Knepley MPI_Comm ncomm; 1297a75e6a4aSMatthew G. Knepley 1298a75e6a4aSMatthew G. Knepley PetscFunctionBegin; 1299908793a3SLisandro 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);} 1300a75e6a4aSMatthew G. Knepley if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) { 1301908793a3SLisandro Dalcin ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,NULL); 1302908793a3SLisandro Dalcin if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 1303a75e6a4aSMatthew G. Knepley } 130447435625SJed Brown ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg); 1305908793a3SLisandro Dalcin if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 1306a75e6a4aSMatthew G. Knepley if (!flg) { /* PetscViewer not yet created */ 1307a75e6a4aSMatthew G. Knepley ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg); 13082cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1309a75e6a4aSMatthew G. Knepley if (!flg) { 1310a75e6a4aSMatthew G. Knepley ierr = PetscStrcpy(fname,"output.h5"); 13112cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1312a75e6a4aSMatthew G. Knepley } 1313a75e6a4aSMatthew G. Knepley ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer); 13142cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1315a75e6a4aSMatthew G. Knepley ierr = PetscObjectRegisterDestroy((PetscObject)viewer); 13162cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 131747435625SJed Brown ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer); 1318908793a3SLisandro Dalcin if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 1319a75e6a4aSMatthew G. Knepley } 1320a75e6a4aSMatthew G. Knepley ierr = PetscCommDestroy(&ncomm); 13212cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1322a75e6a4aSMatthew G. Knepley PetscFunctionReturn(viewer); 1323a75e6a4aSMatthew G. Knepley } 1324