Actual source code: ex9.c
1: static char help[] = "Demonstrates HDF5 vector input/ouput\n\n";
3: #include <petscsys.h>
4: #include <petscdm.h>
5: #include <petscdmda.h>
6: #include <petscviewerhdf5.h>
8: int main(int argc,char **argv)
9: {
11: PetscViewer viewer;
12: DM da;
13: Vec global,local,global2;
14: PetscMPIInt rank;
15: PetscBool flg;
16: PetscInt ndof;
18: /*
19: Every PETSc routine should begin with the PetscInitialize() routine.
20: argc, argv - These command line arguments are taken to extract the options
21: supplied to PETSc and options supplied to MPI.
22: help - When PETSc executable is invoked with the option -help,
23: it prints the various options that can be applied at
24: runtime. The user can use the "help" variable place
25: additional help messages in this printout.
26: */
27: PetscInitialize(&argc,&argv,(char*)0,help);
28: /* Get number of DOF's from command line */
29: PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"DMDA VecView/VecLoad example","");
30: {
31: ndof = 1;
32: PetscOptionsBoundedInt("-ndof","Number of DOF's in DMDA","",ndof,&ndof,NULL,1);
33: }
34: PetscOptionsEnd();
36: /* Create a DMDA and an associated vector */
37: DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,100,90,PETSC_DECIDE,PETSC_DECIDE,ndof,1,NULL,NULL,&da);
38: DMSetFromOptions(da);
39: DMSetUp(da);
40: DMCreateGlobalVector(da,&global);
41: DMCreateLocalVector(da,&local);
42: VecSet(global,-1.0);
43: DMGlobalToLocalBegin(da,global,INSERT_VALUES,local);
44: DMGlobalToLocalEnd(da,global,INSERT_VALUES,local);
45: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
46: VecScale(local,rank+1);
47: DMLocalToGlobalBegin(da,local,ADD_VALUES,global);
48: DMLocalToGlobalEnd(da,local,ADD_VALUES,global);
50: /* Create the HDF5 viewer for writing */
51: PetscViewerHDF5Open(PETSC_COMM_WORLD,"hdf5output.h5",FILE_MODE_WRITE,&viewer);
52: PetscViewerSetFromOptions(viewer);
54: /* Write the Vec without one extra dimension for BS */
55: PetscViewerHDF5SetBaseDimension2(viewer, PETSC_FALSE);
56: PetscObjectSetName((PetscObject) global, "noBsDim");
57: VecView(global,viewer);
59: /* Write the Vec with one extra, 1-sized, dimension for BS */
60: PetscViewerHDF5SetBaseDimension2(viewer, PETSC_TRUE);
61: PetscObjectSetName((PetscObject) global, "bsDim");
62: VecView(global,viewer);
64: PetscViewerDestroy(&viewer);
65: MPI_Barrier(PETSC_COMM_WORLD);
66: VecDuplicate(global,&global2);
68: /* Create the HDF5 viewer for reading */
69: PetscViewerHDF5Open(PETSC_COMM_WORLD,"hdf5output.h5",FILE_MODE_READ,&viewer);
70: PetscViewerSetFromOptions(viewer);
72: /* Load the Vec without the BS dim and compare */
73: PetscObjectSetName((PetscObject) global2, "noBsDim");
74: VecLoad(global2,viewer);
76: VecEqual(global,global2,&flg);
77: if (!flg) {
78: PetscPrintf(PETSC_COMM_WORLD,"Error: Vectors are not equal\n");
79: }
81: /* Load the Vec with one extra, 1-sized, BS dim and compare */
82: PetscObjectSetName((PetscObject) global2, "bsDim");
83: VecLoad(global2,viewer);
84: PetscViewerDestroy(&viewer);
86: VecEqual(global,global2,&flg);
87: if (!flg) {
88: PetscPrintf(PETSC_COMM_WORLD,"Error: Vectors are not equal\n");
89: }
91: /* clean up and exit */
92: VecDestroy(&local);
93: VecDestroy(&global);
94: VecDestroy(&global2);
95: DMDestroy(&da);
97: PetscFinalize();
98: return 0;
99: }
101: /*TEST
103: build:
104: requires: hdf5
106: test:
107: nsize: 4
109: test:
110: nsize: 4
111: suffix: 2
112: args: -ndof 2
113: output_file: output/ex9_1.out
115: TEST*/