15c6c1daeSBarry Smith #include <petsc-private/viewerimpl.h> /*I "petscsys.h" I*/ 2d70abbfaSBarry Smith #include <petscviewerhdf5.h> /*I "petscviewerhdf5.h" I*/ 35c6c1daeSBarry Smith 45c6c1daeSBarry Smith typedef struct GroupList { 55c6c1daeSBarry Smith const char *name; 65c6c1daeSBarry Smith struct GroupList *next; 75c6c1daeSBarry Smith } GroupList; 85c6c1daeSBarry Smith 95c6c1daeSBarry Smith typedef struct { 105c6c1daeSBarry Smith char *filename; 115c6c1daeSBarry Smith PetscFileMode btype; 125c6c1daeSBarry Smith hid_t file_id; 135c6c1daeSBarry Smith PetscInt timestep; 145c6c1daeSBarry Smith GroupList *groups; 15*82ea9b62SBarry Smith PetscBool basedimension2; /* save vectors and DMDA vectors with a dimension of at least 2 even if the bs/dof is 1 */ 165c6c1daeSBarry Smith } PetscViewer_HDF5; 175c6c1daeSBarry Smith 185c6c1daeSBarry Smith #undef __FUNCT__ 19*82ea9b62SBarry Smith #define __FUNCT__ "PetscViewerSetFromOptions_HDF5" 20*82ea9b62SBarry Smith static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptions *PetscOptionsObject,PetscViewer v) 21*82ea9b62SBarry Smith { 22*82ea9b62SBarry Smith PetscErrorCode ierr; 23*82ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data; 24*82ea9b62SBarry Smith 25*82ea9b62SBarry Smith PetscFunctionBegin; 26*82ea9b62SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");CHKERRQ(ierr); 27*82ea9b62SBarry Smith ierr = PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);CHKERRQ(ierr); 28*82ea9b62SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 29*82ea9b62SBarry Smith PetscFunctionReturn(0); 30*82ea9b62SBarry Smith } 31*82ea9b62SBarry Smith 32*82ea9b62SBarry Smith #undef __FUNCT__ 335c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerFileClose_HDF5" 345c6c1daeSBarry Smith static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer) 355c6c1daeSBarry Smith { 365c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data; 375c6c1daeSBarry Smith PetscErrorCode ierr; 385c6c1daeSBarry Smith 395c6c1daeSBarry Smith PetscFunctionBegin; 405c6c1daeSBarry Smith ierr = PetscFree(hdf5->filename);CHKERRQ(ierr); 41729ab6d9SBarry Smith if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 425c6c1daeSBarry Smith PetscFunctionReturn(0); 435c6c1daeSBarry Smith } 445c6c1daeSBarry Smith 455c6c1daeSBarry Smith #undef __FUNCT__ 465c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerDestroy_HDF5" 475c6c1daeSBarry Smith PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer) 485c6c1daeSBarry Smith { 495c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 505c6c1daeSBarry Smith PetscErrorCode ierr; 515c6c1daeSBarry Smith 525c6c1daeSBarry Smith PetscFunctionBegin; 535c6c1daeSBarry Smith ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr); 545c6c1daeSBarry Smith while (hdf5->groups) { 555c6c1daeSBarry Smith GroupList *tmp = hdf5->groups->next; 565c6c1daeSBarry Smith 575c6c1daeSBarry Smith ierr = PetscFree(hdf5->groups->name);CHKERRQ(ierr); 585c6c1daeSBarry Smith ierr = PetscFree(hdf5->groups);CHKERRQ(ierr); 595c6c1daeSBarry Smith hdf5->groups = tmp; 605c6c1daeSBarry Smith } 615c6c1daeSBarry Smith ierr = PetscFree(hdf5);CHKERRQ(ierr); 620b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr); 630b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr); 645c6c1daeSBarry Smith PetscFunctionReturn(0); 655c6c1daeSBarry Smith } 665c6c1daeSBarry Smith 675c6c1daeSBarry Smith #undef __FUNCT__ 685c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerFileSetMode_HDF5" 695c6c1daeSBarry Smith PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type) 705c6c1daeSBarry Smith { 715c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 725c6c1daeSBarry Smith 735c6c1daeSBarry Smith PetscFunctionBegin; 745c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 755c6c1daeSBarry Smith hdf5->btype = type; 765c6c1daeSBarry Smith PetscFunctionReturn(0); 775c6c1daeSBarry Smith } 785c6c1daeSBarry Smith 795c6c1daeSBarry Smith #undef __FUNCT__ 80*82ea9b62SBarry Smith #define __FUNCT__ "PetscViewerHDF5SetBaseDimension2_HDF5" 81*82ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg) 82*82ea9b62SBarry Smith { 83*82ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 84*82ea9b62SBarry Smith 85*82ea9b62SBarry Smith PetscFunctionBegin; 86*82ea9b62SBarry Smith hdf5->basedimension2 = flg; 87*82ea9b62SBarry Smith PetscFunctionReturn(0); 88*82ea9b62SBarry Smith } 89*82ea9b62SBarry Smith 90*82ea9b62SBarry Smith #undef __FUNCT__ 91*82ea9b62SBarry Smith #define __FUNCT__ "PetscViewerHDF5SetBaseDimension2" 92*82ea9b62SBarry Smith /*@C 93*82ea9b62SBarry Smith PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 94*82ea9b62SBarry Smith dimension of 2. 95*82ea9b62SBarry Smith 96*82ea9b62SBarry Smith Logically Collective on PetscViewer 97*82ea9b62SBarry Smith 98*82ea9b62SBarry Smith Input Parameters: 99*82ea9b62SBarry Smith + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 100*82ea9b62SBarry 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 101*82ea9b62SBarry Smith 102*82ea9b62SBarry Smith Options Database: 103*82ea9b62SBarry 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 104*82ea9b62SBarry Smith 105*82ea9b62SBarry Smith 106*82ea9b62SBarry Smith Notes: Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof 107*82ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 108*82ea9b62SBarry Smith 109*82ea9b62SBarry Smith Level: intermediate 110*82ea9b62SBarry Smith 111*82ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 112*82ea9b62SBarry Smith 113*82ea9b62SBarry Smith @*/ 114*82ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg) 115*82ea9b62SBarry Smith { 116*82ea9b62SBarry Smith PetscErrorCode ierr; 117*82ea9b62SBarry Smith 118*82ea9b62SBarry Smith PetscFunctionBegin; 119*82ea9b62SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 120*82ea9b62SBarry Smith ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 121*82ea9b62SBarry Smith PetscFunctionReturn(0); 122*82ea9b62SBarry Smith } 123*82ea9b62SBarry Smith 124*82ea9b62SBarry Smith #undef __FUNCT__ 125*82ea9b62SBarry Smith #define __FUNCT__ "PetscViewerHDF5GetBaseDimension2" 126*82ea9b62SBarry Smith /*@C 127*82ea9b62SBarry Smith PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 128*82ea9b62SBarry Smith dimension of 2. 129*82ea9b62SBarry Smith 130*82ea9b62SBarry Smith Logically Collective on PetscViewer 131*82ea9b62SBarry Smith 132*82ea9b62SBarry Smith Input Parameter: 133*82ea9b62SBarry Smith . viewer - the PetscViewer, must be of type HDF5 134*82ea9b62SBarry Smith 135*82ea9b62SBarry Smith Output Parameter: 136*82ea9b62SBarry 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 137*82ea9b62SBarry Smith 138*82ea9b62SBarry Smith Notes: Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof 139*82ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 140*82ea9b62SBarry Smith 141*82ea9b62SBarry Smith Level: intermediate 142*82ea9b62SBarry Smith 143*82ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 144*82ea9b62SBarry Smith 145*82ea9b62SBarry Smith @*/ 146*82ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg) 147*82ea9b62SBarry Smith { 148*82ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 149*82ea9b62SBarry Smith 150*82ea9b62SBarry Smith PetscFunctionBegin; 151*82ea9b62SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 152*82ea9b62SBarry Smith *flg = hdf5->basedimension2; 153*82ea9b62SBarry Smith PetscFunctionReturn(0); 154*82ea9b62SBarry Smith } 155*82ea9b62SBarry Smith 156*82ea9b62SBarry Smith #undef __FUNCT__ 1575c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerFileSetName_HDF5" 1585c6c1daeSBarry Smith PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[]) 1595c6c1daeSBarry Smith { 1605c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 1615c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 1625c6c1daeSBarry Smith MPI_Info info = MPI_INFO_NULL; 1635c6c1daeSBarry Smith #endif 1645c6c1daeSBarry Smith hid_t plist_id; 1655c6c1daeSBarry Smith PetscErrorCode ierr; 1665c6c1daeSBarry Smith 1675c6c1daeSBarry Smith PetscFunctionBegin; 1685c6c1daeSBarry Smith ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr); 1695c6c1daeSBarry Smith /* Set up file access property list with parallel I/O access */ 170729ab6d9SBarry Smith PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS)); 1715c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 172729ab6d9SBarry Smith PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), info)); 1735c6c1daeSBarry Smith #endif 1745c6c1daeSBarry Smith /* Create or open the file collectively */ 1755c6c1daeSBarry Smith switch (hdf5->btype) { 1765c6c1daeSBarry Smith case FILE_MODE_READ: 177729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id)); 1785c6c1daeSBarry Smith break; 1795c6c1daeSBarry Smith case FILE_MODE_APPEND: 180729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id)); 1815c6c1daeSBarry Smith break; 1825c6c1daeSBarry Smith case FILE_MODE_WRITE: 183729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id)); 1845c6c1daeSBarry Smith break; 1855c6c1daeSBarry Smith default: 1865c6c1daeSBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 1875c6c1daeSBarry Smith } 1885c6c1daeSBarry Smith if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name); 189729ab6d9SBarry Smith PetscStackCallHDF5(H5Pclose,(plist_id)); 1905c6c1daeSBarry Smith PetscFunctionReturn(0); 1915c6c1daeSBarry Smith } 1925c6c1daeSBarry Smith 1935c6c1daeSBarry Smith #undef __FUNCT__ 1945c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerCreate_HDF5" 1958cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v) 1965c6c1daeSBarry Smith { 1975c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5; 1985c6c1daeSBarry Smith PetscErrorCode ierr; 1995c6c1daeSBarry Smith 2005c6c1daeSBarry Smith PetscFunctionBegin; 201b00a9115SJed Brown ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr); 2025c6c1daeSBarry Smith 2035c6c1daeSBarry Smith v->data = (void*) hdf5; 2045c6c1daeSBarry Smith v->ops->destroy = PetscViewerDestroy_HDF5; 205*82ea9b62SBarry Smith v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5; 2065c6c1daeSBarry Smith v->ops->flush = 0; 2075c6c1daeSBarry Smith hdf5->btype = (PetscFileMode) -1; 2085c6c1daeSBarry Smith hdf5->filename = 0; 2095c6c1daeSBarry Smith hdf5->timestep = -1; 2100298fd71SBarry Smith hdf5->groups = NULL; 2115c6c1daeSBarry Smith 2120b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr); 2130b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr); 214*82ea9b62SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr); 2155c6c1daeSBarry Smith PetscFunctionReturn(0); 2165c6c1daeSBarry Smith } 2175c6c1daeSBarry Smith 2185c6c1daeSBarry Smith #undef __FUNCT__ 2195c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5Open" 2205c6c1daeSBarry Smith /*@C 2215c6c1daeSBarry Smith PetscViewerHDF5Open - Opens a file for HDF5 input/output. 2225c6c1daeSBarry Smith 2235c6c1daeSBarry Smith Collective on MPI_Comm 2245c6c1daeSBarry Smith 2255c6c1daeSBarry Smith Input Parameters: 2265c6c1daeSBarry Smith + comm - MPI communicator 2275c6c1daeSBarry Smith . name - name of file 2285c6c1daeSBarry Smith - type - type of file 2295c6c1daeSBarry Smith $ FILE_MODE_WRITE - create new file for binary output 2305c6c1daeSBarry Smith $ FILE_MODE_READ - open existing file for binary input 2315c6c1daeSBarry Smith $ FILE_MODE_APPEND - open existing file for binary output 2325c6c1daeSBarry Smith 2335c6c1daeSBarry Smith Output Parameter: 2345c6c1daeSBarry Smith . hdf5v - PetscViewer for HDF5 input/output to use with the specified file 2355c6c1daeSBarry Smith 236*82ea9b62SBarry Smith Options Database: 237*82ea9b62SBarry 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 238*82ea9b62SBarry Smith 2395c6c1daeSBarry Smith Level: beginner 2405c6c1daeSBarry Smith 2415c6c1daeSBarry Smith Note: 2425c6c1daeSBarry Smith This PetscViewer should be destroyed with PetscViewerDestroy(). 2435c6c1daeSBarry Smith 2445c6c1daeSBarry Smith Concepts: HDF5 files 2455c6c1daeSBarry Smith Concepts: PetscViewerHDF5^creating 2465c6c1daeSBarry Smith 247*82ea9b62SBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerSetFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(), PetscViewerHDF5GetBaseDimension2() 2485c6c1daeSBarry Smith VecView(), MatView(), VecLoad(), MatLoad(), 2495c6c1daeSBarry Smith PetscFileMode, PetscViewer 2505c6c1daeSBarry Smith @*/ 2515c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v) 2525c6c1daeSBarry Smith { 2535c6c1daeSBarry Smith PetscErrorCode ierr; 2545c6c1daeSBarry Smith 2555c6c1daeSBarry Smith PetscFunctionBegin; 2565c6c1daeSBarry Smith ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr); 2575c6c1daeSBarry Smith ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr); 2585c6c1daeSBarry Smith ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr); 2595c6c1daeSBarry Smith ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr); 2605c6c1daeSBarry Smith PetscFunctionReturn(0); 2615c6c1daeSBarry Smith } 2625c6c1daeSBarry Smith 2635c6c1daeSBarry Smith #undef __FUNCT__ 2645c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5GetFileId" 2655c6c1daeSBarry Smith /*@C 2665c6c1daeSBarry Smith PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls 2675c6c1daeSBarry Smith 2685c6c1daeSBarry Smith Not collective 2695c6c1daeSBarry Smith 2705c6c1daeSBarry Smith Input Parameter: 2715c6c1daeSBarry Smith . viewer - the PetscViewer 2725c6c1daeSBarry Smith 2735c6c1daeSBarry Smith Output Parameter: 2745c6c1daeSBarry Smith . file_id - The file id 2755c6c1daeSBarry Smith 2765c6c1daeSBarry Smith Level: intermediate 2775c6c1daeSBarry Smith 2785c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open() 2795c6c1daeSBarry Smith @*/ 2805c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id) 2815c6c1daeSBarry Smith { 2825c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 2835c6c1daeSBarry Smith 2845c6c1daeSBarry Smith PetscFunctionBegin; 2855c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 2865c6c1daeSBarry Smith if (file_id) *file_id = hdf5->file_id; 2875c6c1daeSBarry Smith PetscFunctionReturn(0); 2885c6c1daeSBarry Smith } 2895c6c1daeSBarry Smith 2905c6c1daeSBarry Smith #undef __FUNCT__ 2915c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5PushGroup" 2925c6c1daeSBarry Smith /*@C 2935c6c1daeSBarry Smith PetscViewerHDF5PushGroup - Set the current HDF5 group for output 2945c6c1daeSBarry Smith 2955c6c1daeSBarry Smith Not collective 2965c6c1daeSBarry Smith 2975c6c1daeSBarry Smith Input Parameters: 2985c6c1daeSBarry Smith + viewer - the PetscViewer 2995c6c1daeSBarry Smith - name - The group name 3005c6c1daeSBarry Smith 3015c6c1daeSBarry Smith Level: intermediate 3025c6c1daeSBarry Smith 3035c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 3045c6c1daeSBarry Smith @*/ 3055c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name) 3065c6c1daeSBarry Smith { 3075c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 3085c6c1daeSBarry Smith GroupList *groupNode; 3095c6c1daeSBarry Smith PetscErrorCode ierr; 3105c6c1daeSBarry Smith 3115c6c1daeSBarry Smith PetscFunctionBegin; 3125c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 3135c6c1daeSBarry Smith PetscValidCharPointer(name,2); 3145c6c1daeSBarry Smith ierr = PetscMalloc(sizeof(GroupList), &groupNode);CHKERRQ(ierr); 3155c6c1daeSBarry Smith ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr); 316a297a907SKarl Rupp 3175c6c1daeSBarry Smith groupNode->next = hdf5->groups; 3185c6c1daeSBarry Smith hdf5->groups = groupNode; 3195c6c1daeSBarry Smith PetscFunctionReturn(0); 3205c6c1daeSBarry Smith } 3215c6c1daeSBarry Smith 3225c6c1daeSBarry Smith #undef __FUNCT__ 3235c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5PopGroup" 3243ef9c667SSatish Balay /*@ 3255c6c1daeSBarry Smith PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value 3265c6c1daeSBarry Smith 3275c6c1daeSBarry Smith Not collective 3285c6c1daeSBarry Smith 3295c6c1daeSBarry Smith Input Parameter: 3305c6c1daeSBarry Smith . viewer - the PetscViewer 3315c6c1daeSBarry Smith 3325c6c1daeSBarry Smith Level: intermediate 3335c6c1daeSBarry Smith 3345c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup() 3355c6c1daeSBarry Smith @*/ 3365c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer) 3375c6c1daeSBarry Smith { 3385c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 3395c6c1daeSBarry Smith GroupList *groupNode; 3405c6c1daeSBarry Smith PetscErrorCode ierr; 3415c6c1daeSBarry Smith 3425c6c1daeSBarry Smith PetscFunctionBegin; 3435c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 34482f516ccSBarry Smith if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop"); 3455c6c1daeSBarry Smith groupNode = hdf5->groups; 3465c6c1daeSBarry Smith hdf5->groups = hdf5->groups->next; 3475c6c1daeSBarry Smith ierr = PetscFree(groupNode->name);CHKERRQ(ierr); 3485c6c1daeSBarry Smith ierr = PetscFree(groupNode);CHKERRQ(ierr); 3495c6c1daeSBarry Smith PetscFunctionReturn(0); 3505c6c1daeSBarry Smith } 3515c6c1daeSBarry Smith 3525c6c1daeSBarry Smith #undef __FUNCT__ 3535c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5GetGroup" 3545c6c1daeSBarry Smith /*@C 3550298fd71SBarry Smith PetscViewerHDF5GetGroup - Get the current HDF5 group for output. If none has been assigned, returns NULL. 3565c6c1daeSBarry Smith 3575c6c1daeSBarry Smith Not collective 3585c6c1daeSBarry Smith 3595c6c1daeSBarry Smith Input Parameter: 3605c6c1daeSBarry Smith . viewer - the PetscViewer 3615c6c1daeSBarry Smith 3625c6c1daeSBarry Smith Output Parameter: 3635c6c1daeSBarry Smith . name - The group name 3645c6c1daeSBarry Smith 3655c6c1daeSBarry Smith Level: intermediate 3665c6c1daeSBarry Smith 3675c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup() 3685c6c1daeSBarry Smith @*/ 3695c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name) 3705c6c1daeSBarry Smith { 3715c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 3725c6c1daeSBarry Smith 3735c6c1daeSBarry Smith PetscFunctionBegin; 3745c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 375c959eef4SJed Brown PetscValidPointer(name,2); 376a297a907SKarl Rupp if (hdf5->groups) *name = hdf5->groups->name; 3770298fd71SBarry Smith else *name = NULL; 3785c6c1daeSBarry Smith PetscFunctionReturn(0); 3795c6c1daeSBarry Smith } 3805c6c1daeSBarry Smith 3815c6c1daeSBarry Smith #undef __FUNCT__ 3825c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5IncrementTimestep" 3833ef9c667SSatish Balay /*@ 3845c6c1daeSBarry Smith PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time. 3855c6c1daeSBarry Smith 3865c6c1daeSBarry Smith Not collective 3875c6c1daeSBarry Smith 3885c6c1daeSBarry Smith Input Parameter: 3895c6c1daeSBarry Smith . viewer - the PetscViewer 3905c6c1daeSBarry Smith 3915c6c1daeSBarry Smith Level: intermediate 3925c6c1daeSBarry Smith 3935c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep() 3945c6c1daeSBarry Smith @*/ 3955c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer) 3965c6c1daeSBarry Smith { 3975c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 3985c6c1daeSBarry Smith 3995c6c1daeSBarry Smith PetscFunctionBegin; 4005c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 4015c6c1daeSBarry Smith ++hdf5->timestep; 4025c6c1daeSBarry Smith PetscFunctionReturn(0); 4035c6c1daeSBarry Smith } 4045c6c1daeSBarry Smith 4055c6c1daeSBarry Smith #undef __FUNCT__ 4065c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5SetTimestep" 4073ef9c667SSatish Balay /*@ 4085c6c1daeSBarry Smith PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep 4095c6c1daeSBarry Smith of -1 disables blocking with timesteps. 4105c6c1daeSBarry Smith 4115c6c1daeSBarry Smith Not collective 4125c6c1daeSBarry Smith 4135c6c1daeSBarry Smith Input Parameters: 4145c6c1daeSBarry Smith + viewer - the PetscViewer 4155c6c1daeSBarry Smith - timestep - The timestep number 4165c6c1daeSBarry Smith 4175c6c1daeSBarry Smith Level: intermediate 4185c6c1daeSBarry Smith 4195c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 4205c6c1daeSBarry Smith @*/ 4215c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep) 4225c6c1daeSBarry Smith { 4235c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 4245c6c1daeSBarry Smith 4255c6c1daeSBarry Smith PetscFunctionBegin; 4265c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 4275c6c1daeSBarry Smith hdf5->timestep = timestep; 4285c6c1daeSBarry Smith PetscFunctionReturn(0); 4295c6c1daeSBarry Smith } 4305c6c1daeSBarry Smith 4315c6c1daeSBarry Smith #undef __FUNCT__ 4325c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5GetTimestep" 4333ef9c667SSatish Balay /*@ 4345c6c1daeSBarry Smith PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time. 4355c6c1daeSBarry Smith 4365c6c1daeSBarry Smith Not collective 4375c6c1daeSBarry Smith 4385c6c1daeSBarry Smith Input Parameter: 4395c6c1daeSBarry Smith . viewer - the PetscViewer 4405c6c1daeSBarry Smith 4415c6c1daeSBarry Smith Output Parameter: 4425c6c1daeSBarry Smith . timestep - The timestep number 4435c6c1daeSBarry Smith 4445c6c1daeSBarry Smith Level: intermediate 4455c6c1daeSBarry Smith 4465c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep() 4475c6c1daeSBarry Smith @*/ 4485c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep) 4495c6c1daeSBarry Smith { 4505c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 4515c6c1daeSBarry Smith 4525c6c1daeSBarry Smith PetscFunctionBegin; 4535c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 4545c6c1daeSBarry Smith PetscValidPointer(timestep,2); 4555c6c1daeSBarry Smith *timestep = hdf5->timestep; 4565c6c1daeSBarry Smith PetscFunctionReturn(0); 4575c6c1daeSBarry Smith } 4585c6c1daeSBarry Smith 45936bce6e8SMatthew G. Knepley #undef __FUNCT__ 46036bce6e8SMatthew G. Knepley #define __FUNCT__ "PetscDataTypeToHDF5DataType" 46136bce6e8SMatthew G. Knepley /*@C 46236bce6e8SMatthew G. Knepley PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name. 46336bce6e8SMatthew G. Knepley 46436bce6e8SMatthew G. Knepley Not collective 46536bce6e8SMatthew G. Knepley 46636bce6e8SMatthew G. Knepley Input Parameter: 46736bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 46836bce6e8SMatthew G. Knepley 46936bce6e8SMatthew G. Knepley Output Parameter: 47036bce6e8SMatthew G. Knepley . mtype - the MPI datatype (for example MPI_DOUBLE, ...) 47136bce6e8SMatthew G. Knepley 47236bce6e8SMatthew G. Knepley Level: advanced 47336bce6e8SMatthew G. Knepley 47436bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 47536bce6e8SMatthew G. Knepley @*/ 47636bce6e8SMatthew G. Knepley PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype) 47736bce6e8SMatthew G. Knepley { 47836bce6e8SMatthew G. Knepley PetscFunctionBegin; 47936bce6e8SMatthew G. Knepley if (ptype == PETSC_INT) 48036bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 48136bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_LLONG; 48236bce6e8SMatthew G. Knepley #else 48336bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_INT; 48436bce6e8SMatthew G. Knepley #endif 48536bce6e8SMatthew G. Knepley else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE; 48636bce6e8SMatthew G. Knepley else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG; 48736bce6e8SMatthew G. Knepley else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT; 48836bce6e8SMatthew G. Knepley else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_DOUBLE; 48936bce6e8SMatthew G. Knepley else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_DOUBLE; 49036bce6e8SMatthew G. Knepley else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT; 49136bce6e8SMatthew G. Knepley else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR; 49236bce6e8SMatthew G. Knepley else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR; 4937e97c476SMatthew G. Knepley else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1); 49436bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype"); 49536bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 49636bce6e8SMatthew G. Knepley } 49736bce6e8SMatthew G. Knepley 49836bce6e8SMatthew G. Knepley #undef __FUNCT__ 49936bce6e8SMatthew G. Knepley #define __FUNCT__ "PetscHDF5DataTypeToPetscDataType" 50036bce6e8SMatthew G. Knepley /*@C 50136bce6e8SMatthew G. Knepley PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name 50236bce6e8SMatthew G. Knepley 50336bce6e8SMatthew G. Knepley Not collective 50436bce6e8SMatthew G. Knepley 50536bce6e8SMatthew G. Knepley Input Parameter: 50636bce6e8SMatthew G. Knepley . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...) 50736bce6e8SMatthew G. Knepley 50836bce6e8SMatthew G. Knepley Output Parameter: 50936bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 51036bce6e8SMatthew G. Knepley 51136bce6e8SMatthew G. Knepley Level: advanced 51236bce6e8SMatthew G. Knepley 51336bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 51436bce6e8SMatthew G. Knepley @*/ 51536bce6e8SMatthew G. Knepley PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype) 51636bce6e8SMatthew G. Knepley { 51736bce6e8SMatthew G. Knepley PetscFunctionBegin; 51836bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 51936bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG; 52036bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT; 52136bce6e8SMatthew G. Knepley #else 52236bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT; 52336bce6e8SMatthew G. Knepley #endif 52436bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE; 52536bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG; 52636bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT; 52736bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT; 52836bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR; 52936bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR; 5307e97c476SMatthew G. Knepley else if (htype == H5T_C_S1) *ptype = PETSC_STRING; 53136bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype"); 53236bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 53336bce6e8SMatthew G. Knepley } 53436bce6e8SMatthew G. Knepley 53536bce6e8SMatthew G. Knepley #undef __FUNCT__ 53636bce6e8SMatthew G. Knepley #define __FUNCT__ "PetscViewerHDF5WriteAttribute" 53736bce6e8SMatthew G. Knepley /*@ 53836bce6e8SMatthew G. Knepley PetscViewerHDF5WriteAttribute - Write a scalar attribute 53936bce6e8SMatthew G. Knepley 54036bce6e8SMatthew G. Knepley Input Parameters: 54136bce6e8SMatthew G. Knepley + viewer - The HDF5 viewer 54236bce6e8SMatthew G. Knepley . parent - The parent name 54336bce6e8SMatthew G. Knepley . name - The attribute name 54436bce6e8SMatthew G. Knepley . datatype - The attribute type 54536bce6e8SMatthew G. Knepley - value - The attribute value 54636bce6e8SMatthew G. Knepley 54736bce6e8SMatthew G. Knepley Level: advanced 54836bce6e8SMatthew G. Knepley 549e7bdbf8eSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute() 55036bce6e8SMatthew G. Knepley @*/ 55136bce6e8SMatthew G. Knepley PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value) 55236bce6e8SMatthew G. Knepley { 553729ab6d9SBarry Smith hid_t h5, dataspace, dataset, attribute, dtype; 55436bce6e8SMatthew G. Knepley PetscErrorCode ierr; 55536bce6e8SMatthew G. Knepley 55636bce6e8SMatthew G. Knepley PetscFunctionBegin; 5575cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 55836bce6e8SMatthew G. Knepley PetscValidPointer(parent, 2); 55936bce6e8SMatthew G. Knepley PetscValidPointer(name, 3); 56036bce6e8SMatthew G. Knepley PetscValidPointer(value, 4); 56136bce6e8SMatthew G. Knepley ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 5627e97c476SMatthew G. Knepley if (datatype == PETSC_STRING) { 5637e97c476SMatthew G. Knepley size_t len; 5647e97c476SMatthew G. Knepley ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr); 565729ab6d9SBarry Smith PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 5667e97c476SMatthew G. Knepley } 56736bce6e8SMatthew G. Knepley ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 568729ab6d9SBarry Smith PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR)); 56936bce6e8SMatthew G. Knepley #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800) 570729ab6d9SBarry Smith PetscStackCallHDF5Return(dataset,H5Dopen2,(h5, parent, H5P_DEFAULT)); 571729ab6d9SBarry Smith PetscStackCallHDF5Return(attribute,H5Acreate2,(dataset, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT)); 57236bce6e8SMatthew G. Knepley #else 573729ab6d9SBarry Smith PetscStackCallHDF5Return(dataset,H5Dopen,(h5, parent)); 574729ab6d9SBarry Smith PetscStackCallHDF5Return(attribute,H5Acreate,(dataset, name, dtype, dataspace, H5P_DEFAULT)); 57536bce6e8SMatthew G. Knepley #endif 576729ab6d9SBarry Smith PetscStackCallHDF5(H5Awrite,(attribute, dtype, value)); 577729ab6d9SBarry Smith if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype)); 578729ab6d9SBarry Smith PetscStackCallHDF5(H5Aclose,(attribute)); 579729ab6d9SBarry Smith PetscStackCallHDF5(H5Dclose,(dataset)); 580729ab6d9SBarry Smith PetscStackCallHDF5(H5Sclose,(dataspace)); 58136bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 58236bce6e8SMatthew G. Knepley } 58336bce6e8SMatthew G. Knepley 584ad0c61feSMatthew G. Knepley #undef __FUNCT__ 585ad0c61feSMatthew G. Knepley #define __FUNCT__ "PetscViewerHDF5ReadAttribute" 586ad0c61feSMatthew G. Knepley /*@ 587ad0c61feSMatthew G. Knepley PetscViewerHDF5ReadAttribute - Read a scalar attribute 588ad0c61feSMatthew G. Knepley 589ad0c61feSMatthew G. Knepley Input Parameters: 590ad0c61feSMatthew G. Knepley + viewer - The HDF5 viewer 591ad0c61feSMatthew G. Knepley . parent - The parent name 592ad0c61feSMatthew G. Knepley . name - The attribute name 593ad0c61feSMatthew G. Knepley - datatype - The attribute type 594ad0c61feSMatthew G. Knepley 595ad0c61feSMatthew G. Knepley Output Parameter: 596ad0c61feSMatthew G. Knepley . value - The attribute value 597ad0c61feSMatthew G. Knepley 598ad0c61feSMatthew G. Knepley Level: advanced 599ad0c61feSMatthew G. Knepley 600e7bdbf8eSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute() 601ad0c61feSMatthew G. Knepley @*/ 602ad0c61feSMatthew G. Knepley PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, void *value) 603ad0c61feSMatthew G. Knepley { 604729ab6d9SBarry Smith hid_t h5, dataspace, dataset, attribute, dtype; 605ad0c61feSMatthew G. Knepley PetscErrorCode ierr; 606ad0c61feSMatthew G. Knepley 607ad0c61feSMatthew G. Knepley PetscFunctionBegin; 6085cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 609ad0c61feSMatthew G. Knepley PetscValidPointer(parent, 2); 610ad0c61feSMatthew G. Knepley PetscValidPointer(name, 3); 611ad0c61feSMatthew G. Knepley PetscValidPointer(value, 4); 612ad0c61feSMatthew G. Knepley ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 613ad0c61feSMatthew G. Knepley ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 614729ab6d9SBarry Smith PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR)); 615ad0c61feSMatthew G. Knepley #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800) 616729ab6d9SBarry Smith PetscStackCallHDF5Return(dataset,H5Dopen2,(h5, parent, H5P_DEFAULT)); 617ad0c61feSMatthew G. Knepley #else 618729ab6d9SBarry Smith PetscStackCallHDF5Return(dataset,H5Dopen,(h5, parent)); 619ad0c61feSMatthew G. Knepley #endif 620729ab6d9SBarry Smith PetscStackCallHDF5Return(attribute,H5Aopen_name,(dataset, name)); 62170efba33SBarry Smith PetscStackCallHDF5(H5Aread,(attribute, dtype, value)); 622729ab6d9SBarry Smith PetscStackCallHDF5(H5Aclose,(attribute)); 623729ab6d9SBarry Smith PetscStackCallHDF5(H5Dclose,(dataset)); 624729ab6d9SBarry Smith PetscStackCallHDF5(H5Sclose,(dataspace)); 625ad0c61feSMatthew G. Knepley PetscFunctionReturn(0); 626ad0c61feSMatthew G. Knepley } 627ad0c61feSMatthew G. Knepley 628e7bdbf8eSMatthew G. Knepley #undef __FUNCT__ 6295cdeb108SMatthew G. Knepley #define __FUNCT__ "PetscViewerHDF5HasObject" 6305cdeb108SMatthew G. Knepley static PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, const char name[], H5O_type_t otype, PetscBool *has) 6315cdeb108SMatthew G. Knepley { 6325cdeb108SMatthew G. Knepley hid_t h5; 6335cdeb108SMatthew G. Knepley PetscErrorCode ierr; 6345cdeb108SMatthew G. Knepley 6355cdeb108SMatthew G. Knepley PetscFunctionBegin; 6365cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 6375cdeb108SMatthew G. Knepley PetscValidPointer(name, 2); 6385cdeb108SMatthew G. Knepley PetscValidPointer(has, 3); 6395cdeb108SMatthew G. Knepley *has = PETSC_FALSE; 640c0cd0301SJed Brown ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 6415cdeb108SMatthew G. Knepley if (H5Lexists(h5, name, H5P_DEFAULT)) { 6425cdeb108SMatthew G. Knepley H5O_info_t info; 6435cdeb108SMatthew G. Knepley hid_t obj; 6445cdeb108SMatthew G. Knepley 645729ab6d9SBarry Smith PetscStackCallHDF5Return(obj,H5Oopen,(h5, name, H5P_DEFAULT)); 646729ab6d9SBarry Smith PetscStackCallHDF5(H5Oget_info,(obj, &info)); 6475cdeb108SMatthew G. Knepley if (otype == info.type) *has = PETSC_TRUE; 648729ab6d9SBarry Smith PetscStackCallHDF5(H5Oclose,(obj)); 6495cdeb108SMatthew G. Knepley } 6505cdeb108SMatthew G. Knepley PetscFunctionReturn(0); 6515cdeb108SMatthew G. Knepley } 6525cdeb108SMatthew G. Knepley 6535cdeb108SMatthew G. Knepley #undef __FUNCT__ 654e7bdbf8eSMatthew G. Knepley #define __FUNCT__ "PetscViewerHDF5HasAttribute" 655e7bdbf8eSMatthew G. Knepley /*@ 656e7bdbf8eSMatthew G. Knepley PetscViewerHDF5HasAttribute - Check whether a scalar attribute exists 657e7bdbf8eSMatthew G. Knepley 658e7bdbf8eSMatthew G. Knepley Input Parameters: 659e7bdbf8eSMatthew G. Knepley + viewer - The HDF5 viewer 660e7bdbf8eSMatthew G. Knepley . parent - The parent name 661e7bdbf8eSMatthew G. Knepley - name - The attribute name 662e7bdbf8eSMatthew G. Knepley 663e7bdbf8eSMatthew G. Knepley Output Parameter: 664e7bdbf8eSMatthew G. Knepley . has - Flag for attribute existence 665e7bdbf8eSMatthew G. Knepley 666e7bdbf8eSMatthew G. Knepley Level: advanced 667e7bdbf8eSMatthew G. Knepley 668e7bdbf8eSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute() 669e7bdbf8eSMatthew G. Knepley @*/ 670e7bdbf8eSMatthew G. Knepley PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 671e7bdbf8eSMatthew G. Knepley { 672729ab6d9SBarry Smith hid_t h5, dataset; 6735cdeb108SMatthew G. Knepley htri_t hhas; 6745cdeb108SMatthew G. Knepley PetscBool exists; 675e7bdbf8eSMatthew G. Knepley PetscErrorCode ierr; 676e7bdbf8eSMatthew G. Knepley 677e7bdbf8eSMatthew G. Knepley PetscFunctionBegin; 6785cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 679e7bdbf8eSMatthew G. Knepley PetscValidPointer(parent, 2); 680e7bdbf8eSMatthew G. Knepley PetscValidPointer(name, 3); 681e7bdbf8eSMatthew G. Knepley PetscValidPointer(has, 4); 682e7bdbf8eSMatthew G. Knepley *has = PETSC_FALSE; 683e7bdbf8eSMatthew G. Knepley ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 6845cdeb108SMatthew G. Knepley ierr = PetscViewerHDF5HasObject(viewer, parent, H5O_TYPE_DATASET, &exists);CHKERRQ(ierr); 6855cdeb108SMatthew G. Knepley if (exists) { 686e7bdbf8eSMatthew G. Knepley #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800) 687729ab6d9SBarry Smith PetscStackCall("H5Dopen2",dataset = H5Dopen2(h5, parent, H5P_DEFAULT)); 688e7bdbf8eSMatthew G. Knepley #else 689729ab6d9SBarry Smith PetscStackCall("H5Dopen",dataset = H5Dopen(h5, parent)); 690e7bdbf8eSMatthew G. Knepley #endif 691729ab6d9SBarry Smith if (dataset < 0) PetscFunctionReturn(0); 692729ab6d9SBarry Smith PetscStackCall("H5Aexists",hhas = H5Aexists(dataset, name)); 693729ab6d9SBarry Smith if (hhas < 0) { 694729ab6d9SBarry Smith PetscStackCallHDF5(H5Dclose,(dataset)); 695729ab6d9SBarry Smith PetscFunctionReturn(0); 696729ab6d9SBarry Smith } 697729ab6d9SBarry Smith PetscStackCallHDF5(H5Dclose,(dataset)); 6985cdeb108SMatthew G. Knepley *has = hhas ? PETSC_TRUE : PETSC_FALSE; 6995cdeb108SMatthew G. Knepley } 700e7bdbf8eSMatthew G. Knepley PetscFunctionReturn(0); 701e7bdbf8eSMatthew G. Knepley } 702e7bdbf8eSMatthew G. Knepley 703a75e6a4aSMatthew G. Knepley /* 704a75e6a4aSMatthew G. Knepley The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that 705a75e6a4aSMatthew G. Knepley is attached to a communicator, in this case the attribute is a PetscViewer. 706a75e6a4aSMatthew G. Knepley */ 707a75e6a4aSMatthew G. Knepley static int Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID; 708a75e6a4aSMatthew G. Knepley 709a75e6a4aSMatthew G. Knepley #undef __FUNCT__ 710a75e6a4aSMatthew G. Knepley #define __FUNCT__ "PETSC_VIEWER_HDF5_" 711a75e6a4aSMatthew G. Knepley /*@C 712a75e6a4aSMatthew G. Knepley PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator. 713a75e6a4aSMatthew G. Knepley 714a75e6a4aSMatthew G. Knepley Collective on MPI_Comm 715a75e6a4aSMatthew G. Knepley 716a75e6a4aSMatthew G. Knepley Input Parameter: 717a75e6a4aSMatthew G. Knepley . comm - the MPI communicator to share the HDF5 PetscViewer 718a75e6a4aSMatthew G. Knepley 719a75e6a4aSMatthew G. Knepley Level: intermediate 720a75e6a4aSMatthew G. Knepley 721a75e6a4aSMatthew G. Knepley Options Database Keys: 722a75e6a4aSMatthew G. Knepley . -viewer_hdf5_filename <name> 723a75e6a4aSMatthew G. Knepley 724a75e6a4aSMatthew G. Knepley Environmental variables: 725a75e6a4aSMatthew G. Knepley . PETSC_VIEWER_HDF5_FILENAME 726a75e6a4aSMatthew G. Knepley 727a75e6a4aSMatthew G. Knepley Notes: 728a75e6a4aSMatthew G. Knepley Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return 729a75e6a4aSMatthew G. Knepley an error code. The HDF5 PetscViewer is usually used in the form 730a75e6a4aSMatthew G. Knepley $ XXXView(XXX object, PETSC_VIEWER_HDF5_(comm)); 731a75e6a4aSMatthew G. Knepley 732a75e6a4aSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy() 733a75e6a4aSMatthew G. Knepley @*/ 734a75e6a4aSMatthew G. Knepley PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm) 735a75e6a4aSMatthew G. Knepley { 736a75e6a4aSMatthew G. Knepley PetscErrorCode ierr; 737a75e6a4aSMatthew G. Knepley PetscBool flg; 738a75e6a4aSMatthew G. Knepley PetscViewer viewer; 739a75e6a4aSMatthew G. Knepley char fname[PETSC_MAX_PATH_LEN]; 740a75e6a4aSMatthew G. Knepley MPI_Comm ncomm; 741a75e6a4aSMatthew G. Knepley 742a75e6a4aSMatthew G. Knepley PetscFunctionBegin; 743a75e6a4aSMatthew 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);} 744a75e6a4aSMatthew G. Knepley if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) { 745a75e6a4aSMatthew G. Knepley ierr = MPI_Keyval_create(MPI_NULL_COPY_FN,MPI_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0); 746a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 747a75e6a4aSMatthew G. Knepley } 748a75e6a4aSMatthew G. Knepley ierr = MPI_Attr_get(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg); 749a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 750a75e6a4aSMatthew G. Knepley if (!flg) { /* PetscViewer not yet created */ 751a75e6a4aSMatthew G. Knepley ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg); 752a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 753a75e6a4aSMatthew G. Knepley if (!flg) { 754a75e6a4aSMatthew G. Knepley ierr = PetscStrcpy(fname,"output.h5"); 755a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 756a75e6a4aSMatthew G. Knepley } 757a75e6a4aSMatthew G. Knepley ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer); 758a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 759a75e6a4aSMatthew G. Knepley ierr = PetscObjectRegisterDestroy((PetscObject)viewer); 760a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 761a75e6a4aSMatthew G. Knepley ierr = MPI_Attr_put(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer); 762a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 763a75e6a4aSMatthew G. Knepley } 764a75e6a4aSMatthew G. Knepley ierr = PetscCommDestroy(&ncomm); 765a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 766a75e6a4aSMatthew G. Knepley PetscFunctionReturn(viewer); 767a75e6a4aSMatthew G. Knepley } 768