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