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 */ 278ee8b9732SVaclav Hapla #if H5_VERSION_GE(1,10,3) 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; 286ee8b9732SVaclav Hapla ierr = PetscPrintf(PetscObjectComm((PetscObject)viewer), "Warning: PetscViewerHDF5SetCollective(viewer,PETSC_TRUE) is ignored for HDF5 versions prior to 1.10.3\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. 306ee8b9732SVaclav Hapla However, this works correctly only since HDF5 1.10.3; 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 { 331ee8b9732SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 332ee8b9732SVaclav Hapla H5FD_mpio_xfer_t mode; 333ee8b9732SVaclav Hapla 334ee8b9732SVaclav Hapla PetscFunctionBegin; 335ee8b9732SVaclav Hapla PetscStackCallHDF5(H5Pget_dxpl_mpio,(hdf5->dxpl_id, &mode)); 336ee8b9732SVaclav Hapla *flg = (mode == H5FD_MPIO_COLLECTIVE) ? PETSC_TRUE : PETSC_FALSE; 337ee8b9732SVaclav Hapla PetscFunctionReturn(0); 338ee8b9732SVaclav Hapla } 339ee8b9732SVaclav Hapla 340ee8b9732SVaclav Hapla /*@ 341ee8b9732SVaclav Hapla PetscViewerHDF5GetCollective - Return flag whether collective MPI-IO transfer mode is used for HDF5 reads and writes. 342ee8b9732SVaclav Hapla 343ee8b9732SVaclav Hapla Not Collective 344ee8b9732SVaclav Hapla 345ee8b9732SVaclav Hapla Input Parameters: 346ee8b9732SVaclav Hapla . viewer - the HDF5 PetscViewer 347ee8b9732SVaclav Hapla 348ee8b9732SVaclav Hapla Output Parameters: 349ee8b9732SVaclav Hapla . flg - the flag 350ee8b9732SVaclav Hapla 351ee8b9732SVaclav Hapla Level: intermediate 352ee8b9732SVaclav Hapla 353ee8b9732SVaclav Hapla Notes: 354ee8b9732SVaclav Hapla This setting works correctly only since HDF5 1.10.3. For older versions, PETSC_FALSE will be always returned. 355ee8b9732SVaclav Hapla For more details, see PetscViewerHDF5SetCollective(). 356ee8b9732SVaclav Hapla 357ee8b9732SVaclav Hapla .seealso: PetscViewerHDF5SetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open() 358ee8b9732SVaclav Hapla 359ee8b9732SVaclav Hapla @*/ 360ee8b9732SVaclav Hapla PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer viewer,PetscBool *flg) 361ee8b9732SVaclav Hapla { 362ee8b9732SVaclav Hapla PetscErrorCode ierr; 363ee8b9732SVaclav Hapla 364ee8b9732SVaclav Hapla PetscFunctionBegin; 365ee8b9732SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 366534a8f05SLisandro Dalcin PetscValidBoolPointer(flg,2); 367ee8b9732SVaclav Hapla 368ee8b9732SVaclav Hapla ierr = PetscUseMethod(viewer,"PetscViewerHDF5GetCollective_C",(PetscViewer,PetscBool*),(viewer,flg));CHKERRQ(ierr); 369ee8b9732SVaclav Hapla PetscFunctionReturn(0); 370ee8b9732SVaclav Hapla } 371ee8b9732SVaclav Hapla 3729b1403beSVaclav Hapla static PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[]) 3735c6c1daeSBarry Smith { 3745c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 3755c6c1daeSBarry Smith hid_t plist_id; 3765c6c1daeSBarry Smith PetscErrorCode ierr; 3775c6c1daeSBarry Smith 3785c6c1daeSBarry Smith PetscFunctionBegin; 379f683cc58SToby Isaac if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 380f683cc58SToby Isaac if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);} 3815c6c1daeSBarry Smith ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr); 3825c6c1daeSBarry Smith /* Set up file access property list with parallel I/O access */ 383729ab6d9SBarry Smith PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS)); 384d28bfa55SVaclav Hapla PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), MPI_INFO_NULL)); 3855c6c1daeSBarry Smith /* Create or open the file collectively */ 3865c6c1daeSBarry Smith switch (hdf5->btype) { 3875c6c1daeSBarry Smith case FILE_MODE_READ: 388729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id)); 3895c6c1daeSBarry Smith break; 3905c6c1daeSBarry Smith case FILE_MODE_APPEND: 391729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id)); 3925c6c1daeSBarry Smith break; 3935c6c1daeSBarry Smith case FILE_MODE_WRITE: 394729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id)); 3955c6c1daeSBarry Smith break; 3965c6c1daeSBarry Smith default: 3975c6c1daeSBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 3985c6c1daeSBarry Smith } 3995c6c1daeSBarry Smith if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name); 400729ab6d9SBarry Smith PetscStackCallHDF5(H5Pclose,(plist_id)); 4015c6c1daeSBarry Smith PetscFunctionReturn(0); 4025c6c1daeSBarry Smith } 4035c6c1daeSBarry Smith 404d1232d7fSToby Isaac static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name) 405d1232d7fSToby Isaac { 406d1232d7fSToby Isaac PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data; 407d1232d7fSToby Isaac 408d1232d7fSToby Isaac PetscFunctionBegin; 409d1232d7fSToby Isaac *name = vhdf5->filename; 410d1232d7fSToby Isaac PetscFunctionReturn(0); 411d1232d7fSToby Isaac } 412d1232d7fSToby Isaac 413b723ab35SVaclav Hapla static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer) 414b723ab35SVaclav Hapla { 4155dc64a97SVaclav Hapla /* 416b723ab35SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 417b723ab35SVaclav Hapla PetscErrorCode ierr; 4185dc64a97SVaclav Hapla */ 419b723ab35SVaclav Hapla 420b723ab35SVaclav Hapla PetscFunctionBegin; 421b723ab35SVaclav Hapla PetscFunctionReturn(0); 422b723ab35SVaclav Hapla } 423b723ab35SVaclav Hapla 4248556b5ebSBarry Smith /*MC 4258556b5ebSBarry Smith PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file 4268556b5ebSBarry Smith 4278556b5ebSBarry Smith 4288556b5ebSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET, 4298556b5ebSBarry Smith PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING, 4308556b5ebSBarry Smith PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB, 4318556b5ebSBarry Smith PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType() 4328556b5ebSBarry Smith 4331b266c99SBarry Smith Level: beginner 4348556b5ebSBarry Smith M*/ 435d1232d7fSToby Isaac 4368cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v) 4375c6c1daeSBarry Smith { 4385c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5; 4395c6c1daeSBarry Smith PetscErrorCode ierr; 4405c6c1daeSBarry Smith 4415c6c1daeSBarry Smith PetscFunctionBegin; 442b00a9115SJed Brown ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr); 4435c6c1daeSBarry Smith 4445c6c1daeSBarry Smith v->data = (void*) hdf5; 4455c6c1daeSBarry Smith v->ops->destroy = PetscViewerDestroy_HDF5; 44682ea9b62SBarry Smith v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5; 447b723ab35SVaclav Hapla v->ops->setup = PetscViewerSetUp_HDF5; 4481b793a25SVaclav Hapla v->ops->view = PetscViewerView_HDF5; 4495c6c1daeSBarry Smith v->ops->flush = 0; 4505c6c1daeSBarry Smith hdf5->btype = (PetscFileMode) -1; 4515c6c1daeSBarry Smith hdf5->filename = 0; 4525c6c1daeSBarry Smith hdf5->timestep = -1; 4530298fd71SBarry Smith hdf5->groups = NULL; 4545c6c1daeSBarry Smith 4559c5072fbSVaclav Hapla PetscStackCallHDF5Return(hdf5->dxpl_id,H5Pcreate,(H5P_DATASET_XFER)); 4569c5072fbSVaclav Hapla 4570b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr); 458d1232d7fSToby Isaac ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr); 4590b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr); 4600b62783dSVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_HDF5);CHKERRQ(ierr); 46182ea9b62SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr); 4629a0502c6SHåkon Strandenes ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr); 463ee8b9732SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetCollective_C",PetscViewerHDF5SetCollective_HDF5);CHKERRQ(ierr); 464ee8b9732SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetCollective_C",PetscViewerHDF5GetCollective_HDF5);CHKERRQ(ierr); 4655c6c1daeSBarry Smith PetscFunctionReturn(0); 4665c6c1daeSBarry Smith } 4675c6c1daeSBarry Smith 4685c6c1daeSBarry Smith /*@C 4695c6c1daeSBarry Smith PetscViewerHDF5Open - Opens a file for HDF5 input/output. 4705c6c1daeSBarry Smith 471d083f849SBarry Smith Collective 4725c6c1daeSBarry Smith 4735c6c1daeSBarry Smith Input Parameters: 4745c6c1daeSBarry Smith + comm - MPI communicator 4755c6c1daeSBarry Smith . name - name of file 4765c6c1daeSBarry Smith - type - type of file 4775c6c1daeSBarry Smith $ FILE_MODE_WRITE - create new file for binary output 4785c6c1daeSBarry Smith $ FILE_MODE_READ - open existing file for binary input 4795c6c1daeSBarry Smith $ FILE_MODE_APPEND - open existing file for binary output 4805c6c1daeSBarry Smith 4815c6c1daeSBarry Smith Output Parameter: 4825c6c1daeSBarry Smith . hdf5v - PetscViewer for HDF5 input/output to use with the specified file 4835c6c1daeSBarry Smith 48482ea9b62SBarry Smith Options Database: 485a2b725a8SWilliam 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 486a2b725a8SWilliam Gropp - -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal 48782ea9b62SBarry Smith 4885c6c1daeSBarry Smith Level: beginner 4895c6c1daeSBarry Smith 4905c6c1daeSBarry Smith Note: 4915c6c1daeSBarry Smith This PetscViewer should be destroyed with PetscViewerDestroy(). 4925c6c1daeSBarry Smith 4935c6c1daeSBarry Smith 4946a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(), 4959a0502c6SHåkon Strandenes PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(), 496a56f64adSBarry Smith MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName() 4975c6c1daeSBarry Smith @*/ 4985c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v) 4995c6c1daeSBarry Smith { 5005c6c1daeSBarry Smith PetscErrorCode ierr; 5015c6c1daeSBarry Smith 5025c6c1daeSBarry Smith PetscFunctionBegin; 5035c6c1daeSBarry Smith ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr); 5045c6c1daeSBarry Smith ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr); 5055c6c1daeSBarry Smith ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr); 5065c6c1daeSBarry Smith ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr); 507b663d2d2SVaclav Hapla ierr = PetscViewerSetFromOptions(*hdf5v);CHKERRQ(ierr); 5085c6c1daeSBarry Smith PetscFunctionReturn(0); 5095c6c1daeSBarry Smith } 5105c6c1daeSBarry Smith 5115c6c1daeSBarry Smith /*@C 5125c6c1daeSBarry Smith PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls 5135c6c1daeSBarry Smith 5145c6c1daeSBarry Smith Not collective 5155c6c1daeSBarry Smith 5165c6c1daeSBarry Smith Input Parameter: 5175c6c1daeSBarry Smith . viewer - the PetscViewer 5185c6c1daeSBarry Smith 5195c6c1daeSBarry Smith Output Parameter: 5205c6c1daeSBarry Smith . file_id - The file id 5215c6c1daeSBarry Smith 5225c6c1daeSBarry Smith Level: intermediate 5235c6c1daeSBarry Smith 5245c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open() 5255c6c1daeSBarry Smith @*/ 5265c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id) 5275c6c1daeSBarry Smith { 5285c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 5295c6c1daeSBarry Smith 5305c6c1daeSBarry Smith PetscFunctionBegin; 5315c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 5325c6c1daeSBarry Smith if (file_id) *file_id = hdf5->file_id; 5335c6c1daeSBarry Smith PetscFunctionReturn(0); 5345c6c1daeSBarry Smith } 5355c6c1daeSBarry Smith 5365c6c1daeSBarry Smith /*@C 5375c6c1daeSBarry Smith PetscViewerHDF5PushGroup - Set the current HDF5 group for output 5385c6c1daeSBarry Smith 5395c6c1daeSBarry Smith Not collective 5405c6c1daeSBarry Smith 5415c6c1daeSBarry Smith Input Parameters: 5425c6c1daeSBarry Smith + viewer - the PetscViewer 5435c6c1daeSBarry Smith - name - The group name 5445c6c1daeSBarry Smith 5455c6c1daeSBarry Smith Level: intermediate 5465c6c1daeSBarry Smith 547*820fc2d1SVaclav Hapla Note: The group name being NULL, empty string, or a sequence of all slashes (e.g. "///") is always internally stored as NULL and interpreted as "/". 548*820fc2d1SVaclav Hapla 549874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 5505c6c1daeSBarry Smith @*/ 551be7aa430SVaclav Hapla PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char name[]) 5525c6c1daeSBarry Smith { 5535c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 5547d964842SVaclav Hapla PetscViewerHDF5GroupList *groupNode; 5555c6c1daeSBarry Smith PetscErrorCode ierr; 5565c6c1daeSBarry Smith 5575c6c1daeSBarry Smith PetscFunctionBegin; 5585c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 559*820fc2d1SVaclav Hapla if (name) PetscValidCharPointer(name,2); 560*820fc2d1SVaclav Hapla if (name && name[0]) { 561*820fc2d1SVaclav Hapla size_t i,len; 562*820fc2d1SVaclav Hapla ierr = PetscStrlen(name, &len);CHKERRQ(ierr); 563*820fc2d1SVaclav Hapla for (i=0; i<len; i++) if (name[i] != '/') break; 564*820fc2d1SVaclav Hapla if (i == len) name = NULL; 565*820fc2d1SVaclav Hapla } else name = NULL; 56695dccacaSBarry Smith ierr = PetscNew(&groupNode);CHKERRQ(ierr); 5675c6c1daeSBarry Smith ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr); 5685c6c1daeSBarry Smith groupNode->next = hdf5->groups; 5695c6c1daeSBarry Smith hdf5->groups = groupNode; 5705c6c1daeSBarry Smith PetscFunctionReturn(0); 5715c6c1daeSBarry Smith } 5725c6c1daeSBarry Smith 5733ef9c667SSatish Balay /*@ 5745c6c1daeSBarry Smith PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value 5755c6c1daeSBarry Smith 5765c6c1daeSBarry Smith Not collective 5775c6c1daeSBarry Smith 5785c6c1daeSBarry Smith Input Parameter: 5795c6c1daeSBarry Smith . viewer - the PetscViewer 5805c6c1daeSBarry Smith 5815c6c1daeSBarry Smith Level: intermediate 5825c6c1daeSBarry Smith 583874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 5845c6c1daeSBarry Smith @*/ 5855c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer) 5865c6c1daeSBarry Smith { 5875c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 5887d964842SVaclav Hapla PetscViewerHDF5GroupList *groupNode; 5895c6c1daeSBarry Smith PetscErrorCode ierr; 5905c6c1daeSBarry Smith 5915c6c1daeSBarry Smith PetscFunctionBegin; 5925c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 59382f516ccSBarry Smith if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop"); 5945c6c1daeSBarry Smith groupNode = hdf5->groups; 5955c6c1daeSBarry Smith hdf5->groups = hdf5->groups->next; 5965c6c1daeSBarry Smith ierr = PetscFree(groupNode->name);CHKERRQ(ierr); 5975c6c1daeSBarry Smith ierr = PetscFree(groupNode);CHKERRQ(ierr); 5985c6c1daeSBarry Smith PetscFunctionReturn(0); 5995c6c1daeSBarry Smith } 6005c6c1daeSBarry Smith 6015c6c1daeSBarry Smith /*@C 602874270d9SVaclav Hapla PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup(). 603874270d9SVaclav Hapla If none has been assigned, returns NULL. 6045c6c1daeSBarry Smith 6055c6c1daeSBarry Smith Not collective 6065c6c1daeSBarry Smith 6075c6c1daeSBarry Smith Input Parameter: 6085c6c1daeSBarry Smith . viewer - the PetscViewer 6095c6c1daeSBarry Smith 6105c6c1daeSBarry Smith Output Parameter: 6115c6c1daeSBarry Smith . name - The group name 6125c6c1daeSBarry Smith 6135c6c1daeSBarry Smith Level: intermediate 6145c6c1daeSBarry Smith 615874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup() 6165c6c1daeSBarry Smith @*/ 617be7aa430SVaclav Hapla PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char *name[]) 6185c6c1daeSBarry Smith { 6195c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 6205c6c1daeSBarry Smith 6215c6c1daeSBarry Smith PetscFunctionBegin; 6225c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 623c959eef4SJed Brown PetscValidPointer(name,2); 624a297a907SKarl Rupp if (hdf5->groups) *name = hdf5->groups->name; 6250298fd71SBarry Smith else *name = NULL; 6265c6c1daeSBarry Smith PetscFunctionReturn(0); 6275c6c1daeSBarry Smith } 6285c6c1daeSBarry Smith 6298c557b5aSVaclav Hapla /*@ 630874270d9SVaclav Hapla PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(), 631874270d9SVaclav Hapla and return this group's ID and file ID. 632874270d9SVaclav Hapla If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID. 633874270d9SVaclav Hapla 634874270d9SVaclav Hapla Not collective 635874270d9SVaclav Hapla 636874270d9SVaclav Hapla Input Parameter: 637874270d9SVaclav Hapla . viewer - the PetscViewer 638874270d9SVaclav Hapla 639874270d9SVaclav Hapla Output Parameter: 640874270d9SVaclav Hapla + fileId - The HDF5 file ID 641874270d9SVaclav Hapla - groupId - The HDF5 group ID 642874270d9SVaclav Hapla 643e74616beSVaclav Hapla Notes: 644e74616beSVaclav Hapla If the viewer is writable, the group is created if it doesn't exist yet. 645e74616beSVaclav Hapla 646874270d9SVaclav Hapla Level: intermediate 647874270d9SVaclav Hapla 648874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 649874270d9SVaclav Hapla @*/ 65054dbf706SVaclav Hapla PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId) 65154dbf706SVaclav Hapla { 65290e03537SVaclav Hapla hid_t file_id; 65390e03537SVaclav Hapla H5O_type_t type; 65454dbf706SVaclav Hapla const char *groupName = NULL; 655e74616beSVaclav Hapla PetscBool create; 65654dbf706SVaclav Hapla PetscErrorCode ierr; 65754dbf706SVaclav Hapla 65854dbf706SVaclav Hapla PetscFunctionBegin; 659e74616beSVaclav Hapla ierr = PetscViewerWritable(viewer, &create);CHKERRQ(ierr); 66054dbf706SVaclav Hapla ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr); 66154dbf706SVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr); 662e74616beSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, groupName, create, NULL, &type);CHKERRQ(ierr); 66390e03537SVaclav 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); 66490e03537SVaclav Hapla PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT)); 66554dbf706SVaclav Hapla *fileId = file_id; 66654dbf706SVaclav Hapla PetscFunctionReturn(0); 66754dbf706SVaclav Hapla } 66854dbf706SVaclav Hapla 6693ef9c667SSatish Balay /*@ 6705c6c1daeSBarry Smith PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time. 6715c6c1daeSBarry Smith 6725c6c1daeSBarry Smith Not collective 6735c6c1daeSBarry Smith 6745c6c1daeSBarry Smith Input Parameter: 6755c6c1daeSBarry Smith . viewer - the PetscViewer 6765c6c1daeSBarry Smith 6775c6c1daeSBarry Smith Level: intermediate 6785c6c1daeSBarry Smith 6795c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep() 6805c6c1daeSBarry Smith @*/ 6815c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer) 6825c6c1daeSBarry Smith { 6835c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 6845c6c1daeSBarry Smith 6855c6c1daeSBarry Smith PetscFunctionBegin; 6865c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 6875c6c1daeSBarry Smith ++hdf5->timestep; 6885c6c1daeSBarry Smith PetscFunctionReturn(0); 6895c6c1daeSBarry Smith } 6905c6c1daeSBarry Smith 6913ef9c667SSatish Balay /*@ 6925c6c1daeSBarry Smith PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep 6935c6c1daeSBarry Smith of -1 disables blocking with timesteps. 6945c6c1daeSBarry Smith 6955c6c1daeSBarry Smith Not collective 6965c6c1daeSBarry Smith 6975c6c1daeSBarry Smith Input Parameters: 6985c6c1daeSBarry Smith + viewer - the PetscViewer 6995c6c1daeSBarry Smith - timestep - The timestep number 7005c6c1daeSBarry Smith 7015c6c1daeSBarry Smith Level: intermediate 7025c6c1daeSBarry Smith 7035c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 7045c6c1daeSBarry Smith @*/ 7055c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep) 7065c6c1daeSBarry Smith { 7075c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 7085c6c1daeSBarry Smith 7095c6c1daeSBarry Smith PetscFunctionBegin; 7105c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 7115c6c1daeSBarry Smith hdf5->timestep = timestep; 7125c6c1daeSBarry Smith PetscFunctionReturn(0); 7135c6c1daeSBarry Smith } 7145c6c1daeSBarry Smith 7153ef9c667SSatish Balay /*@ 7165c6c1daeSBarry Smith PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time. 7175c6c1daeSBarry Smith 7185c6c1daeSBarry Smith Not collective 7195c6c1daeSBarry Smith 7205c6c1daeSBarry Smith Input Parameter: 7215c6c1daeSBarry Smith . viewer - the PetscViewer 7225c6c1daeSBarry Smith 7235c6c1daeSBarry Smith Output Parameter: 7245c6c1daeSBarry Smith . timestep - The timestep number 7255c6c1daeSBarry Smith 7265c6c1daeSBarry Smith Level: intermediate 7275c6c1daeSBarry Smith 7285c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep() 7295c6c1daeSBarry Smith @*/ 7305c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep) 7315c6c1daeSBarry Smith { 7325c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 7335c6c1daeSBarry Smith 7345c6c1daeSBarry Smith PetscFunctionBegin; 7355c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 7365c6c1daeSBarry Smith PetscValidPointer(timestep,2); 7375c6c1daeSBarry Smith *timestep = hdf5->timestep; 7385c6c1daeSBarry Smith PetscFunctionReturn(0); 7395c6c1daeSBarry Smith } 7405c6c1daeSBarry Smith 74136bce6e8SMatthew G. Knepley /*@C 74236bce6e8SMatthew G. Knepley PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name. 74336bce6e8SMatthew G. Knepley 74436bce6e8SMatthew G. Knepley Not collective 74536bce6e8SMatthew G. Knepley 74636bce6e8SMatthew G. Knepley Input Parameter: 74736bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 74836bce6e8SMatthew G. Knepley 74936bce6e8SMatthew G. Knepley Output Parameter: 75036bce6e8SMatthew G. Knepley . mtype - the MPI datatype (for example MPI_DOUBLE, ...) 75136bce6e8SMatthew G. Knepley 75236bce6e8SMatthew G. Knepley Level: advanced 75336bce6e8SMatthew G. Knepley 75436bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 75536bce6e8SMatthew G. Knepley @*/ 75636bce6e8SMatthew G. Knepley PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype) 75736bce6e8SMatthew G. Knepley { 75836bce6e8SMatthew G. Knepley PetscFunctionBegin; 75936bce6e8SMatthew G. Knepley if (ptype == PETSC_INT) 76036bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 76136bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_LLONG; 76236bce6e8SMatthew G. Knepley #else 76336bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_INT; 76436bce6e8SMatthew G. Knepley #endif 76536bce6e8SMatthew G. Knepley else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE; 76636bce6e8SMatthew G. Knepley else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG; 76736bce6e8SMatthew G. Knepley else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT; 768c9450970SBarry Smith else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_INT; 769de3ba262SVaclav Hapla else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_INT; 77036bce6e8SMatthew G. Knepley else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT; 77136bce6e8SMatthew G. Knepley else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR; 77236bce6e8SMatthew G. Knepley else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR; 7737e97c476SMatthew G. Knepley else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1); 77436bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype"); 77536bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 77636bce6e8SMatthew G. Knepley } 77736bce6e8SMatthew G. Knepley 77836bce6e8SMatthew G. Knepley /*@C 77936bce6e8SMatthew G. Knepley PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name 78036bce6e8SMatthew G. Knepley 78136bce6e8SMatthew G. Knepley Not collective 78236bce6e8SMatthew G. Knepley 78336bce6e8SMatthew G. Knepley Input Parameter: 78436bce6e8SMatthew G. Knepley . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...) 78536bce6e8SMatthew G. Knepley 78636bce6e8SMatthew G. Knepley Output Parameter: 78736bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 78836bce6e8SMatthew G. Knepley 78936bce6e8SMatthew G. Knepley Level: advanced 79036bce6e8SMatthew G. Knepley 79136bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 79236bce6e8SMatthew G. Knepley @*/ 79336bce6e8SMatthew G. Knepley PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype) 79436bce6e8SMatthew G. Knepley { 79536bce6e8SMatthew G. Knepley PetscFunctionBegin; 79636bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 79736bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG; 79836bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT; 79936bce6e8SMatthew G. Knepley #else 80036bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT; 80136bce6e8SMatthew G. Knepley #endif 80236bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE; 80336bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG; 80436bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT; 80536bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT; 80636bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR; 80736bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR; 8087e97c476SMatthew G. Knepley else if (htype == H5T_C_S1) *ptype = PETSC_STRING; 80936bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype"); 81036bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 81136bce6e8SMatthew G. Knepley } 81236bce6e8SMatthew G. Knepley 813df863907SAlex Fikl /*@C 814b8ee85a1SVaclav Hapla PetscViewerHDF5WriteAttribute - Write an attribute 81536bce6e8SMatthew G. Knepley 81636bce6e8SMatthew G. Knepley Input Parameters: 81736bce6e8SMatthew G. Knepley + viewer - The HDF5 viewer 81857d22f7dSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 81936bce6e8SMatthew G. Knepley . name - The attribute name 82036bce6e8SMatthew G. Knepley . datatype - The attribute type 82136bce6e8SMatthew G. Knepley - value - The attribute value 82236bce6e8SMatthew G. Knepley 82336bce6e8SMatthew G. Knepley Level: advanced 82436bce6e8SMatthew G. Knepley 825577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 82636bce6e8SMatthew G. Knepley @*/ 82757d22f7dSVaclav Hapla PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, const void *value) 82836bce6e8SMatthew G. Knepley { 82957d22f7dSVaclav Hapla char *parent; 83060568592SMatthew G. Knepley hid_t h5, dataspace, obj, attribute, dtype; 831080f144cSVaclav Hapla PetscBool has; 83236bce6e8SMatthew G. Knepley PetscErrorCode ierr; 83336bce6e8SMatthew G. Knepley 83436bce6e8SMatthew G. Knepley PetscFunctionBegin; 8355cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 83657d22f7dSVaclav Hapla if (dataset) PetscValidCharPointer(dataset, 2); 837c57a1a9aSVaclav Hapla PetscValidCharPointer(name, 3); 838b7e81f0fSVaclav Hapla PetscValidPointer(value, 5); 83957d22f7dSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 840bb286ee1SVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_TRUE, NULL, NULL);CHKERRQ(ierr); 841b67695ecSVaclav Hapla ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr); 84236bce6e8SMatthew G. Knepley ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 8437e97c476SMatthew G. Knepley if (datatype == PETSC_STRING) { 8447e97c476SMatthew G. Knepley size_t len; 8457e97c476SMatthew G. Knepley ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr); 846729ab6d9SBarry Smith PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 8477e97c476SMatthew G. Knepley } 84836bce6e8SMatthew G. Knepley ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 849729ab6d9SBarry Smith PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR)); 85060568592SMatthew G. Knepley PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT)); 851080f144cSVaclav Hapla if (has) { 852080f144cSVaclav Hapla PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 853080f144cSVaclav Hapla } else { 85460568592SMatthew G. Knepley PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT)); 855080f144cSVaclav Hapla } 856729ab6d9SBarry Smith PetscStackCallHDF5(H5Awrite,(attribute, dtype, value)); 857729ab6d9SBarry Smith if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype)); 858729ab6d9SBarry Smith PetscStackCallHDF5(H5Aclose,(attribute)); 85960568592SMatthew G. Knepley PetscStackCallHDF5(H5Oclose,(obj)); 860729ab6d9SBarry Smith PetscStackCallHDF5(H5Sclose,(dataspace)); 86157d22f7dSVaclav Hapla ierr = PetscFree(parent);CHKERRQ(ierr); 86236bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 86336bce6e8SMatthew G. Knepley } 86436bce6e8SMatthew G. Knepley 865df863907SAlex Fikl /*@C 866577e0e04SVaclav Hapla PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name 867577e0e04SVaclav Hapla 868577e0e04SVaclav Hapla Input Parameters: 869577e0e04SVaclav Hapla + viewer - The HDF5 viewer 870577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 871577e0e04SVaclav Hapla . name - The attribute name 872577e0e04SVaclav Hapla . datatype - The attribute type 873577e0e04SVaclav Hapla - value - The attribute value 874577e0e04SVaclav Hapla 875577e0e04SVaclav Hapla Notes: 876577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 877577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 878577e0e04SVaclav Hapla 879577e0e04SVaclav Hapla Level: advanced 880577e0e04SVaclav Hapla 881577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 882577e0e04SVaclav Hapla @*/ 883577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value) 884577e0e04SVaclav Hapla { 885577e0e04SVaclav Hapla PetscErrorCode ierr; 886577e0e04SVaclav Hapla 887577e0e04SVaclav Hapla PetscFunctionBegin; 888577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 889577e0e04SVaclav Hapla PetscValidHeader(obj,2); 890577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 891b7e81f0fSVaclav Hapla PetscValidPointer(value,5); 892577e0e04SVaclav Hapla ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 893577e0e04SVaclav Hapla ierr = PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr); 894577e0e04SVaclav Hapla PetscFunctionReturn(0); 895577e0e04SVaclav Hapla } 896577e0e04SVaclav Hapla 897577e0e04SVaclav Hapla /*@C 898b8ee85a1SVaclav Hapla PetscViewerHDF5ReadAttribute - Read an attribute 899ad0c61feSMatthew G. Knepley 900ad0c61feSMatthew G. Knepley Input Parameters: 901ad0c61feSMatthew G. Knepley + viewer - The HDF5 viewer 90257d22f7dSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 903ad0c61feSMatthew G. Knepley . name - The attribute name 904ad0c61feSMatthew G. Knepley - datatype - The attribute type 905ad0c61feSMatthew G. Knepley 906ad0c61feSMatthew G. Knepley Output Parameter: 907ad0c61feSMatthew G. Knepley . value - The attribute value 908ad0c61feSMatthew G. Knepley 909ad0c61feSMatthew G. Knepley Level: advanced 910ad0c61feSMatthew G. Knepley 911577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 912ad0c61feSMatthew G. Knepley @*/ 91357d22f7dSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, void *value) 914ad0c61feSMatthew G. Knepley { 91557d22f7dSVaclav Hapla char *parent; 916f0b58479SMatthew G. Knepley hid_t h5, obj, attribute, atype, dtype; 917969748fdSVaclav Hapla PetscBool has; 918ad0c61feSMatthew G. Knepley PetscErrorCode ierr; 919ad0c61feSMatthew G. Knepley 920ad0c61feSMatthew G. Knepley PetscFunctionBegin; 9215cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 92257d22f7dSVaclav Hapla if (dataset) PetscValidCharPointer(dataset, 2); 923c57a1a9aSVaclav Hapla PetscValidCharPointer(name, 3); 924b7e81f0fSVaclav Hapla PetscValidPointer(value, 5); 92557d22f7dSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 926969748fdSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, &has, NULL);CHKERRQ(ierr); 927969748fdSVaclav Hapla if (has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr);} 928969748fdSVaclav Hapla if (!has) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist", parent, name); 929ad0c61feSMatthew G. Knepley ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 930ad0c61feSMatthew G. Knepley ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 93160568592SMatthew G. Knepley PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT)); 93260568592SMatthew G. Knepley PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 933f0b58479SMatthew G. Knepley if (datatype == PETSC_STRING) { 934f0b58479SMatthew G. Knepley size_t len; 935e90c5075SMarek Pecha PetscStackCallHDF5Return(atype,H5Aget_type,(attribute)); 93615b861d2SVaclav Hapla PetscStackCall("H5Tget_size",len = H5Tget_size(atype)); 937f0b58479SMatthew G. Knepley PetscStackCallHDF5(H5Tclose,(atype)); 938f0b58479SMatthew G. Knepley ierr = PetscMalloc((len+1) * sizeof(char *), &value);CHKERRQ(ierr); 939f0b58479SMatthew G. Knepley } 94070efba33SBarry Smith PetscStackCallHDF5(H5Aread,(attribute, dtype, value)); 941729ab6d9SBarry Smith PetscStackCallHDF5(H5Aclose,(attribute)); 942e90c5075SMarek Pecha /* H5Oclose can be used to close groups, datasets, or committed datatypes */ 943e90c5075SMarek Pecha PetscStackCallHDF5(H5Oclose,(obj)); 94457d22f7dSVaclav Hapla ierr = PetscFree(parent);CHKERRQ(ierr); 945ad0c61feSMatthew G. Knepley PetscFunctionReturn(0); 946ad0c61feSMatthew G. Knepley } 947ad0c61feSMatthew G. Knepley 948577e0e04SVaclav Hapla /*@C 949577e0e04SVaclav Hapla PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name 950577e0e04SVaclav Hapla 951577e0e04SVaclav Hapla Input Parameters: 952577e0e04SVaclav Hapla + viewer - The HDF5 viewer 953577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 954577e0e04SVaclav Hapla . name - The attribute name 955577e0e04SVaclav Hapla - datatype - The attribute type 956577e0e04SVaclav Hapla 957577e0e04SVaclav Hapla Output Parameter: 958577e0e04SVaclav Hapla . value - The attribute value 959577e0e04SVaclav Hapla 960577e0e04SVaclav Hapla Notes: 961577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 962577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 963577e0e04SVaclav Hapla 964577e0e04SVaclav Hapla Level: advanced 965577e0e04SVaclav Hapla 966577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute() PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 967577e0e04SVaclav Hapla @*/ 968577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *value) 969577e0e04SVaclav Hapla { 970577e0e04SVaclav Hapla PetscErrorCode ierr; 971577e0e04SVaclav Hapla 972577e0e04SVaclav Hapla PetscFunctionBegin; 973577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 974577e0e04SVaclav Hapla PetscValidHeader(obj,2); 975577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 976b7e81f0fSVaclav Hapla PetscValidPointer(value, 5); 977577e0e04SVaclav Hapla ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 978577e0e04SVaclav Hapla ierr = PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr); 979577e0e04SVaclav Hapla PetscFunctionReturn(0); 980577e0e04SVaclav Hapla } 981577e0e04SVaclav Hapla 982a75c3fd4SVaclav Hapla PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_) 983a75c3fd4SVaclav Hapla { 984a75c3fd4SVaclav Hapla htri_t exists; 985a75c3fd4SVaclav Hapla hid_t group; 986a75c3fd4SVaclav Hapla 987a75c3fd4SVaclav Hapla PetscFunctionBegin; 988a75c3fd4SVaclav Hapla PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT)); 989a75c3fd4SVaclav Hapla if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT)); 990a75c3fd4SVaclav Hapla if (!exists && createGroup) { 991a75c3fd4SVaclav Hapla PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT)); 992a75c3fd4SVaclav Hapla PetscStackCallHDF5(H5Gclose,(group)); 993a75c3fd4SVaclav Hapla exists = PETSC_TRUE; 994a75c3fd4SVaclav Hapla } 995a75c3fd4SVaclav Hapla *exists_ = (PetscBool) exists; 996a75c3fd4SVaclav Hapla PetscFunctionReturn(0); 997a75c3fd4SVaclav Hapla } 998a75c3fd4SVaclav Hapla 999bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype) 10005cdeb108SMatthew G. Knepley { 100190e03537SVaclav Hapla const char rootGroupName[] = "/"; 10025cdeb108SMatthew G. Knepley hid_t h5; 1003e5bf9ebcSVaclav Hapla PetscBool exists=PETSC_FALSE; 100415b861d2SVaclav Hapla PetscInt i; 100515b861d2SVaclav Hapla int n; 100685688ae2SVaclav Hapla char **hierarchy; 100785688ae2SVaclav Hapla char buf[PETSC_MAX_PATH_LEN]=""; 10085cdeb108SMatthew G. Knepley PetscErrorCode ierr; 10095cdeb108SMatthew G. Knepley 10105cdeb108SMatthew G. Knepley PetscFunctionBegin; 10115cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 101290e03537SVaclav Hapla if (name) PetscValidCharPointer(name, 2); 101390e03537SVaclav Hapla else name = rootGroupName; 1014ccd66a83SVaclav Hapla if (has) { 101556cc0592SVaclav Hapla PetscValidIntPointer(has, 3); 10165cdeb108SMatthew G. Knepley *has = PETSC_FALSE; 1017ccd66a83SVaclav Hapla } 1018ccd66a83SVaclav Hapla if (otype) { 1019ccd66a83SVaclav Hapla PetscValidIntPointer(otype, 4); 102056cc0592SVaclav Hapla *otype = H5O_TYPE_UNKNOWN; 1021ccd66a83SVaclav Hapla } 1022c0cd0301SJed Brown ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 102385688ae2SVaclav Hapla 102485688ae2SVaclav Hapla /* 102585688ae2SVaclav Hapla Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing. 102685688ae2SVaclav Hapla Hence, each of them needs to be tested separately: 102785688ae2SVaclav Hapla 1) whether it's a valid link 102885688ae2SVaclav Hapla 2) whether this link resolves to an object 102985688ae2SVaclav Hapla See H5Oexists_by_name() documentation. 103085688ae2SVaclav Hapla */ 103185688ae2SVaclav Hapla ierr = PetscStrToArray(name,'/',&n,&hierarchy);CHKERRQ(ierr); 103285688ae2SVaclav Hapla if (!n) { 103385688ae2SVaclav Hapla /* Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */ 1034ccd66a83SVaclav Hapla if (has) *has = PETSC_TRUE; 1035ccd66a83SVaclav Hapla if (otype) *otype = H5O_TYPE_GROUP; 103685688ae2SVaclav Hapla ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 103785688ae2SVaclav Hapla PetscFunctionReturn(0); 103885688ae2SVaclav Hapla } 103985688ae2SVaclav Hapla for (i=0; i<n; i++) { 104085688ae2SVaclav Hapla ierr = PetscStrcat(buf,"/");CHKERRQ(ierr); 104185688ae2SVaclav Hapla ierr = PetscStrcat(buf,hierarchy[i]);CHKERRQ(ierr); 1042a75c3fd4SVaclav Hapla ierr = PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);CHKERRQ(ierr); 1043a75c3fd4SVaclav Hapla if (!exists) break; 104485688ae2SVaclav Hapla } 104585688ae2SVaclav Hapla ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 104685688ae2SVaclav Hapla 104785688ae2SVaclav Hapla /* If the object exists, get its type */ 1048ccd66a83SVaclav Hapla if (exists && otype) { 10495cdeb108SMatthew G. Knepley H5O_info_t info; 10505cdeb108SMatthew G. Knepley 105176276748SVaclav Hapla /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */ 105204633f7fSVaclav Hapla PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT)); 105356cc0592SVaclav Hapla *otype = info.type; 10545cdeb108SMatthew G. Knepley } 1055ccd66a83SVaclav Hapla if (has) *has = exists; 10565cdeb108SMatthew G. Knepley PetscFunctionReturn(0); 10575cdeb108SMatthew G. Knepley } 10585cdeb108SMatthew G. Knepley 1059ecce1506SVaclav Hapla /*@ 10600a9f38efSVaclav Hapla PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file 10610a9f38efSVaclav Hapla 10620a9f38efSVaclav Hapla Input Parameters: 10630a9f38efSVaclav Hapla . viewer - The HDF5 viewer 10640a9f38efSVaclav Hapla 10650a9f38efSVaclav Hapla Output Parameter: 10660a9f38efSVaclav Hapla . has - Flag for group existence 10670a9f38efSVaclav Hapla 10680a9f38efSVaclav Hapla Notes: 10690a9f38efSVaclav Hapla If the path exists but is not a group, this returns PETSC_FALSE as well. 10700a9f38efSVaclav Hapla 10710a9f38efSVaclav Hapla Level: advanced 10720a9f38efSVaclav Hapla 10730a9f38efSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5OpenGroup() 10740a9f38efSVaclav Hapla @*/ 10750a9f38efSVaclav Hapla PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, PetscBool *has) 10760a9f38efSVaclav Hapla { 10770a9f38efSVaclav Hapla H5O_type_t type; 10780a9f38efSVaclav Hapla const char *name; 10790a9f38efSVaclav Hapla PetscErrorCode ierr; 10800a9f38efSVaclav Hapla 10810a9f38efSVaclav Hapla PetscFunctionBegin; 10820a9f38efSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 10830a9f38efSVaclav Hapla PetscValidIntPointer(has,2); 10840a9f38efSVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &name);CHKERRQ(ierr); 10850a9f38efSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, name, PETSC_FALSE, has, &type);CHKERRQ(ierr); 10860a9f38efSVaclav Hapla *has = (type == H5O_TYPE_GROUP) ? PETSC_TRUE : PETSC_FALSE; 10870a9f38efSVaclav Hapla PetscFunctionReturn(0); 10880a9f38efSVaclav Hapla } 10890a9f38efSVaclav Hapla 10900a9f38efSVaclav Hapla /*@ 1091e3f143f7SVaclav Hapla PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group 1092ecce1506SVaclav Hapla 1093ecce1506SVaclav Hapla Input Parameters: 1094ecce1506SVaclav Hapla + viewer - The HDF5 viewer 1095ecce1506SVaclav Hapla - obj - The named object 1096ecce1506SVaclav Hapla 1097ecce1506SVaclav Hapla Output Parameter: 1098ecce1506SVaclav Hapla . has - Flag for dataset existence; PETSC_FALSE for unnamed object 1099ecce1506SVaclav Hapla 1100e3f143f7SVaclav Hapla Notes: 1101e3f143f7SVaclav Hapla If the path exists but is not a dataset, this returns PETSC_FALSE as well. 1102e3f143f7SVaclav Hapla 1103ecce1506SVaclav Hapla Level: advanced 1104ecce1506SVaclav Hapla 1105e3f143f7SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1106ecce1506SVaclav Hapla @*/ 1107ecce1506SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has) 1108ecce1506SVaclav Hapla { 110956cc0592SVaclav Hapla H5O_type_t type; 11106c132bc1SVaclav Hapla char *path; 1111ecce1506SVaclav Hapla PetscErrorCode ierr; 1112ecce1506SVaclav Hapla 1113ecce1506SVaclav Hapla PetscFunctionBegin; 1114c57a1a9aSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1115c57a1a9aSVaclav Hapla PetscValidHeader(obj,2); 1116c57a1a9aSVaclav Hapla PetscValidIntPointer(has,3); 1117ecce1506SVaclav Hapla *has = PETSC_FALSE; 1118e3f143f7SVaclav Hapla if (!obj->name) PetscFunctionReturn(0); 11196c132bc1SVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, obj->name, &path);CHKERRQ(ierr); 11206c132bc1SVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, path, PETSC_FALSE, has, &type);CHKERRQ(ierr); 112156cc0592SVaclav Hapla *has = (type == H5O_TYPE_DATASET) ? PETSC_TRUE : PETSC_FALSE; 11226c132bc1SVaclav Hapla ierr = PetscFree(path);CHKERRQ(ierr); 1123ecce1506SVaclav Hapla PetscFunctionReturn(0); 1124ecce1506SVaclav Hapla } 1125ecce1506SVaclav Hapla 1126df863907SAlex Fikl /*@C 1127ece83b6cSVaclav Hapla PetscViewerHDF5HasAttribute - Check whether an attribute exists 1128e7bdbf8eSMatthew G. Knepley 1129e7bdbf8eSMatthew G. Knepley Input Parameters: 1130e7bdbf8eSMatthew G. Knepley + viewer - The HDF5 viewer 113110e69e8fSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 1132e7bdbf8eSMatthew G. Knepley - name - The attribute name 1133e7bdbf8eSMatthew G. Knepley 1134e7bdbf8eSMatthew G. Knepley Output Parameter: 1135e7bdbf8eSMatthew G. Knepley . has - Flag for attribute existence 1136e7bdbf8eSMatthew G. Knepley 1137e7bdbf8eSMatthew G. Knepley Level: advanced 1138e7bdbf8eSMatthew G. Knepley 1139577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1140e7bdbf8eSMatthew G. Knepley @*/ 114110e69e8fSVaclav Hapla PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscBool *has) 1142e7bdbf8eSMatthew G. Knepley { 114310e69e8fSVaclav Hapla char *parent; 1144e7bdbf8eSMatthew G. Knepley PetscErrorCode ierr; 1145e7bdbf8eSMatthew G. Knepley 1146e7bdbf8eSMatthew G. Knepley PetscFunctionBegin; 11475cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 114810e69e8fSVaclav Hapla if (dataset) PetscValidCharPointer(dataset,2); 1149c57a1a9aSVaclav Hapla PetscValidCharPointer(name,3); 1150c57a1a9aSVaclav Hapla PetscValidIntPointer(has,4); 115110e69e8fSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 1152bb286ee1SVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, has, NULL);CHKERRQ(ierr); 115310e69e8fSVaclav Hapla if (*has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, has);CHKERRQ(ierr);} 115410e69e8fSVaclav Hapla ierr = PetscFree(parent);CHKERRQ(ierr); 115506db490cSVaclav Hapla PetscFunctionReturn(0); 115606db490cSVaclav Hapla } 115706db490cSVaclav Hapla 1158577e0e04SVaclav Hapla /*@C 1159577e0e04SVaclav Hapla PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name 1160577e0e04SVaclav Hapla 1161577e0e04SVaclav Hapla Input Parameters: 1162577e0e04SVaclav Hapla + viewer - The HDF5 viewer 1163577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1164577e0e04SVaclav Hapla - name - The attribute name 1165577e0e04SVaclav Hapla 1166577e0e04SVaclav Hapla Output Parameter: 1167577e0e04SVaclav Hapla . has - Flag for attribute existence 1168577e0e04SVaclav Hapla 1169577e0e04SVaclav Hapla Notes: 1170577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1171577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 1172577e0e04SVaclav Hapla 1173577e0e04SVaclav Hapla Level: advanced 1174577e0e04SVaclav Hapla 1175577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1176577e0e04SVaclav Hapla @*/ 1177577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has) 1178577e0e04SVaclav Hapla { 1179577e0e04SVaclav Hapla PetscErrorCode ierr; 1180577e0e04SVaclav Hapla 1181577e0e04SVaclav Hapla PetscFunctionBegin; 1182577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1183577e0e04SVaclav Hapla PetscValidHeader(obj,2); 1184577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 1185577e0e04SVaclav Hapla PetscValidIntPointer(has,4); 1186577e0e04SVaclav Hapla ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 1187577e0e04SVaclav Hapla ierr = PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);CHKERRQ(ierr); 1188577e0e04SVaclav Hapla PetscFunctionReturn(0); 1189577e0e04SVaclav Hapla } 1190577e0e04SVaclav Hapla 119106db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 119206db490cSVaclav Hapla { 119306db490cSVaclav Hapla hid_t h5; 119406db490cSVaclav Hapla htri_t hhas; 119506db490cSVaclav Hapla PetscErrorCode ierr; 119606db490cSVaclav Hapla 119706db490cSVaclav Hapla PetscFunctionBegin; 119806db490cSVaclav Hapla ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 11992f936e54SVaclav Hapla PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT)); 12005cdeb108SMatthew G. Knepley *has = hhas ? PETSC_TRUE : PETSC_FALSE; 1201e7bdbf8eSMatthew G. Knepley PetscFunctionReturn(0); 1202e7bdbf8eSMatthew G. Knepley } 1203e7bdbf8eSMatthew G. Knepley 1204a75e6a4aSMatthew G. Knepley /* 1205a75e6a4aSMatthew G. Knepley The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that 1206a75e6a4aSMatthew G. Knepley is attached to a communicator, in this case the attribute is a PetscViewer. 1207a75e6a4aSMatthew G. Knepley */ 1208d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID; 1209a75e6a4aSMatthew G. Knepley 1210a75e6a4aSMatthew G. Knepley /*@C 1211a75e6a4aSMatthew G. Knepley PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator. 1212a75e6a4aSMatthew G. Knepley 1213d083f849SBarry Smith Collective 1214a75e6a4aSMatthew G. Knepley 1215a75e6a4aSMatthew G. Knepley Input Parameter: 1216a75e6a4aSMatthew G. Knepley . comm - the MPI communicator to share the HDF5 PetscViewer 1217a75e6a4aSMatthew G. Knepley 1218a75e6a4aSMatthew G. Knepley Level: intermediate 1219a75e6a4aSMatthew G. Knepley 1220a75e6a4aSMatthew G. Knepley Options Database Keys: 1221a75e6a4aSMatthew G. Knepley . -viewer_hdf5_filename <name> 1222a75e6a4aSMatthew G. Knepley 1223a75e6a4aSMatthew G. Knepley Environmental variables: 1224a75e6a4aSMatthew G. Knepley . PETSC_VIEWER_HDF5_FILENAME 1225a75e6a4aSMatthew G. Knepley 1226a75e6a4aSMatthew G. Knepley Notes: 1227a75e6a4aSMatthew G. Knepley Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return 1228a75e6a4aSMatthew G. Knepley an error code. The HDF5 PetscViewer is usually used in the form 1229a75e6a4aSMatthew G. Knepley $ XXXView(XXX object, PETSC_VIEWER_HDF5_(comm)); 1230a75e6a4aSMatthew G. Knepley 1231a75e6a4aSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy() 1232a75e6a4aSMatthew G. Knepley @*/ 1233a75e6a4aSMatthew G. Knepley PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm) 1234a75e6a4aSMatthew G. Knepley { 1235a75e6a4aSMatthew G. Knepley PetscErrorCode ierr; 1236a75e6a4aSMatthew G. Knepley PetscBool flg; 1237a75e6a4aSMatthew G. Knepley PetscViewer viewer; 1238a75e6a4aSMatthew G. Knepley char fname[PETSC_MAX_PATH_LEN]; 1239a75e6a4aSMatthew G. Knepley MPI_Comm ncomm; 1240a75e6a4aSMatthew G. Knepley 1241a75e6a4aSMatthew G. Knepley PetscFunctionBegin; 1242a75e6a4aSMatthew G. Knepley ierr = PetscCommDuplicate(comm,&ncomm,NULL);if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1243a75e6a4aSMatthew G. Knepley if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) { 124412801b39SBarry Smith ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0); 1245a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1246a75e6a4aSMatthew G. Knepley } 124747435625SJed Brown ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg); 1248a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1249a75e6a4aSMatthew G. Knepley if (!flg) { /* PetscViewer not yet created */ 1250a75e6a4aSMatthew G. Knepley ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg); 1251a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1252a75e6a4aSMatthew G. Knepley if (!flg) { 1253a75e6a4aSMatthew G. Knepley ierr = PetscStrcpy(fname,"output.h5"); 1254a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1255a75e6a4aSMatthew G. Knepley } 1256a75e6a4aSMatthew G. Knepley ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer); 1257a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1258a75e6a4aSMatthew G. Knepley ierr = PetscObjectRegisterDestroy((PetscObject)viewer); 1259a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 126047435625SJed Brown ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer); 1261a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1262a75e6a4aSMatthew G. Knepley } 1263a75e6a4aSMatthew G. Knepley ierr = PetscCommDestroy(&ncomm); 1264a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1265a75e6a4aSMatthew G. Knepley PetscFunctionReturn(viewer); 1266a75e6a4aSMatthew G. Knepley } 1267