xref: /petsc/src/dm/impls/plex/tests/ex15.c (revision 030f984af8d8bb4c203755d35bded3c05b3d83ce)
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       numFields   = 1;
14   PetscInt       numBC       = 0;
15   PetscInt       numComp[1]  = {2};
16   PetscInt       numDof[3]   = {2, 0, 0};
17   PetscInt       bcFields[1] = {0};
18   IS             bcPoints[1] = {NULL};
19   PetscSection   section;
20   PetscReal      norm;
21   PetscInt       dim;
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 = PetscOptionsEnd();CHKERRQ(ierr);
31 
32   ierr = DMCreate(comm, &dm);CHKERRQ(ierr);
33   ierr = DMSetType(dm, DMPLEX);CHKERRQ(ierr);
34   ierr = DMSetFromOptions(dm);CHKERRQ(ierr);
35   ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr);
36   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
37   numDof[0] = dim;
38   ierr = DMSetNumFields(dm, numFields);CHKERRQ(ierr);
39   ierr = DMPlexCreateSection(dm, NULL, numComp, numDof, numBC, bcFields, bcPoints, NULL, NULL, &section);CHKERRQ(ierr);
40   ierr = DMSetLocalSection(dm, section);CHKERRQ(ierr);
41   ierr = PetscSectionDestroy(&section);CHKERRQ(ierr);
42   ierr = DMSetUseNatural(dm, PETSC_TRUE);CHKERRQ(ierr);
43   {
44     PetscPartitioner part;
45     DM               dmDist;
46 
47     ierr = DMPlexGetPartitioner(dm,&part);CHKERRQ(ierr);
48     ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
49     ierr = DMPlexDistribute(dm, 0, NULL, &dmDist);CHKERRQ(ierr);
50     if (dmDist) {
51       ierr = DMDestroy(&dm);CHKERRQ(ierr);
52       dm   = dmDist;
53     }
54   }
55 
56   ierr = DMCreateGlobalVector(dm, &v);CHKERRQ(ierr);
57   ierr = PetscObjectSetName((PetscObject) v, "V");CHKERRQ(ierr);
58   ierr = DMGetCoordinates(dm, &coord);CHKERRQ(ierr);
59   ierr = VecCopy(coord, v);CHKERRQ(ierr);
60 
61   if (verbose) {
62     PetscInt size, bs;
63 
64     ierr = VecGetSize(v, &size);CHKERRQ(ierr);
65     ierr = VecGetBlockSize(v, &bs);CHKERRQ(ierr);
66     ierr = PetscPrintf(PETSC_COMM_WORLD, "==== original V in global ordering. size==%d\tblock size=%d\n", size, bs);CHKERRQ(ierr);
67     ierr = VecView(v, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
68   }
69 
70   ierr = DMCreateGlobalVector(dm, &nv);CHKERRQ(ierr);
71   ierr = PetscObjectSetName((PetscObject) nv, "NV");CHKERRQ(ierr);
72   ierr = DMPlexGlobalToNaturalBegin(dm, v, nv);CHKERRQ(ierr);
73   ierr = DMPlexGlobalToNaturalEnd(dm, v, nv);CHKERRQ(ierr);
74 
75   if (verbose) {
76     PetscInt size, bs;
77 
78     ierr = VecGetSize(nv, &size);CHKERRQ(ierr);
79     ierr = VecGetBlockSize(nv, &bs);CHKERRQ(ierr);
80     ierr = PetscPrintf(PETSC_COMM_WORLD, "====  V in natural ordering. size==%d\tblock size=%d\n", size, bs);CHKERRQ(ierr);
81     ierr = VecView(nv, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
82   }
83 
84   ierr = VecViewFromOptions(v, NULL, "-global_vec_view");CHKERRQ(ierr);
85 
86   if (test_read) {
87     ierr = DMCreateGlobalVector(dm, &rv);CHKERRQ(ierr);
88     ierr = PetscObjectSetName((PetscObject) rv, "V");CHKERRQ(ierr);
89     /* Test native read */
90     ierr = PetscViewerHDF5Open(comm, "V.h5", FILE_MODE_READ, &hdf5Viewer);CHKERRQ(ierr);
91     ierr = PetscViewerPushFormat(hdf5Viewer, PETSC_VIEWER_NATIVE);CHKERRQ(ierr);
92     ierr = VecLoad(rv, hdf5Viewer);CHKERRQ(ierr);
93     ierr = PetscViewerPopFormat(hdf5Viewer);CHKERRQ(ierr);
94     ierr = PetscViewerDestroy(&hdf5Viewer);CHKERRQ(ierr);
95     if (verbose) {
96       PetscInt size, bs;
97 
98       ierr = VecGetSize(rv, &size);CHKERRQ(ierr);
99       ierr = VecGetBlockSize(rv, &bs);CHKERRQ(ierr);
100       ierr = PetscPrintf(PETSC_COMM_WORLD, "==== Vector from file. size==%d\tblock size=%d\n", size, bs);CHKERRQ(ierr);
101       ierr = VecView(rv, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
102     }
103     ierr = VecEqual(rv, v, &flg);CHKERRQ(ierr);
104     if (flg) {
105       ierr = PetscPrintf(PETSC_COMM_WORLD, "V and RV are equal\n");CHKERRQ(ierr);
106     } else {
107       ierr = PetscPrintf(PETSC_COMM_WORLD, "V and RV are not equal\n\n");CHKERRQ(ierr);
108       ierr = VecAXPY(rv, -1.0, v);CHKERRQ(ierr);
109       ierr = VecNorm(rv, NORM_INFINITY, &norm);CHKERRQ(ierr);
110       ierr = PetscPrintf(PETSC_COMM_WORLD, "diff norm is = %g\n", (double) norm);CHKERRQ(ierr);
111       ierr = VecView(rv, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
112     }
113     /* Test raw read */
114     ierr = PetscViewerHDF5Open(comm, "V.h5", FILE_MODE_READ, &hdf5Viewer);CHKERRQ(ierr);
115     ierr = VecLoad(rv, hdf5Viewer);CHKERRQ(ierr);
116     ierr = PetscViewerDestroy(&hdf5Viewer);CHKERRQ(ierr);
117     if (verbose) {
118       PetscInt size, bs;
119 
120       ierr = VecGetSize(rv, &size);CHKERRQ(ierr);
121       ierr = VecGetBlockSize(rv, &bs);CHKERRQ(ierr);
122       ierr = PetscPrintf(PETSC_COMM_WORLD, "==== Vector from file. size==%d\tblock size=%d\n", size, bs);CHKERRQ(ierr);
123       ierr = VecView(rv, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
124     }
125     ierr = VecEqual(rv, nv, &flg);CHKERRQ(ierr);
126     if (flg) {
127       ierr = PetscPrintf(PETSC_COMM_WORLD, "NV and RV are equal\n");CHKERRQ(ierr);
128     } else {
129       ierr = PetscPrintf(PETSC_COMM_WORLD, "NV and RV are not equal\n\n");CHKERRQ(ierr);
130       ierr = VecAXPY(rv, -1.0, v);CHKERRQ(ierr);
131       ierr = VecNorm(rv, NORM_INFINITY, &norm);CHKERRQ(ierr);
132       ierr = PetscPrintf(PETSC_COMM_WORLD, "diff norm is = %g\n", (double) norm);CHKERRQ(ierr);
133       ierr = VecView(rv, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
134     }
135     ierr = VecDestroy(&rv);CHKERRQ(ierr);
136   }
137   ierr = VecDestroy(&nv);CHKERRQ(ierr);
138   ierr = VecDestroy(&v);CHKERRQ(ierr);
139   ierr = DMDestroy(&dm);CHKERRQ(ierr);
140   ierr = PetscFinalize();
141   return ierr;
142 }
143 
144 /*TEST
145   build:
146     requires: triangle hdf5
147   test:
148     suffix: 0
149     requires: triangle hdf5
150     nsize: 2
151     args: -petscpartitioner_type simple -verbose -globaltonatural_sf_view
152   test:
153     suffix: 1
154     requires: triangle hdf5
155     nsize: 2
156     args: -petscpartitioner_type simple -verbose -global_vec_view hdf5:V.h5:native -test_read
157 
158 TEST*/
159