xref: /petsc/src/dm/tutorials/ex9.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1c4762a1bSJed Brown static char help[] = "Demonstrates HDF5 vector input/ouput\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscsys.h>
4c4762a1bSJed Brown #include <petscdm.h>
5c4762a1bSJed Brown #include <petscdmda.h>
6c4762a1bSJed Brown #include <petscviewerhdf5.h>
7c4762a1bSJed Brown 
8c4762a1bSJed Brown int main(int argc,char **argv)
9c4762a1bSJed Brown {
10c4762a1bSJed Brown   PetscViewer    viewer;
11c4762a1bSJed Brown   DM             da;
12c4762a1bSJed Brown   Vec            global,local,global2;
13c4762a1bSJed Brown   PetscMPIInt    rank;
14c4762a1bSJed Brown   PetscBool      flg;
15c4762a1bSJed Brown   PetscInt       ndof;
16c4762a1bSJed Brown 
17c4762a1bSJed Brown   /*
18c4762a1bSJed Brown     Every PETSc routine should begin with the PetscInitialize() routine.
19c4762a1bSJed Brown     argc, argv - These command line arguments are taken to extract the options
20c4762a1bSJed Brown                  supplied to PETSc and options supplied to MPI.
21c4762a1bSJed Brown     help       - When PETSc executable is invoked with the option -help,
22c4762a1bSJed Brown                  it prints the various options that can be applied at
23c4762a1bSJed Brown                  runtime.  The user can use the "help" variable place
24c4762a1bSJed Brown                  additional help messages in this printout.
25c4762a1bSJed Brown   */
26*327415f7SBarry Smith   PetscFunctionBeginUser;
279566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
28c4762a1bSJed Brown   /* Get number of DOF's from command line */
29d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"DMDA VecView/VecLoad example","");
30c4762a1bSJed Brown   {
31c4762a1bSJed Brown     ndof = 1;
32c4762a1bSJed Brown     PetscOptionsBoundedInt("-ndof","Number of DOF's in DMDA","",ndof,&ndof,NULL,1);
33c4762a1bSJed Brown   }
34d0609cedSBarry Smith   PetscOptionsEnd();
35c4762a1bSJed Brown 
36c4762a1bSJed Brown   /* Create a DMDA and an associated vector */
379566063dSJacob Faibussowitsch   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));
389566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
399566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
409566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da,&global));
419566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(da,&local));
429566063dSJacob Faibussowitsch   PetscCall(VecSet(global,-1.0));
439566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,global,INSERT_VALUES,local));
449566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,global,INSERT_VALUES,local));
459566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
469566063dSJacob Faibussowitsch   PetscCall(VecScale(local,rank+1));
479566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(da,local,ADD_VALUES,global));
489566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(da,local,ADD_VALUES,global));
49c4762a1bSJed Brown 
50c4762a1bSJed Brown   /* Create the HDF5 viewer for writing */
519566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Open(PETSC_COMM_WORLD,"hdf5output.h5",FILE_MODE_WRITE,&viewer));
529566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetFromOptions(viewer));
53c4762a1bSJed Brown 
54c4762a1bSJed Brown   /* Write the Vec without one extra dimension for BS */
559566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5SetBaseDimension2(viewer, PETSC_FALSE));
569566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) global, "noBsDim"));
579566063dSJacob Faibussowitsch   PetscCall(VecView(global,viewer));
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   /* Write the Vec with one extra, 1-sized, dimension for BS */
609566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5SetBaseDimension2(viewer, PETSC_TRUE));
619566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) global, "bsDim"));
629566063dSJacob Faibussowitsch   PetscCall(VecView(global,viewer));
63c4762a1bSJed Brown 
649566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
659566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD));
669566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(global,&global2));
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   /* Create the HDF5 viewer for reading */
699566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Open(PETSC_COMM_WORLD,"hdf5output.h5",FILE_MODE_READ,&viewer));
709566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetFromOptions(viewer));
71c4762a1bSJed Brown 
72c4762a1bSJed Brown   /* Load the Vec without the BS dim and compare */
739566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) global2, "noBsDim"));
749566063dSJacob Faibussowitsch   PetscCall(VecLoad(global2,viewer));
75c4762a1bSJed Brown 
769566063dSJacob Faibussowitsch   PetscCall(VecEqual(global,global2,&flg));
77c4762a1bSJed Brown   if (!flg) {
789566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error: Vectors are not equal\n"));
79c4762a1bSJed Brown   }
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   /* Load the Vec with one extra, 1-sized, BS dim and compare */
829566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) global2, "bsDim"));
839566063dSJacob Faibussowitsch   PetscCall(VecLoad(global2,viewer));
849566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
85c4762a1bSJed Brown 
869566063dSJacob Faibussowitsch   PetscCall(VecEqual(global,global2,&flg));
87c4762a1bSJed Brown   if (!flg) {
889566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error: Vectors are not equal\n"));
89c4762a1bSJed Brown   }
90c4762a1bSJed Brown 
91c4762a1bSJed Brown   /* clean up and exit */
929566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&local));
939566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&global));
949566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&global2));
959566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
96c4762a1bSJed Brown 
979566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
98b122ec5aSJacob Faibussowitsch   return 0;
99c4762a1bSJed Brown }
100c4762a1bSJed Brown 
101c4762a1bSJed Brown /*TEST
102c4762a1bSJed Brown 
103c4762a1bSJed Brown       build:
104c4762a1bSJed Brown          requires: hdf5
105c4762a1bSJed Brown 
106c4762a1bSJed Brown       test:
107c4762a1bSJed Brown          nsize: 4
108c4762a1bSJed Brown 
109c4762a1bSJed Brown       test:
110c4762a1bSJed Brown          nsize: 4
111c4762a1bSJed Brown          suffix: 2
112c4762a1bSJed Brown          args: -ndof 2
113c4762a1bSJed Brown          output_file: output/ex9_1.out
114c4762a1bSJed Brown 
115c4762a1bSJed Brown TEST*/
116