Actual source code: ex19.c

petsc-3.11.4 2019-09-28
Report Typos and Errors
  1: static char help[] = "Parallel HDF5 Vec Viewing.\n\n";

  3: /*T
  4:    Concepts: vectors^viewing
  5:    Concepts: viewers^hdf5
  6:    Processors: n
  7: T*/

  9:  #include <petscvec.h>
 10: #include <petscviewerhdf5.h>

 12: int main(int argc,char **argv)
 13: {
 14:   Vec            x1, x2, y1, y2, y3, y4;
 15:   PetscViewer    viewer;
 16:   PetscMPIInt    rank;
 17:   PetscInt       i, nlocal, n = 6;
 18:   PetscScalar    *array;
 19:   PetscBool      equal;

 23:   PetscInitialize(&argc, &argv, (char*) 0, help);if (ierr) return ierr;
 24:   MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
 25:   PetscOptionsGetInt(NULL,NULL, "-n", &n, NULL);

 27:   VecCreate(PETSC_COMM_WORLD, &x1);
 28:   PetscObjectSetName((PetscObject) x1, "TestVec");
 29:   VecSetSizes(x1, PETSC_DECIDE, n);
 30:   VecSetFromOptions(x1);

 32:   /* initialize x1 */
 33:   VecGetLocalSize(x1, &nlocal);
 34:   VecGetArray(x1, &array);
 35:   for (i = 0; i < nlocal; i++) array[i] = rank + 1;
 36:   VecRestoreArray(x1, &array);

 38:   VecCreate(PETSC_COMM_WORLD, &x2);
 39:   PetscObjectSetName((PetscObject) x2, "TestVec2");
 40:   VecSetSizes(x2, PETSC_DECIDE, n);
 41:   VecSetBlockSize(x2, 2);
 42:   VecSetFromOptions(x2);

 44:   /* initialize x2 */
 45:   VecGetLocalSize(x2, &nlocal);
 46:   VecGetArray(x2, &array);
 47:   for (i = 0; i < nlocal; i++) array[i] = rank + 1;
 48:   VecRestoreArray(x2, &array);

 50:   PetscViewerHDF5Open(PETSC_COMM_WORLD, "ex19.h5", FILE_MODE_WRITE, &viewer);
 51:   VecView(x1, viewer);
 52:   PetscViewerHDF5PushGroup(viewer, "/testBlockSize");
 53:   VecView(x2, viewer);
 54:   PetscViewerHDF5PushGroup(viewer, "/testTimestep");
 55:   PetscViewerHDF5SetTimestep(viewer, 0);
 56:   VecView(x2, viewer);
 57:   PetscViewerHDF5SetTimestep(viewer, 1);
 58:   VecView(x2, viewer);
 59:   PetscViewerHDF5PopGroup(viewer);
 60:   PetscViewerHDF5PopGroup(viewer);
 61:   PetscViewerDestroy(&viewer);

 63:   VecCreate(PETSC_COMM_WORLD, &y1);
 64:   PetscObjectSetName((PetscObject) y1, "TestVec");
 65:   VecSetSizes(y1, PETSC_DECIDE, n);
 66:   VecSetFromOptions(y1);

 68:   VecCreate(PETSC_COMM_WORLD,&y2);
 69:   VecSetBlockSize(y2, 2);
 70:   PetscObjectSetName((PetscObject) y2, "TestVec2");

 72:   VecCreate(PETSC_COMM_WORLD,&y3);
 73:   VecSetBlockSize(y3, 2);
 74:   PetscObjectSetName((PetscObject) y3, "TestVec2");

 76:   VecCreate(PETSC_COMM_WORLD,&y4);
 77:   VecSetBlockSize(y4, 2);
 78:   PetscObjectSetName((PetscObject) y4, "TestVec2");

 80:   PetscViewerHDF5Open(PETSC_COMM_WORLD, "ex19.h5", FILE_MODE_READ, &viewer);
 81:   VecLoad(y1, viewer);

 83:   PetscViewerHDF5PushGroup(viewer, "/testBlockSize");
 84:   VecLoad(y2, viewer);
 85:   PetscViewerHDF5PushGroup(viewer, "/testTimestep");
 86:   PetscViewerHDF5SetTimestep(viewer, 0);
 87:   VecLoad(y3, viewer);
 88:   PetscViewerHDF5SetTimestep(viewer, 1);
 89:   VecLoad(y4, viewer);
 90:   PetscViewerHDF5PopGroup(viewer);
 91:   PetscViewerHDF5PopGroup(viewer);
 92:   PetscViewerDestroy(&viewer);

 94:   VecEqual(x1, y1, &equal);
 95:   if (!equal) {
 96:     VecView(x1, PETSC_VIEWER_STDOUT_WORLD);
 97:     VecView(y1, PETSC_VIEWER_STDOUT_WORLD);
 98:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB, "Error in HDF5 viewer");
 99:   }
100:   VecEqual(x2, y2, &equal);
101:   if (!equal) {
102:     VecView(x2, PETSC_VIEWER_STDOUT_WORLD);
103:     VecView(y2, PETSC_VIEWER_STDOUT_WORLD);
104:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB, "Error in HDF5 viewer");
105:   }

107:   VecDestroy(&x1);
108:   VecDestroy(&x2);
109:   VecDestroy(&y1);
110:   VecDestroy(&y2);
111:   VecDestroy(&y3);
112:   VecDestroy(&y4);
113:   PetscFinalize();
114:   return ierr;
115: }

117: /*TEST

119:      build:
120:        requires: hdf5

122:      test:
123:        nsize: 1

125:      test:
126:        nsize: 2
127:        suffix: 2

129:      test:
130:        nsize: 3
131:        suffix: 3

133:      test:
134:        nsize: 4
135:        suffix: 4

137: TEST*/