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