Actual source code: ex19.c
petsc-3.7.5 2017-01-01
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>
14: int main(int argc,char **argv)
15: {
16: Vec x1, x2, y1, y2, y3, y4;
17: PetscViewer viewer;
18: PetscMPIInt rank;
19: PetscInt i, nlocal, n = 6;
20: PetscScalar *array;
21: PetscBool equal;
25: PetscInitialize(&argc, &argv, (char*) 0, help);
26: MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
27: PetscOptionsGetInt(NULL,NULL, "-n", &n, NULL);
29: VecCreate(PETSC_COMM_WORLD, &x1);
30: PetscObjectSetName((PetscObject) x1, "TestVec");
31: VecSetSizes(x1, PETSC_DECIDE, n);
32: VecSetFromOptions(x1);
34: VecGetLocalSize(x1, &nlocal);
35: VecGetArray(x1, &array);
36: for (i = 0; i < nlocal; i++) array[i] = rank + 1;
37: VecRestoreArray(x1, &array);
38: VecAssemblyBegin(x1);
39: VecAssemblyEnd(x1);
41: VecCreate(PETSC_COMM_WORLD, &x2);
42: PetscObjectSetName((PetscObject) x2, "TestVec2");
43: VecSetSizes(x2, PETSC_DECIDE, n);
44: VecSetBlockSize(x2, 2);
45: VecSetFromOptions(x2);
46: VecCopy(x1, x2);
48: PetscViewerHDF5Open(PETSC_COMM_WORLD, "ex19.h5", FILE_MODE_WRITE, &viewer);
49: PetscViewerHDF5PushGroup(viewer, "/");
50: VecView(x1, viewer);
51: PetscViewerHDF5PushGroup(viewer, "/testBlockSize");
52: VecView(x2, viewer);
53: PetscViewerHDF5PushGroup(viewer, "/testTimestep");
54: PetscViewerHDF5SetTimestep(viewer, 0);
55: VecView(x2, viewer);
56: PetscViewerHDF5SetTimestep(viewer, 1);
57: VecView(x2, viewer);
58: PetscViewerHDF5PopGroup(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: PetscViewerHDF5PushGroup(viewer, "/");
82: VecLoad(y1, viewer);
84: PetscViewerHDF5PushGroup(viewer, "/testBlockSize");
85: VecLoad(y2, viewer);
86: PetscViewerHDF5PushGroup(viewer, "/testTimestep");
87: PetscViewerHDF5SetTimestep(viewer, 0);
88: VecLoad(y3, viewer);
89: PetscViewerHDF5SetTimestep(viewer, 1);
90: VecLoad(y4, viewer);
91: PetscViewerHDF5PopGroup(viewer);
92: PetscViewerHDF5PopGroup(viewer);
93: PetscViewerHDF5PopGroup(viewer);
94: PetscViewerDestroy(&viewer);
96: VecEqual(x1, y1, &equal);
97: if (!equal) {
98: VecView(x1, PETSC_VIEWER_STDOUT_WORLD);
99: VecView(y1, PETSC_VIEWER_STDOUT_WORLD);
100: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB, "Error in HDF5 viewer");
101: }
102: VecEqual(x2, y2, &equal);
103: if (!equal) {
104: VecView(x2, PETSC_VIEWER_STDOUT_WORLD);
105: VecView(y2, PETSC_VIEWER_STDOUT_WORLD);
106: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB, "Error in HDF5 viewer");
107: }
109: VecDestroy(&x1);
110: VecDestroy(&x2);
111: VecDestroy(&y1);
112: VecDestroy(&y2);
113: VecDestroy(&y3);
114: VecDestroy(&y4);
115: PetscFinalize();
116: return(0);
117: }