xref: /petsc/src/dm/impls/plex/tests/ex15.c (revision 732aec7a18f2199fb53bb9a2f3aef439a834ce31)
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 
main(int argc,char ** argv)6d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
7d71ae5a4SJacob Faibussowitsch {
8c4762a1bSJed Brown   MPI_Comm     comm;
9e2739ba6SAlexis Marboeuf   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 
23327415f7SBarry Smith   PetscFunctionBeginUser;
24*c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
25c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
26c4762a1bSJed Brown 
27d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, "", "Test Options", "none");
289566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_read", "Test reading from the HDF5 file", "", PETSC_FALSE, &test_read, NULL));
299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-verbose", "print the Vecs", "", PETSC_FALSE, &verbose, NULL));
30d0609cedSBarry Smith   PetscOptionsEnd();
31c4762a1bSJed Brown 
329566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, &dm));
339566063dSJacob Faibussowitsch   PetscCall(DMSetType(dm, DMPLEX));
349566063dSJacob Faibussowitsch   PetscCall(DMPlexDistributeSetDefault(dm, PETSC_FALSE));
359566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
369566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
379566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
38c4762a1bSJed Brown   numDof[0] = dim;
399566063dSJacob Faibussowitsch   PetscCall(DMSetNumFields(dm, numFields));
409566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateSection(dm, NULL, numComp, numDof, numBC, bcFields, bcPoints, NULL, NULL, &section));
419566063dSJacob Faibussowitsch   PetscCall(DMSetLocalSection(dm, section));
429566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&section));
439566063dSJacob Faibussowitsch   PetscCall(DMSetUseNatural(dm, PETSC_TRUE));
44c4762a1bSJed Brown   {
45c4762a1bSJed Brown     PetscPartitioner part;
46e2739ba6SAlexis Marboeuf     DM               dmDist;
47c4762a1bSJed Brown 
489566063dSJacob Faibussowitsch     PetscCall(DMPlexGetPartitioner(dm, &part));
499566063dSJacob Faibussowitsch     PetscCall(PetscPartitionerSetFromOptions(part));
50e2739ba6SAlexis Marboeuf     PetscCall(DMPlexDistribute(dm, 0, NULL, &dmDist));
51e2739ba6SAlexis Marboeuf     if (dmDist) {
52e2739ba6SAlexis Marboeuf       PetscCall(DMDestroy(&dm));
53e2739ba6SAlexis Marboeuf       dm = dmDist;
54e2739ba6SAlexis Marboeuf     }
55c4762a1bSJed Brown   }
56c4762a1bSJed Brown 
57e2739ba6SAlexis Marboeuf   PetscCall(DMCreateGlobalVector(dm, &v));
589566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)v, "V"));
59e2739ba6SAlexis Marboeuf   PetscCall(DMGetCoordinates(dm, &coord));
609566063dSJacob Faibussowitsch   PetscCall(VecCopy(coord, v));
61c4762a1bSJed Brown 
62c4762a1bSJed Brown   if (verbose) {
63c4762a1bSJed Brown     PetscInt size, bs;
64c4762a1bSJed Brown 
659566063dSJacob Faibussowitsch     PetscCall(VecGetSize(v, &size));
669566063dSJacob Faibussowitsch     PetscCall(VecGetBlockSize(v, &bs));
67f59864b2SJames Wright     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "==== original V in global ordering. size==%" PetscInt_FMT "\tblock size=%" PetscInt_FMT "\n", size, bs));
689566063dSJacob Faibussowitsch     PetscCall(VecView(v, PETSC_VIEWER_STDOUT_WORLD));
69c4762a1bSJed Brown   }
70c4762a1bSJed Brown 
71e2739ba6SAlexis Marboeuf   PetscCall(DMPlexCreateNaturalVector(dm, &nv));
729566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)nv, "NV"));
73e2739ba6SAlexis Marboeuf   PetscCall(DMPlexGlobalToNaturalBegin(dm, v, nv));
74e2739ba6SAlexis Marboeuf   PetscCall(DMPlexGlobalToNaturalEnd(dm, v, nv));
75c4762a1bSJed Brown 
76c4762a1bSJed Brown   if (verbose) {
77c4762a1bSJed Brown     PetscInt size, bs;
78c4762a1bSJed Brown 
799566063dSJacob Faibussowitsch     PetscCall(VecGetSize(nv, &size));
809566063dSJacob Faibussowitsch     PetscCall(VecGetBlockSize(nv, &bs));
81f59864b2SJames Wright     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "====  V in natural ordering. size==%" PetscInt_FMT "\tblock size=%" PetscInt_FMT "\n", size, bs));
829566063dSJacob Faibussowitsch     PetscCall(VecView(nv, PETSC_VIEWER_STDOUT_WORLD));
83c4762a1bSJed Brown   }
84c4762a1bSJed Brown 
859566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(v, NULL, "-global_vec_view"));
86c4762a1bSJed Brown 
8768b974d2SMatthew G. Knepley   if (test_read) {
88e2739ba6SAlexis Marboeuf     PetscCall(DMCreateGlobalVector(dm, &rv));
899566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)rv, "V"));
90c4762a1bSJed Brown     /* Test native read */
919566063dSJacob Faibussowitsch     PetscCall(PetscViewerHDF5Open(comm, "V.h5", FILE_MODE_READ, &hdf5Viewer));
929566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(hdf5Viewer, PETSC_VIEWER_NATIVE));
939566063dSJacob Faibussowitsch     PetscCall(VecLoad(rv, hdf5Viewer));
949566063dSJacob Faibussowitsch     PetscCall(PetscViewerPopFormat(hdf5Viewer));
959566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&hdf5Viewer));
96c4762a1bSJed Brown     if (verbose) {
97c4762a1bSJed Brown       PetscInt size, bs;
98c4762a1bSJed Brown 
999566063dSJacob Faibussowitsch       PetscCall(VecGetSize(rv, &size));
1009566063dSJacob Faibussowitsch       PetscCall(VecGetBlockSize(rv, &bs));
101f59864b2SJames Wright       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "==== Vector from file. size==%" PetscInt_FMT "\tblock size=%" PetscInt_FMT "\n", size, bs));
1029566063dSJacob Faibussowitsch       PetscCall(VecView(rv, PETSC_VIEWER_STDOUT_WORLD));
103c4762a1bSJed Brown     }
1049566063dSJacob Faibussowitsch     PetscCall(VecEqual(rv, v, &flg));
105c4762a1bSJed Brown     if (flg) {
1069566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "V and RV are equal\n"));
107c4762a1bSJed Brown     } else {
1089566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "V and RV are not equal\n\n"));
1099566063dSJacob Faibussowitsch       PetscCall(VecAXPY(rv, -1.0, v));
1109566063dSJacob Faibussowitsch       PetscCall(VecNorm(rv, NORM_INFINITY, &norm));
1119566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "diff norm is = %g\n", (double)norm));
1129566063dSJacob Faibussowitsch       PetscCall(VecView(rv, PETSC_VIEWER_STDOUT_WORLD));
113c4762a1bSJed Brown     }
114c4762a1bSJed Brown     /* Test raw read */
1159566063dSJacob Faibussowitsch     PetscCall(PetscViewerHDF5Open(comm, "V.h5", FILE_MODE_READ, &hdf5Viewer));
1169566063dSJacob Faibussowitsch     PetscCall(VecLoad(rv, hdf5Viewer));
1179566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&hdf5Viewer));
118c4762a1bSJed Brown     if (verbose) {
119c4762a1bSJed Brown       PetscInt size, bs;
120c4762a1bSJed Brown 
1219566063dSJacob Faibussowitsch       PetscCall(VecGetSize(rv, &size));
1229566063dSJacob Faibussowitsch       PetscCall(VecGetBlockSize(rv, &bs));
123f59864b2SJames Wright       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "==== Vector from file. size==%" PetscInt_FMT "\tblock size=%" PetscInt_FMT "\n", size, bs));
1249566063dSJacob Faibussowitsch       PetscCall(VecView(rv, PETSC_VIEWER_STDOUT_WORLD));
125c4762a1bSJed Brown     }
1269566063dSJacob Faibussowitsch     PetscCall(VecEqual(rv, nv, &flg));
127c4762a1bSJed Brown     if (flg) {
1289566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "NV and RV are equal\n"));
129c4762a1bSJed Brown     } else {
1309566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "NV and RV are not equal\n\n"));
1319566063dSJacob Faibussowitsch       PetscCall(VecAXPY(rv, -1.0, v));
1329566063dSJacob Faibussowitsch       PetscCall(VecNorm(rv, NORM_INFINITY, &norm));
1339566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "diff norm is = %g\n", (double)norm));
1349566063dSJacob Faibussowitsch       PetscCall(VecView(rv, PETSC_VIEWER_STDOUT_WORLD));
135c4762a1bSJed Brown     }
1369566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&rv));
137c4762a1bSJed Brown   }
1389566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&nv));
1399566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&v));
1409566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
1419566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
142b122ec5aSJacob Faibussowitsch   return 0;
143c4762a1bSJed Brown }
144c4762a1bSJed Brown 
145c4762a1bSJed Brown /*TEST
146e2739ba6SAlexis Marboeuf   build:
147e2739ba6SAlexis Marboeuf     requires: triangle hdf5
148e2739ba6SAlexis Marboeuf   test:
149e2739ba6SAlexis Marboeuf     suffix: 0
150e2739ba6SAlexis Marboeuf     requires: triangle hdf5
151e2739ba6SAlexis Marboeuf     nsize: 2
152e2739ba6SAlexis Marboeuf     args: -petscpartitioner_type simple -verbose -globaltonatural_sf_view
153e2739ba6SAlexis Marboeuf   test:
154e2739ba6SAlexis Marboeuf     suffix: 1
155e2739ba6SAlexis Marboeuf     requires: triangle hdf5
156e2739ba6SAlexis Marboeuf     nsize: 2
157e2739ba6SAlexis Marboeuf     args: -petscpartitioner_type simple -verbose -global_vec_view hdf5:V.h5:native -test_read
158c4762a1bSJed Brown 
159c4762a1bSJed Brown TEST*/
160