xref: /petsc/src/dm/tutorials/ex9.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown static char help[] = "Demonstrates HDF5 vector input/ouput\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*T
4c4762a1bSJed Brown    Concepts: viewers
5c4762a1bSJed Brown    Concepts: HDF5
6c4762a1bSJed Brown    Processors: n
7c4762a1bSJed Brown T*/
8c4762a1bSJed Brown #include <petscsys.h>
9c4762a1bSJed Brown #include <petscdm.h>
10c4762a1bSJed Brown #include <petscdmda.h>
11c4762a1bSJed Brown #include <petscviewerhdf5.h>
12c4762a1bSJed Brown 
13c4762a1bSJed Brown int main(int argc,char **argv)
14c4762a1bSJed Brown {
15c4762a1bSJed Brown   PetscErrorCode ierr;
16c4762a1bSJed Brown   PetscViewer    viewer;
17c4762a1bSJed Brown   DM             da;
18c4762a1bSJed Brown   Vec            global,local,global2;
19c4762a1bSJed Brown   PetscMPIInt    rank;
20c4762a1bSJed Brown   PetscBool      flg;
21c4762a1bSJed Brown   PetscInt       ndof;
22c4762a1bSJed Brown 
23c4762a1bSJed Brown   /*
24c4762a1bSJed Brown     Every PETSc routine should begin with the PetscInitialize() routine.
25c4762a1bSJed Brown     argc, argv - These command line arguments are taken to extract the options
26c4762a1bSJed Brown                  supplied to PETSc and options supplied to MPI.
27c4762a1bSJed Brown     help       - When PETSc executable is invoked with the option -help,
28c4762a1bSJed Brown                  it prints the various options that can be applied at
29c4762a1bSJed Brown                  runtime.  The user can use the "help" variable place
30c4762a1bSJed Brown                  additional help messages in this printout.
31c4762a1bSJed Brown   */
32*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
33c4762a1bSJed Brown   /* Get number of DOF's from command line */
34c4762a1bSJed Brown   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"DMDA VecView/VecLoad example","");CHKERRQ(ierr);
35c4762a1bSJed Brown   {
36c4762a1bSJed Brown     ndof = 1;
37c4762a1bSJed Brown     PetscOptionsBoundedInt("-ndof","Number of DOF's in DMDA","",ndof,&ndof,NULL,1);
38c4762a1bSJed Brown   }
39c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
40c4762a1bSJed Brown 
41c4762a1bSJed Brown   /* Create a DMDA and an associated vector */
425f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,100,90,PETSC_DECIDE,PETSC_DECIDE,ndof,1,NULL,NULL,&da));
435f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da));
445f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
455f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(da,&global));
465f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateLocalVector(da,&local));
475f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(global,-1.0));
485f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,global,INSERT_VALUES,local));
495f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,global,INSERT_VALUES,local));
505f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
515f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScale(local,rank+1));
525f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLocalToGlobalBegin(da,local,ADD_VALUES,global));
535f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLocalToGlobalEnd(da,local,ADD_VALUES,global));
54c4762a1bSJed Brown 
55c4762a1bSJed Brown   /* Create the HDF5 viewer for writing */
565f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerHDF5Open(PETSC_COMM_WORLD,"hdf5output.h5",FILE_MODE_WRITE,&viewer));
575f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerSetFromOptions(viewer));
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   /* Write the Vec without one extra dimension for BS */
605f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerHDF5SetBaseDimension2(viewer, PETSC_FALSE));
615f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) global, "noBsDim"));
625f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(global,viewer));
63c4762a1bSJed Brown 
64c4762a1bSJed Brown   /* Write the Vec with one extra, 1-sized, dimension for BS */
655f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerHDF5SetBaseDimension2(viewer, PETSC_TRUE));
665f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) global, "bsDim"));
675f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(global,viewer));
68c4762a1bSJed Brown 
695f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&viewer));
705f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Barrier(PETSC_COMM_WORLD));
715f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(global,&global2));
72c4762a1bSJed Brown 
73c4762a1bSJed Brown   /* Create the HDF5 viewer for reading */
745f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerHDF5Open(PETSC_COMM_WORLD,"hdf5output.h5",FILE_MODE_READ,&viewer));
755f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerSetFromOptions(viewer));
76c4762a1bSJed Brown 
77c4762a1bSJed Brown   /* Load the Vec without the BS dim and compare */
785f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) global2, "noBsDim"));
795f80ce2aSJacob Faibussowitsch   CHKERRQ(VecLoad(global2,viewer));
80c4762a1bSJed Brown 
815f80ce2aSJacob Faibussowitsch   CHKERRQ(VecEqual(global,global2,&flg));
82c4762a1bSJed Brown   if (!flg) {
835f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Error: Vectors are not equal\n"));
84c4762a1bSJed Brown   }
85c4762a1bSJed Brown 
86c4762a1bSJed Brown   /* Load the Vec with one extra, 1-sized, BS dim and compare */
875f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) global2, "bsDim"));
885f80ce2aSJacob Faibussowitsch   CHKERRQ(VecLoad(global2,viewer));
895f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&viewer));
90c4762a1bSJed Brown 
915f80ce2aSJacob Faibussowitsch   CHKERRQ(VecEqual(global,global2,&flg));
92c4762a1bSJed Brown   if (!flg) {
935f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Error: Vectors are not equal\n"));
94c4762a1bSJed Brown   }
95c4762a1bSJed Brown 
96c4762a1bSJed Brown   /* clean up and exit */
975f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&local));
985f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&global));
995f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&global2));
1005f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
101c4762a1bSJed Brown 
102*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
103*b122ec5aSJacob Faibussowitsch   return 0;
104c4762a1bSJed Brown }
105c4762a1bSJed Brown 
106c4762a1bSJed Brown /*TEST
107c4762a1bSJed Brown 
108c4762a1bSJed Brown       build:
109c4762a1bSJed Brown          requires: hdf5
110c4762a1bSJed Brown 
111c4762a1bSJed Brown       test:
112c4762a1bSJed Brown          nsize: 4
113c4762a1bSJed Brown 
114c4762a1bSJed Brown       test:
115c4762a1bSJed Brown          nsize: 4
116c4762a1bSJed Brown          suffix: 2
117c4762a1bSJed Brown          args: -ndof 2
118c4762a1bSJed Brown          output_file: output/ex9_1.out
119c4762a1bSJed Brown 
120c4762a1bSJed Brown TEST*/
121