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