xref: /petsc/src/dm/impls/plex/tests/ex15.c (revision 6a98f8dc3f2c9149905a87dc2e9d0fedaf64e09a)
1 static char help[] = "An example of writing a global Vec from a DMPlex with HDF5 format.\n";
2 
3 #include <petscdmplex.h>
4 #include <petscviewerhdf5.h>
5 
6 int main(int argc, char **argv)
7 {
8   MPI_Comm       comm;
9   DM             dm;
10   Vec            v, nv, rv, coord;
11   PetscBool      test_read = PETSC_FALSE, verbose = PETSC_FALSE, flg;
12   PetscViewer    hdf5Viewer;
13   PetscInt       dim         = 2;
14   PetscInt       numFields   = 1;
15   PetscInt       numBC       = 0;
16   PetscInt       numComp[1]  = {2};
17   PetscInt       numDof[3]   = {2, 0, 0};
18   PetscInt       bcFields[1] = {0};
19   IS             bcPoints[1] = {NULL};
20   PetscSection   section;
21   PetscReal      norm;
22   PetscErrorCode ierr;
23 
24   ierr = PetscInitialize(&argc, &argv, (char *) 0, help);if (ierr) return ierr;
25   comm = PETSC_COMM_WORLD;
26 
27   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,"","Test Options","none");CHKERRQ(ierr);
28   ierr = PetscOptionsBool("-test_read","Test reading from the HDF5 file","",PETSC_FALSE,&test_read,NULL);CHKERRQ(ierr);
29   ierr = PetscOptionsBool("-verbose","print the Vecs","",PETSC_FALSE,&verbose,NULL);CHKERRQ(ierr);
30   ierr = PetscOptionsRangeInt("-dim","the dimension of the problem","",2,&dim,NULL,1,3);CHKERRQ(ierr);
31   ierr = PetscOptionsEnd();
32 
33   ierr = DMPlexCreateBoxMesh(comm, dim, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, &dm);CHKERRQ(ierr);
34   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
35   numDof[0] = dim;
36   ierr = DMSetNumFields(dm, numFields);CHKERRQ(ierr);
37   ierr = DMPlexCreateSection(dm, NULL, numComp, numDof, numBC, bcFields, bcPoints, NULL, NULL, &section);CHKERRQ(ierr);
38   ierr = DMSetLocalSection(dm, section);CHKERRQ(ierr);
39   ierr = PetscSectionDestroy(&section);CHKERRQ(ierr);
40   ierr = DMSetUseNatural(dm, PETSC_TRUE);CHKERRQ(ierr);
41   {
42     PetscPartitioner part;
43     DM               dmDist;
44 
45     ierr = DMPlexGetPartitioner(dm,&part);CHKERRQ(ierr);
46     ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
47     ierr = DMPlexDistribute(dm, 0, NULL, &dmDist);CHKERRQ(ierr);
48     if (dmDist) {
49       ierr = DMDestroy(&dm);CHKERRQ(ierr);
50       dm   = dmDist;
51     }
52   }
53 
54   ierr = DMCreateGlobalVector(dm, &v);CHKERRQ(ierr);
55   ierr = PetscObjectSetName((PetscObject) v, "V");CHKERRQ(ierr);
56   ierr = DMGetCoordinates(dm, &coord);CHKERRQ(ierr);
57   ierr = VecCopy(coord, v);CHKERRQ(ierr);
58 
59   if (verbose) {
60     PetscInt size, bs;
61 
62     ierr = VecGetSize(v, &size);CHKERRQ(ierr);
63     ierr = VecGetBlockSize(v, &bs);CHKERRQ(ierr);
64     ierr = PetscPrintf(PETSC_COMM_WORLD, "==== original V in global ordering. size==%d\tblock size=%d\n", size, bs);CHKERRQ(ierr);
65     ierr = VecView(v, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
66   }
67 
68   ierr = DMCreateGlobalVector(dm, &nv);CHKERRQ(ierr);
69   ierr = PetscObjectSetName((PetscObject) nv, "NV");CHKERRQ(ierr);
70   ierr = DMPlexGlobalToNaturalBegin(dm, v, nv);CHKERRQ(ierr);
71   ierr = DMPlexGlobalToNaturalEnd(dm, v, nv);CHKERRQ(ierr);
72 
73   if (verbose) {
74     PetscInt size, bs;
75 
76     ierr = VecGetSize(nv, &size);CHKERRQ(ierr);
77     ierr = VecGetBlockSize(nv, &bs);CHKERRQ(ierr);
78     ierr = PetscPrintf(PETSC_COMM_WORLD, "====  V in natural ordering. size==%d\tblock size=%d\n", size, bs);CHKERRQ(ierr);
79     ierr = VecView(nv, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
80   }
81 
82   ierr = VecViewFromOptions(v, NULL, "-global_vec_view");CHKERRQ(ierr);
83 
84   if (test_read) {
85     ierr = DMCreateGlobalVector(dm, &rv);CHKERRQ(ierr);
86     ierr = PetscObjectSetName((PetscObject) rv, "V");CHKERRQ(ierr);
87     /* Test native read */
88     ierr = PetscViewerHDF5Open(comm, "V.h5", FILE_MODE_READ, &hdf5Viewer);CHKERRQ(ierr);
89     ierr = PetscViewerPushFormat(hdf5Viewer, PETSC_VIEWER_NATIVE);CHKERRQ(ierr);
90     ierr = VecLoad(rv, hdf5Viewer);CHKERRQ(ierr);
91     ierr = PetscViewerPopFormat(hdf5Viewer);CHKERRQ(ierr);
92     ierr = PetscViewerDestroy(&hdf5Viewer);CHKERRQ(ierr);
93     if (verbose) {
94       PetscInt size, bs;
95 
96       ierr = VecGetSize(rv, &size);CHKERRQ(ierr);
97       ierr = VecGetBlockSize(rv, &bs);CHKERRQ(ierr);
98       ierr = PetscPrintf(PETSC_COMM_WORLD, "==== Vector from file. size==%d\tblock size=%d\n", size, bs);CHKERRQ(ierr);
99       ierr = VecView(rv, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
100     }
101     ierr = VecEqual(rv, v, &flg);CHKERRQ(ierr);
102     if (flg) {
103       ierr = PetscPrintf(PETSC_COMM_WORLD, "V and RV are equal\n");CHKERRQ(ierr);
104     } else {
105       ierr = PetscPrintf(PETSC_COMM_WORLD, "V and RV are not equal\n\n");CHKERRQ(ierr);
106       ierr = VecAXPY(rv, -1.0, v);CHKERRQ(ierr);
107       ierr = VecNorm(rv, NORM_INFINITY, &norm);CHKERRQ(ierr);
108       ierr = PetscPrintf(PETSC_COMM_WORLD, "diff norm is = %g\n", (double) norm);CHKERRQ(ierr);
109       ierr = VecView(rv, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
110     }
111     /* Test raw read */
112     ierr = PetscViewerHDF5Open(comm, "V.h5", FILE_MODE_READ, &hdf5Viewer);CHKERRQ(ierr);
113     ierr = VecLoad(rv, hdf5Viewer);CHKERRQ(ierr);
114     ierr = PetscViewerDestroy(&hdf5Viewer);CHKERRQ(ierr);
115     if (verbose) {
116       PetscInt size, bs;
117 
118       ierr = VecGetSize(rv, &size);CHKERRQ(ierr);
119       ierr = VecGetBlockSize(rv, &bs);CHKERRQ(ierr);
120       ierr = PetscPrintf(PETSC_COMM_WORLD, "==== Vector from file. size==%d\tblock size=%d\n", size, bs);CHKERRQ(ierr);
121       ierr = VecView(rv, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
122     }
123     ierr = VecEqual(rv, nv, &flg);CHKERRQ(ierr);
124     if (flg) {
125       ierr = PetscPrintf(PETSC_COMM_WORLD, "NV and RV are equal\n");CHKERRQ(ierr);
126     } else {
127       ierr = PetscPrintf(PETSC_COMM_WORLD, "NV and RV are not equal\n\n");CHKERRQ(ierr);
128       ierr = VecAXPY(rv, -1.0, v);CHKERRQ(ierr);
129       ierr = VecNorm(rv, NORM_INFINITY, &norm);CHKERRQ(ierr);
130       ierr = PetscPrintf(PETSC_COMM_WORLD, "diff norm is = %g\n", (double) norm);CHKERRQ(ierr);
131       ierr = VecView(rv, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
132     }
133     ierr = VecDestroy(&rv);CHKERRQ(ierr);
134   }
135   ierr = VecDestroy(&nv);CHKERRQ(ierr);
136   ierr = VecDestroy(&v);CHKERRQ(ierr);
137   ierr = DMDestroy(&dm);CHKERRQ(ierr);
138   ierr = PetscFinalize();
139   return ierr;
140 }
141 
142 /*TEST
143   build:
144     requires: triangle hdf5
145   test:
146     suffix: 0
147     requires: triangle hdf5
148     nsize: 2
149     args: -petscpartitioner_type simple -verbose -globaltonatural_sf_view
150   test:
151     suffix: 1
152     requires: triangle hdf5
153     nsize: 2
154     args: -petscpartitioner_type simple -verbose -global_vec_view hdf5:V.h5:native -test_read
155 
156 TEST*/
157