xref: /petsc/src/dm/tutorials/ex9.c (revision d0609ced746bc51b019815ca91d747429db24893)
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   */
269566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
27c4762a1bSJed Brown   /* Get number of DOF's from command line */
28*d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"DMDA VecView/VecLoad example","");
29c4762a1bSJed Brown   {
30c4762a1bSJed Brown     ndof = 1;
31c4762a1bSJed Brown     PetscOptionsBoundedInt("-ndof","Number of DOF's in DMDA","",ndof,&ndof,NULL,1);
32c4762a1bSJed Brown   }
33*d0609cedSBarry Smith   PetscOptionsEnd();
34c4762a1bSJed Brown 
35c4762a1bSJed Brown   /* Create a DMDA and an associated vector */
369566063dSJacob 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));
379566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
389566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
399566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da,&global));
409566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(da,&local));
419566063dSJacob Faibussowitsch   PetscCall(VecSet(global,-1.0));
429566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,global,INSERT_VALUES,local));
439566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,global,INSERT_VALUES,local));
449566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
459566063dSJacob Faibussowitsch   PetscCall(VecScale(local,rank+1));
469566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(da,local,ADD_VALUES,global));
479566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(da,local,ADD_VALUES,global));
48c4762a1bSJed Brown 
49c4762a1bSJed Brown   /* Create the HDF5 viewer for writing */
509566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Open(PETSC_COMM_WORLD,"hdf5output.h5",FILE_MODE_WRITE,&viewer));
519566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetFromOptions(viewer));
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   /* Write the Vec without one extra dimension for BS */
549566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5SetBaseDimension2(viewer, PETSC_FALSE));
559566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) global, "noBsDim"));
569566063dSJacob Faibussowitsch   PetscCall(VecView(global,viewer));
57c4762a1bSJed Brown 
58c4762a1bSJed Brown   /* Write the Vec with one extra, 1-sized, dimension for BS */
599566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5SetBaseDimension2(viewer, PETSC_TRUE));
609566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) global, "bsDim"));
619566063dSJacob Faibussowitsch   PetscCall(VecView(global,viewer));
62c4762a1bSJed Brown 
639566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
649566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD));
659566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(global,&global2));
66c4762a1bSJed Brown 
67c4762a1bSJed Brown   /* Create the HDF5 viewer for reading */
689566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Open(PETSC_COMM_WORLD,"hdf5output.h5",FILE_MODE_READ,&viewer));
699566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetFromOptions(viewer));
70c4762a1bSJed Brown 
71c4762a1bSJed Brown   /* Load the Vec without the BS dim and compare */
729566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) global2, "noBsDim"));
739566063dSJacob Faibussowitsch   PetscCall(VecLoad(global2,viewer));
74c4762a1bSJed Brown 
759566063dSJacob Faibussowitsch   PetscCall(VecEqual(global,global2,&flg));
76c4762a1bSJed Brown   if (!flg) {
779566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error: Vectors are not equal\n"));
78c4762a1bSJed Brown   }
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   /* Load the Vec with one extra, 1-sized, BS dim and compare */
819566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) global2, "bsDim"));
829566063dSJacob Faibussowitsch   PetscCall(VecLoad(global2,viewer));
839566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
84c4762a1bSJed Brown 
859566063dSJacob Faibussowitsch   PetscCall(VecEqual(global,global2,&flg));
86c4762a1bSJed Brown   if (!flg) {
879566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error: Vectors are not equal\n"));
88c4762a1bSJed Brown   }
89c4762a1bSJed Brown 
90c4762a1bSJed Brown   /* clean up and exit */
919566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&local));
929566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&global));
939566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&global2));
949566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
95c4762a1bSJed Brown 
969566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
97b122ec5aSJacob Faibussowitsch   return 0;
98c4762a1bSJed Brown }
99c4762a1bSJed Brown 
100c4762a1bSJed Brown /*TEST
101c4762a1bSJed Brown 
102c4762a1bSJed Brown       build:
103c4762a1bSJed Brown          requires: hdf5
104c4762a1bSJed Brown 
105c4762a1bSJed Brown       test:
106c4762a1bSJed Brown          nsize: 4
107c4762a1bSJed Brown 
108c4762a1bSJed Brown       test:
109c4762a1bSJed Brown          nsize: 4
110c4762a1bSJed Brown          suffix: 2
111c4762a1bSJed Brown          args: -ndof 2
112c4762a1bSJed Brown          output_file: output/ex9_1.out
113c4762a1bSJed Brown 
114c4762a1bSJed Brown TEST*/
115