Actual source code: ex19.c

  1: static char help[] = "Parallel HDF5 Vec Viewing.\n\n";

  3: #include <petscvec.h>
  4: #include <petscviewerhdf5.h>

  6: int main(int argc,char **argv)
  7: {
  8:   Vec            x1, x2, *x3ts, *x4ts;
  9:   Vec            x1r, x2r, x3r, x4r;
 10:   PetscViewer    viewer;
 11:   PetscRandom    rand;
 12:   PetscMPIInt    rank;
 13:   PetscInt       i, n = 6, n_timesteps = 5;
 14:   PetscBool      equal;
 15:   MPI_Comm       comm;

 17:   PetscInitialize(&argc, &argv, (char*) 0, help);
 18:   comm = PETSC_COMM_WORLD;
 19:   MPI_Comm_rank(comm, &rank);
 20:   PetscOptionsGetInt(NULL,NULL, "-n", &n, NULL);
 21:   PetscOptionsGetInt(NULL,NULL, "-n_timesteps", &n_timesteps, NULL);

 24:   /* create, initialize and write vectors */
 25:   PetscRandomCreate(comm, &rand);
 26:   PetscRandomSetFromOptions(rand);
 27:   PetscViewerHDF5Open(comm, "ex19.h5", FILE_MODE_WRITE, &viewer);

 29:   VecCreate(comm, &x1);
 30:   PetscObjectSetName((PetscObject) x1, "x1");
 31:   VecSetSizes(x1, PETSC_DECIDE, n);
 32:   VecSetFromOptions(x1);
 33:   VecSetRandom(x1, rand);
 34:   VecView(x1, viewer);

 36:   PetscViewerHDF5PushGroup(viewer, "/testBlockSize");
 37:   VecCreate(comm, &x2);
 38:   PetscObjectSetName((PetscObject) x2, "x2");
 39:   VecSetSizes(x2, PETSC_DECIDE, n);
 40:   VecSetBlockSize(x2, 2);
 41:   VecSetFromOptions(x2);
 42:   VecSetRandom(x2, rand);
 43:   VecView(x2, viewer);
 44:   PetscViewerHDF5PopGroup(viewer);

 46:   PetscViewerHDF5PushGroup(viewer, "/testTimestep");
 47:   PetscViewerHDF5PushTimestepping(viewer);

 49:   VecDuplicateVecs(x1, n_timesteps, &x3ts);
 50:   for (i=0; i<n_timesteps; i++) {
 51:     PetscObjectSetName((PetscObject) x3ts[i], "x3");
 52:     VecSetRandom(x3ts[i], rand);
 53:     VecView(x3ts[i], viewer);
 54:     PetscViewerHDF5IncrementTimestep(viewer);
 55:   }

 57:   PetscViewerHDF5PushGroup(viewer, "testBlockSize");
 58:   VecDuplicateVecs(x2, n_timesteps, &x4ts);
 59:   for (i=0; i<n_timesteps; i++) {
 60:     PetscObjectSetName((PetscObject) x4ts[i], "x4");
 61:     VecSetRandom(x4ts[i], rand);
 62:     PetscViewerHDF5SetTimestep(viewer, i);
 63:     VecView(x4ts[i], viewer);
 64:   }
 65:   PetscViewerHDF5PopGroup(viewer);

 67:   PetscViewerHDF5PopTimestepping(viewer);
 68:   PetscViewerHDF5PopGroup(viewer);

 70:   PetscViewerDestroy(&viewer);
 71:   PetscRandomDestroy(&rand);

 73:   /* read and compare */
 74:   PetscViewerHDF5Open(comm, "ex19.h5", FILE_MODE_READ, &viewer);

 76:   VecDuplicate(x1, &x1r);
 77:   PetscObjectSetName((PetscObject) x1r, "x1");
 78:   VecLoad(x1r, viewer);
 79:   VecEqual(x1, x1r, &equal);
 80:   if (!equal) {
 81:     VecView(x1, PETSC_VIEWER_STDOUT_WORLD);
 82:     VecView(x1r, PETSC_VIEWER_STDOUT_WORLD);
 83:     SETERRQ(comm, PETSC_ERR_PLIB, "Error in HDF5 viewer: x1 != x1r");
 84:   }

 86:   PetscViewerHDF5PushGroup(viewer, "/testBlockSize");
 87:   VecDuplicate(x2, &x2r);
 88:   PetscObjectSetName((PetscObject) x2r, "x2");
 89:   VecLoad(x2r, viewer);
 90:   VecEqual(x2, x2r, &equal);
 91:   if (!equal) {
 92:     VecView(x2, PETSC_VIEWER_STDOUT_WORLD);
 93:     VecView(x2r, PETSC_VIEWER_STDOUT_WORLD);
 94:     SETERRQ(comm, PETSC_ERR_PLIB, "Error in HDF5 viewer: x2 != x2r");
 95:   }
 96:   PetscViewerHDF5PopGroup(viewer);

 98:   PetscViewerHDF5PushGroup(viewer, "/testTimestep");
 99:   PetscViewerHDF5PushTimestepping(viewer);

101:   VecDuplicate(x1, &x3r);
102:   PetscObjectSetName((PetscObject) x3r, "x3");
103:   for (i=0; i<n_timesteps; i++) {
104:     PetscViewerHDF5SetTimestep(viewer, i);
105:     VecLoad(x3r, viewer);
106:     VecEqual(x3r, x3ts[i], &equal);
107:     if (!equal) {
108:       VecView(x3r, PETSC_VIEWER_STDOUT_WORLD);
109:       VecView(x3ts[i], PETSC_VIEWER_STDOUT_WORLD);
110:       SETERRQ(comm, PETSC_ERR_PLIB, "Error in HDF5 viewer: x3ts[%" PetscInt_FMT "] != x3r", i);
111:     }
112:   }

114:   PetscViewerHDF5PushGroup(viewer, "testBlockSize");
115:   VecDuplicate(x2, &x4r);
116:   PetscObjectSetName((PetscObject) x4r, "x4");
117:   PetscViewerHDF5SetTimestep(viewer, 0);
118:   for (i=0; i<n_timesteps; i++) {
119:     VecLoad(x4r, viewer);
120:     VecEqual(x4r, x4ts[i], &equal);
121:     if (!equal) {
122:       VecView(x4r, PETSC_VIEWER_STDOUT_WORLD);
123:       VecView(x4ts[i], PETSC_VIEWER_STDOUT_WORLD);
124:       SETERRQ(comm, PETSC_ERR_PLIB, "Error in HDF5 viewer: x4ts[%" PetscInt_FMT "] != x4r", i);
125:     }
126:     PetscViewerHDF5IncrementTimestep(viewer);
127:   }
128:   PetscViewerHDF5PopGroup(viewer);

130:   PetscViewerHDF5PopTimestepping(viewer);
131:   PetscViewerHDF5PopGroup(viewer);

133:   /* cleanup */
134:   PetscViewerDestroy(&viewer);
135:   VecDestroy(&x1);
136:   VecDestroy(&x2);
137:   VecDestroyVecs(n_timesteps, &x3ts);
138:   VecDestroyVecs(n_timesteps, &x4ts);
139:   VecDestroy(&x1r);
140:   VecDestroy(&x2r);
141:   VecDestroy(&x3r);
142:   VecDestroy(&x4r);
143:   PetscFinalize();
144:   return 0;
145: }

147: /*TEST

149:      build:
150:        requires: hdf5

152:      test:
153:        nsize: {{1 2 3 4}}

155: TEST*/