1*5c6c1daeSBarry Smith #include <petsc-private/viewerimpl.h> /*I "petscsys.h" I*/ 2*5c6c1daeSBarry Smith #include <hdf5.h> 3*5c6c1daeSBarry Smith 4*5c6c1daeSBarry Smith typedef struct GroupList { 5*5c6c1daeSBarry Smith const char *name; 6*5c6c1daeSBarry Smith struct GroupList *next; 7*5c6c1daeSBarry Smith } GroupList; 8*5c6c1daeSBarry Smith 9*5c6c1daeSBarry Smith typedef struct { 10*5c6c1daeSBarry Smith char *filename; 11*5c6c1daeSBarry Smith PetscFileMode btype; 12*5c6c1daeSBarry Smith hid_t file_id; 13*5c6c1daeSBarry Smith PetscInt timestep; 14*5c6c1daeSBarry Smith GroupList *groups; 15*5c6c1daeSBarry Smith } PetscViewer_HDF5; 16*5c6c1daeSBarry Smith 17*5c6c1daeSBarry Smith #undef __FUNCT__ 18*5c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerFileClose_HDF5" 19*5c6c1daeSBarry Smith static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer) 20*5c6c1daeSBarry Smith { 21*5c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data; 22*5c6c1daeSBarry Smith PetscErrorCode ierr; 23*5c6c1daeSBarry Smith 24*5c6c1daeSBarry Smith PetscFunctionBegin; 25*5c6c1daeSBarry Smith ierr = PetscFree(hdf5->filename);CHKERRQ(ierr); 26*5c6c1daeSBarry Smith if (hdf5->file_id) { 27*5c6c1daeSBarry Smith H5Fclose(hdf5->file_id); 28*5c6c1daeSBarry Smith } 29*5c6c1daeSBarry Smith PetscFunctionReturn(0); 30*5c6c1daeSBarry Smith } 31*5c6c1daeSBarry Smith 32*5c6c1daeSBarry Smith #undef __FUNCT__ 33*5c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerDestroy_HDF5" 34*5c6c1daeSBarry Smith PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer) 35*5c6c1daeSBarry Smith { 36*5c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 37*5c6c1daeSBarry Smith PetscErrorCode ierr; 38*5c6c1daeSBarry Smith 39*5c6c1daeSBarry Smith PetscFunctionBegin; 40*5c6c1daeSBarry Smith ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr); 41*5c6c1daeSBarry Smith if (hdf5->groups) { 42*5c6c1daeSBarry Smith while(hdf5->groups) { 43*5c6c1daeSBarry Smith GroupList *tmp = hdf5->groups->next; 44*5c6c1daeSBarry Smith 45*5c6c1daeSBarry Smith ierr = PetscFree(hdf5->groups->name);CHKERRQ(ierr); 46*5c6c1daeSBarry Smith ierr = PetscFree(hdf5->groups);CHKERRQ(ierr); 47*5c6c1daeSBarry Smith hdf5->groups = tmp; 48*5c6c1daeSBarry Smith } 49*5c6c1daeSBarry Smith } 50*5c6c1daeSBarry Smith ierr = PetscFree(hdf5);CHKERRQ(ierr); 51*5c6c1daeSBarry Smith PetscFunctionReturn(0); 52*5c6c1daeSBarry Smith } 53*5c6c1daeSBarry Smith 54*5c6c1daeSBarry Smith EXTERN_C_BEGIN 55*5c6c1daeSBarry Smith #undef __FUNCT__ 56*5c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerFileSetMode_HDF5" 57*5c6c1daeSBarry Smith PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type) 58*5c6c1daeSBarry Smith { 59*5c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 60*5c6c1daeSBarry Smith 61*5c6c1daeSBarry Smith PetscFunctionBegin; 62*5c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 63*5c6c1daeSBarry Smith hdf5->btype = type; 64*5c6c1daeSBarry Smith PetscFunctionReturn(0); 65*5c6c1daeSBarry Smith } 66*5c6c1daeSBarry Smith EXTERN_C_END 67*5c6c1daeSBarry Smith 68*5c6c1daeSBarry Smith EXTERN_C_BEGIN 69*5c6c1daeSBarry Smith #undef __FUNCT__ 70*5c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerFileSetName_HDF5" 71*5c6c1daeSBarry Smith PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[]) 72*5c6c1daeSBarry Smith { 73*5c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 74*5c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 75*5c6c1daeSBarry Smith MPI_Info info = MPI_INFO_NULL; 76*5c6c1daeSBarry Smith #endif 77*5c6c1daeSBarry Smith hid_t plist_id; 78*5c6c1daeSBarry Smith herr_t herr; 79*5c6c1daeSBarry Smith PetscErrorCode ierr; 80*5c6c1daeSBarry Smith 81*5c6c1daeSBarry Smith PetscFunctionBegin; 82*5c6c1daeSBarry Smith ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr); 83*5c6c1daeSBarry Smith /* Set up file access property list with parallel I/O access */ 84*5c6c1daeSBarry Smith plist_id = H5Pcreate(H5P_FILE_ACCESS); 85*5c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 86*5c6c1daeSBarry Smith herr = H5Pset_fapl_mpio(plist_id, ((PetscObject)viewer)->comm, info);CHKERRQ(herr); 87*5c6c1daeSBarry Smith #endif 88*5c6c1daeSBarry Smith /* Create or open the file collectively */ 89*5c6c1daeSBarry Smith switch (hdf5->btype) { 90*5c6c1daeSBarry Smith case FILE_MODE_READ: 91*5c6c1daeSBarry Smith hdf5->file_id = H5Fopen(name, H5F_ACC_RDONLY, plist_id); 92*5c6c1daeSBarry Smith break; 93*5c6c1daeSBarry Smith case FILE_MODE_APPEND: 94*5c6c1daeSBarry Smith hdf5->file_id = H5Fopen(name, H5F_ACC_RDWR, plist_id); 95*5c6c1daeSBarry Smith break; 96*5c6c1daeSBarry Smith case FILE_MODE_WRITE: 97*5c6c1daeSBarry Smith hdf5->file_id = H5Fcreate(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id); 98*5c6c1daeSBarry Smith break; 99*5c6c1daeSBarry Smith default: 100*5c6c1daeSBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 101*5c6c1daeSBarry Smith } 102*5c6c1daeSBarry Smith if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name); 103*5c6c1daeSBarry Smith viewer->format = PETSC_VIEWER_NOFORMAT; 104*5c6c1daeSBarry Smith H5Pclose(plist_id); 105*5c6c1daeSBarry Smith PetscFunctionReturn(0); 106*5c6c1daeSBarry Smith } 107*5c6c1daeSBarry Smith EXTERN_C_END 108*5c6c1daeSBarry Smith 109*5c6c1daeSBarry Smith EXTERN_C_BEGIN 110*5c6c1daeSBarry Smith #undef __FUNCT__ 111*5c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerCreate_HDF5" 112*5c6c1daeSBarry Smith PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v) 113*5c6c1daeSBarry Smith { 114*5c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5; 115*5c6c1daeSBarry Smith PetscErrorCode ierr; 116*5c6c1daeSBarry Smith 117*5c6c1daeSBarry Smith PetscFunctionBegin; 118*5c6c1daeSBarry Smith ierr = PetscNewLog(v, PetscViewer_HDF5, &hdf5);CHKERRQ(ierr); 119*5c6c1daeSBarry Smith 120*5c6c1daeSBarry Smith v->data = (void *) hdf5; 121*5c6c1daeSBarry Smith v->ops->destroy = PetscViewerDestroy_HDF5; 122*5c6c1daeSBarry Smith v->ops->flush = 0; 123*5c6c1daeSBarry Smith v->iformat = 0; 124*5c6c1daeSBarry Smith hdf5->btype = (PetscFileMode) -1; 125*5c6c1daeSBarry Smith hdf5->filename = 0; 126*5c6c1daeSBarry Smith hdf5->timestep = -1; 127*5c6c1daeSBarry Smith hdf5->groups = PETSC_NULL; 128*5c6c1daeSBarry Smith 129*5c6c1daeSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)v,"PetscViewerFileSetName_C","PetscViewerFileSetName_HDF5", 130*5c6c1daeSBarry Smith PetscViewerFileSetName_HDF5);CHKERRQ(ierr); 131*5c6c1daeSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)v,"PetscViewerFileSetMode_C","PetscViewerFileSetMode_HDF5", 132*5c6c1daeSBarry Smith PetscViewerFileSetMode_HDF5);CHKERRQ(ierr); 133*5c6c1daeSBarry Smith PetscFunctionReturn(0); 134*5c6c1daeSBarry Smith } 135*5c6c1daeSBarry Smith EXTERN_C_END 136*5c6c1daeSBarry Smith 137*5c6c1daeSBarry Smith #undef __FUNCT__ 138*5c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5Open" 139*5c6c1daeSBarry Smith /*@C 140*5c6c1daeSBarry Smith PetscViewerHDF5Open - Opens a file for HDF5 input/output. 141*5c6c1daeSBarry Smith 142*5c6c1daeSBarry Smith Collective on MPI_Comm 143*5c6c1daeSBarry Smith 144*5c6c1daeSBarry Smith Input Parameters: 145*5c6c1daeSBarry Smith + comm - MPI communicator 146*5c6c1daeSBarry Smith . name - name of file 147*5c6c1daeSBarry Smith - type - type of file 148*5c6c1daeSBarry Smith $ FILE_MODE_WRITE - create new file for binary output 149*5c6c1daeSBarry Smith $ FILE_MODE_READ - open existing file for binary input 150*5c6c1daeSBarry Smith $ FILE_MODE_APPEND - open existing file for binary output 151*5c6c1daeSBarry Smith 152*5c6c1daeSBarry Smith Output Parameter: 153*5c6c1daeSBarry Smith . hdf5v - PetscViewer for HDF5 input/output to use with the specified file 154*5c6c1daeSBarry Smith 155*5c6c1daeSBarry Smith Level: beginner 156*5c6c1daeSBarry Smith 157*5c6c1daeSBarry Smith Note: 158*5c6c1daeSBarry Smith This PetscViewer should be destroyed with PetscViewerDestroy(). 159*5c6c1daeSBarry Smith 160*5c6c1daeSBarry Smith Concepts: HDF5 files 161*5c6c1daeSBarry Smith Concepts: PetscViewerHDF5^creating 162*5c6c1daeSBarry Smith 163*5c6c1daeSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerSetFormat(), PetscViewerDestroy(), 164*5c6c1daeSBarry Smith VecView(), MatView(), VecLoad(), MatLoad(), 165*5c6c1daeSBarry Smith PetscFileMode, PetscViewer 166*5c6c1daeSBarry Smith @*/ 167*5c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v) 168*5c6c1daeSBarry Smith { 169*5c6c1daeSBarry Smith PetscErrorCode ierr; 170*5c6c1daeSBarry Smith 171*5c6c1daeSBarry Smith PetscFunctionBegin; 172*5c6c1daeSBarry Smith ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr); 173*5c6c1daeSBarry Smith ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr); 174*5c6c1daeSBarry Smith ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr); 175*5c6c1daeSBarry Smith ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr); 176*5c6c1daeSBarry Smith PetscFunctionReturn(0); 177*5c6c1daeSBarry Smith } 178*5c6c1daeSBarry Smith 179*5c6c1daeSBarry Smith #undef __FUNCT__ 180*5c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5GetFileId" 181*5c6c1daeSBarry Smith /*@C 182*5c6c1daeSBarry Smith PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls 183*5c6c1daeSBarry Smith 184*5c6c1daeSBarry Smith Not collective 185*5c6c1daeSBarry Smith 186*5c6c1daeSBarry Smith Input Parameter: 187*5c6c1daeSBarry Smith . viewer - the PetscViewer 188*5c6c1daeSBarry Smith 189*5c6c1daeSBarry Smith Output Parameter: 190*5c6c1daeSBarry Smith . file_id - The file id 191*5c6c1daeSBarry Smith 192*5c6c1daeSBarry Smith Level: intermediate 193*5c6c1daeSBarry Smith 194*5c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open() 195*5c6c1daeSBarry Smith @*/ 196*5c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id) 197*5c6c1daeSBarry Smith { 198*5c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 199*5c6c1daeSBarry Smith 200*5c6c1daeSBarry Smith PetscFunctionBegin; 201*5c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 202*5c6c1daeSBarry Smith if (file_id) *file_id = hdf5->file_id; 203*5c6c1daeSBarry Smith PetscFunctionReturn(0); 204*5c6c1daeSBarry Smith } 205*5c6c1daeSBarry Smith 206*5c6c1daeSBarry Smith #undef __FUNCT__ 207*5c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5PushGroup" 208*5c6c1daeSBarry Smith /*@C 209*5c6c1daeSBarry Smith PetscViewerHDF5PushGroup - Set the current HDF5 group for output 210*5c6c1daeSBarry Smith 211*5c6c1daeSBarry Smith Not collective 212*5c6c1daeSBarry Smith 213*5c6c1daeSBarry Smith Input Parameters: 214*5c6c1daeSBarry Smith + viewer - the PetscViewer 215*5c6c1daeSBarry Smith - name - The group name 216*5c6c1daeSBarry Smith 217*5c6c1daeSBarry Smith Level: intermediate 218*5c6c1daeSBarry Smith 219*5c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 220*5c6c1daeSBarry Smith @*/ 221*5c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name) 222*5c6c1daeSBarry Smith { 223*5c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 224*5c6c1daeSBarry Smith GroupList *groupNode; 225*5c6c1daeSBarry Smith PetscErrorCode ierr; 226*5c6c1daeSBarry Smith 227*5c6c1daeSBarry Smith PetscFunctionBegin; 228*5c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 229*5c6c1daeSBarry Smith PetscValidCharPointer(name,2); 230*5c6c1daeSBarry Smith ierr = PetscMalloc(sizeof(GroupList), &groupNode);CHKERRQ(ierr); 231*5c6c1daeSBarry Smith ierr = PetscStrallocpy(name, (char **) &groupNode->name);CHKERRQ(ierr); 232*5c6c1daeSBarry Smith groupNode->next = hdf5->groups; 233*5c6c1daeSBarry Smith hdf5->groups = groupNode; 234*5c6c1daeSBarry Smith PetscFunctionReturn(0); 235*5c6c1daeSBarry Smith } 236*5c6c1daeSBarry Smith 237*5c6c1daeSBarry Smith #undef __FUNCT__ 238*5c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5PopGroup" 239*5c6c1daeSBarry Smith /*@C 240*5c6c1daeSBarry Smith PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value 241*5c6c1daeSBarry Smith 242*5c6c1daeSBarry Smith Not collective 243*5c6c1daeSBarry Smith 244*5c6c1daeSBarry Smith Input Parameter: 245*5c6c1daeSBarry Smith . viewer - the PetscViewer 246*5c6c1daeSBarry Smith 247*5c6c1daeSBarry Smith Level: intermediate 248*5c6c1daeSBarry Smith 249*5c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup() 250*5c6c1daeSBarry Smith @*/ 251*5c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer) 252*5c6c1daeSBarry Smith { 253*5c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 254*5c6c1daeSBarry Smith GroupList *groupNode; 255*5c6c1daeSBarry Smith PetscErrorCode ierr; 256*5c6c1daeSBarry Smith 257*5c6c1daeSBarry Smith PetscFunctionBegin; 258*5c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 259*5c6c1daeSBarry Smith if (!hdf5->groups) SETERRQ(((PetscObject) viewer)->comm, PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop"); 260*5c6c1daeSBarry Smith groupNode = hdf5->groups; 261*5c6c1daeSBarry Smith hdf5->groups = hdf5->groups->next; 262*5c6c1daeSBarry Smith ierr = PetscFree(groupNode->name);CHKERRQ(ierr); 263*5c6c1daeSBarry Smith ierr = PetscFree(groupNode);CHKERRQ(ierr); 264*5c6c1daeSBarry Smith PetscFunctionReturn(0); 265*5c6c1daeSBarry Smith } 266*5c6c1daeSBarry Smith 267*5c6c1daeSBarry Smith #undef __FUNCT__ 268*5c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5GetGroup" 269*5c6c1daeSBarry Smith /*@C 270*5c6c1daeSBarry Smith PetscViewerHDF5GetGroup - Get the current HDF5 group for output. If none has been assigned, returns PETSC_NULL. 271*5c6c1daeSBarry Smith 272*5c6c1daeSBarry Smith Not collective 273*5c6c1daeSBarry Smith 274*5c6c1daeSBarry Smith Input Parameter: 275*5c6c1daeSBarry Smith . viewer - the PetscViewer 276*5c6c1daeSBarry Smith 277*5c6c1daeSBarry Smith Output Parameter: 278*5c6c1daeSBarry Smith . name - The group name 279*5c6c1daeSBarry Smith 280*5c6c1daeSBarry Smith Level: intermediate 281*5c6c1daeSBarry Smith 282*5c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup() 283*5c6c1daeSBarry Smith @*/ 284*5c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name) 285*5c6c1daeSBarry Smith { 286*5c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 287*5c6c1daeSBarry Smith 288*5c6c1daeSBarry Smith PetscFunctionBegin; 289*5c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 290*5c6c1daeSBarry Smith PetscValidCharPointer(name,2); 291*5c6c1daeSBarry Smith if (hdf5->groups) { 292*5c6c1daeSBarry Smith *name = hdf5->groups->name; 293*5c6c1daeSBarry Smith } else { 294*5c6c1daeSBarry Smith *name = PETSC_NULL; 295*5c6c1daeSBarry Smith } 296*5c6c1daeSBarry Smith PetscFunctionReturn(0); 297*5c6c1daeSBarry Smith } 298*5c6c1daeSBarry Smith 299*5c6c1daeSBarry Smith #undef __FUNCT__ 300*5c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5IncrementTimestep" 301*5c6c1daeSBarry Smith /*@C 302*5c6c1daeSBarry Smith PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time. 303*5c6c1daeSBarry Smith 304*5c6c1daeSBarry Smith Not collective 305*5c6c1daeSBarry Smith 306*5c6c1daeSBarry Smith Input Parameter: 307*5c6c1daeSBarry Smith . viewer - the PetscViewer 308*5c6c1daeSBarry Smith 309*5c6c1daeSBarry Smith Level: intermediate 310*5c6c1daeSBarry Smith 311*5c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep() 312*5c6c1daeSBarry Smith @*/ 313*5c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer) 314*5c6c1daeSBarry Smith { 315*5c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 316*5c6c1daeSBarry Smith 317*5c6c1daeSBarry Smith PetscFunctionBegin; 318*5c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 319*5c6c1daeSBarry Smith ++hdf5->timestep; 320*5c6c1daeSBarry Smith PetscFunctionReturn(0); 321*5c6c1daeSBarry Smith } 322*5c6c1daeSBarry Smith 323*5c6c1daeSBarry Smith #undef __FUNCT__ 324*5c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5SetTimestep" 325*5c6c1daeSBarry Smith /*@C 326*5c6c1daeSBarry Smith PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep 327*5c6c1daeSBarry Smith of -1 disables blocking with timesteps. 328*5c6c1daeSBarry Smith 329*5c6c1daeSBarry Smith Not collective 330*5c6c1daeSBarry Smith 331*5c6c1daeSBarry Smith Input Parameters: 332*5c6c1daeSBarry Smith + viewer - the PetscViewer 333*5c6c1daeSBarry Smith - timestep - The timestep number 334*5c6c1daeSBarry Smith 335*5c6c1daeSBarry Smith Level: intermediate 336*5c6c1daeSBarry Smith 337*5c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 338*5c6c1daeSBarry Smith @*/ 339*5c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep) 340*5c6c1daeSBarry Smith { 341*5c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 342*5c6c1daeSBarry Smith 343*5c6c1daeSBarry Smith PetscFunctionBegin; 344*5c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 345*5c6c1daeSBarry Smith hdf5->timestep = timestep; 346*5c6c1daeSBarry Smith PetscFunctionReturn(0); 347*5c6c1daeSBarry Smith } 348*5c6c1daeSBarry Smith 349*5c6c1daeSBarry Smith #undef __FUNCT__ 350*5c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5GetTimestep" 351*5c6c1daeSBarry Smith /*@C 352*5c6c1daeSBarry Smith PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time. 353*5c6c1daeSBarry Smith 354*5c6c1daeSBarry Smith Not collective 355*5c6c1daeSBarry Smith 356*5c6c1daeSBarry Smith Input Parameter: 357*5c6c1daeSBarry Smith . viewer - the PetscViewer 358*5c6c1daeSBarry Smith 359*5c6c1daeSBarry Smith Output Parameter: 360*5c6c1daeSBarry Smith . timestep - The timestep number 361*5c6c1daeSBarry Smith 362*5c6c1daeSBarry Smith Level: intermediate 363*5c6c1daeSBarry Smith 364*5c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep() 365*5c6c1daeSBarry Smith @*/ 366*5c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep) 367*5c6c1daeSBarry Smith { 368*5c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 369*5c6c1daeSBarry Smith 370*5c6c1daeSBarry Smith PetscFunctionBegin; 371*5c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 372*5c6c1daeSBarry Smith PetscValidPointer(timestep,2); 373*5c6c1daeSBarry Smith *timestep = hdf5->timestep; 374*5c6c1daeSBarry Smith PetscFunctionReturn(0); 375*5c6c1daeSBarry Smith } 376*5c6c1daeSBarry Smith 377*5c6c1daeSBarry Smith #if defined(oldhdf4stuff) 378*5c6c1daeSBarry Smith #undef __FUNCT__ 379*5c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5WriteSDS" 380*5c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5WriteSDS(PetscViewer viewer, float *xf, int d, int *dims,int bs) 381*5c6c1daeSBarry Smith { 382*5c6c1daeSBarry Smith int i; 383*5c6c1daeSBarry Smith PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data; 384*5c6c1daeSBarry Smith int32 sds_id,zero32[3],dims32[3]; 385*5c6c1daeSBarry Smith 386*5c6c1daeSBarry Smith PetscFunctionBegin; 387*5c6c1daeSBarry Smith 388*5c6c1daeSBarry Smith for (i = 0; i < d; i++) { 389*5c6c1daeSBarry Smith zero32[i] = 0; 390*5c6c1daeSBarry Smith dims32[i] = dims[i]; 391*5c6c1daeSBarry Smith } 392*5c6c1daeSBarry Smith sds_id = SDcreate(vhdf5->sd_id, "Vec", DFNT_FLOAT32, d, dims32); 393*5c6c1daeSBarry Smith if (sds_id < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"SDcreate failed"); 394*5c6c1daeSBarry Smith SDwritedata(sds_id, zero32, 0, dims32, xf); 395*5c6c1daeSBarry Smith SDendaccess(sds_id); 396*5c6c1daeSBarry Smith PetscFunctionReturn(0); 397*5c6c1daeSBarry Smith } 398*5c6c1daeSBarry Smith 399*5c6c1daeSBarry Smith #endif 400