xref: /petsc/src/dm/impls/plex/tests/ex15.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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;
9445019dbSMatthew G. Knepley   DM             dm, ddm;
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 
23*327415f7SBarry Smith   PetscFunctionBeginUser;
249566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *) 0, 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;
46c4762a1bSJed Brown 
479566063dSJacob Faibussowitsch     PetscCall(DMPlexGetPartitioner(dm,&part));
489566063dSJacob Faibussowitsch     PetscCall(PetscPartitionerSetFromOptions(part));
49445019dbSMatthew G. Knepley     PetscCall(DMPlexDistribute(dm, 0, NULL, &ddm));
50c4762a1bSJed Brown   }
51c4762a1bSJed Brown 
52445019dbSMatthew G. Knepley   PetscCall(DMCreateGlobalVector(ddm, &v));
539566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) v, "V"));
54445019dbSMatthew G. Knepley   PetscCall(DMGetCoordinates(ddm, &coord));
559566063dSJacob Faibussowitsch   PetscCall(VecCopy(coord, v));
56c4762a1bSJed Brown 
57c4762a1bSJed Brown   if (verbose) {
58c4762a1bSJed Brown     PetscInt size, bs;
59c4762a1bSJed Brown 
609566063dSJacob Faibussowitsch     PetscCall(VecGetSize(v, &size));
619566063dSJacob Faibussowitsch     PetscCall(VecGetBlockSize(v, &bs));
629566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "==== original V in global ordering. size==%d\tblock size=%d\n", size, bs));
639566063dSJacob Faibussowitsch     PetscCall(VecView(v, PETSC_VIEWER_STDOUT_WORLD));
64c4762a1bSJed Brown   }
65c4762a1bSJed Brown 
669566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &nv));
679566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) nv, "NV"));
68445019dbSMatthew G. Knepley   PetscCall(DMPlexGlobalToNaturalBegin(ddm, v, nv));
69445019dbSMatthew G. Knepley   PetscCall(DMPlexGlobalToNaturalEnd(ddm, v, nv));
70c4762a1bSJed Brown 
71c4762a1bSJed Brown   if (verbose) {
72c4762a1bSJed Brown     PetscInt size, bs;
73c4762a1bSJed Brown 
749566063dSJacob Faibussowitsch     PetscCall(VecGetSize(nv, &size));
759566063dSJacob Faibussowitsch     PetscCall(VecGetBlockSize(nv, &bs));
769566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "====  V in natural ordering. size==%d\tblock size=%d\n", size, bs));
779566063dSJacob Faibussowitsch     PetscCall(VecView(nv, PETSC_VIEWER_STDOUT_WORLD));
78c4762a1bSJed Brown   }
79c4762a1bSJed Brown 
809566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(v, NULL, "-global_vec_view"));
81c4762a1bSJed Brown 
82445019dbSMatthew G. Knepley   if (0 && test_read) {
83445019dbSMatthew G. Knepley     PetscCall(DMCreateGlobalVector(ddm, &rv));
849566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) rv, "V"));
85c4762a1bSJed Brown     /* Test native read */
869566063dSJacob Faibussowitsch     PetscCall(PetscViewerHDF5Open(comm, "V.h5", FILE_MODE_READ, &hdf5Viewer));
879566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(hdf5Viewer, PETSC_VIEWER_NATIVE));
889566063dSJacob Faibussowitsch     PetscCall(VecLoad(rv, hdf5Viewer));
899566063dSJacob Faibussowitsch     PetscCall(PetscViewerPopFormat(hdf5Viewer));
909566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&hdf5Viewer));
91c4762a1bSJed Brown     if (verbose) {
92c4762a1bSJed Brown       PetscInt size, bs;
93c4762a1bSJed Brown 
949566063dSJacob Faibussowitsch       PetscCall(VecGetSize(rv, &size));
959566063dSJacob Faibussowitsch       PetscCall(VecGetBlockSize(rv, &bs));
969566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "==== Vector from file. size==%d\tblock size=%d\n", size, bs));
979566063dSJacob Faibussowitsch       PetscCall(VecView(rv, PETSC_VIEWER_STDOUT_WORLD));
98c4762a1bSJed Brown     }
999566063dSJacob Faibussowitsch     PetscCall(VecEqual(rv, v, &flg));
100c4762a1bSJed Brown     if (flg) {
1019566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "V and RV are equal\n"));
102c4762a1bSJed Brown     } else {
1039566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "V and RV are not equal\n\n"));
1049566063dSJacob Faibussowitsch       PetscCall(VecAXPY(rv, -1.0, v));
1059566063dSJacob Faibussowitsch       PetscCall(VecNorm(rv, NORM_INFINITY, &norm));
1069566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "diff norm is = %g\n", (double) norm));
1079566063dSJacob Faibussowitsch       PetscCall(VecView(rv, PETSC_VIEWER_STDOUT_WORLD));
108c4762a1bSJed Brown     }
109c4762a1bSJed Brown     /* Test raw read */
1109566063dSJacob Faibussowitsch     PetscCall(PetscViewerHDF5Open(comm, "V.h5", FILE_MODE_READ, &hdf5Viewer));
1119566063dSJacob Faibussowitsch     PetscCall(VecLoad(rv, hdf5Viewer));
1129566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&hdf5Viewer));
113c4762a1bSJed Brown     if (verbose) {
114c4762a1bSJed Brown       PetscInt size, bs;
115c4762a1bSJed Brown 
1169566063dSJacob Faibussowitsch       PetscCall(VecGetSize(rv, &size));
1179566063dSJacob Faibussowitsch       PetscCall(VecGetBlockSize(rv, &bs));
1189566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "==== Vector from file. size==%d\tblock size=%d\n", size, bs));
1199566063dSJacob Faibussowitsch       PetscCall(VecView(rv, PETSC_VIEWER_STDOUT_WORLD));
120c4762a1bSJed Brown     }
1219566063dSJacob Faibussowitsch     PetscCall(VecEqual(rv, nv, &flg));
122c4762a1bSJed Brown     if (flg) {
1239566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "NV and RV are equal\n"));
124c4762a1bSJed Brown     } else {
1259566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "NV and RV are not equal\n\n"));
1269566063dSJacob Faibussowitsch       PetscCall(VecAXPY(rv, -1.0, v));
1279566063dSJacob Faibussowitsch       PetscCall(VecNorm(rv, NORM_INFINITY, &norm));
1289566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "diff norm is = %g\n", (double) norm));
1299566063dSJacob Faibussowitsch       PetscCall(VecView(rv, PETSC_VIEWER_STDOUT_WORLD));
130c4762a1bSJed Brown     }
1319566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&rv));
132c4762a1bSJed Brown   }
1339566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&nv));
1349566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&v));
135445019dbSMatthew G. Knepley   PetscCall(DMDestroy(&ddm));
1369566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
1379566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
138b122ec5aSJacob Faibussowitsch   return 0;
139c4762a1bSJed Brown }
140c4762a1bSJed Brown 
141c4762a1bSJed Brown /*TEST
142445019dbSMatthew G. Knepley   #build:
143445019dbSMatthew G. Knepley   #  requires: triangle hdf5
144445019dbSMatthew G. Knepley   #test:
145445019dbSMatthew G. Knepley   #  suffix: 0
146445019dbSMatthew G. Knepley   #  requires: triangle hdf5
147445019dbSMatthew G. Knepley   #  nsize: 2
148445019dbSMatthew G. Knepley   #  args: -petscpartitioner_type simple -verbose -globaltonatural_sf_view
149445019dbSMatthew G. Knepley   #test:
150445019dbSMatthew G. Knepley   #  suffix: 1
151445019dbSMatthew G. Knepley   #  requires: triangle hdf5
152445019dbSMatthew G. Knepley   #  nsize: 2
153445019dbSMatthew G. Knepley   #  args: -petscpartitioner_type simple -verbose -global_vec_view hdf5:V.h5:native -test_read
154c4762a1bSJed Brown 
155c4762a1bSJed Brown TEST*/
156