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