xref: /petsc/src/dm/impls/plex/tests/ex15.c (revision 445019dbb80171c1132cad6620e62279b7cbfa97)
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;
9*445019dbSMatthew 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 
239566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *) 0, help));
24c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
25c4762a1bSJed Brown 
26d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD,"","Test Options","none");
279566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_read","Test reading from the HDF5 file","",PETSC_FALSE,&test_read,NULL));
289566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-verbose","print the Vecs","",PETSC_FALSE,&verbose,NULL));
29d0609cedSBarry Smith   PetscOptionsEnd();
30c4762a1bSJed Brown 
319566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, &dm));
329566063dSJacob Faibussowitsch   PetscCall(DMSetType(dm, DMPLEX));
339566063dSJacob Faibussowitsch   PetscCall(DMPlexDistributeSetDefault(dm, PETSC_FALSE));
349566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
359566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
369566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
37c4762a1bSJed Brown   numDof[0] = dim;
389566063dSJacob Faibussowitsch   PetscCall(DMSetNumFields(dm, numFields));
399566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateSection(dm, NULL, numComp, numDof, numBC, bcFields, bcPoints, NULL, NULL, &section));
409566063dSJacob Faibussowitsch   PetscCall(DMSetLocalSection(dm, section));
419566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&section));
429566063dSJacob Faibussowitsch   PetscCall(DMSetUseNatural(dm, PETSC_TRUE));
43c4762a1bSJed Brown   {
44c4762a1bSJed Brown     PetscPartitioner part;
45c4762a1bSJed Brown 
469566063dSJacob Faibussowitsch     PetscCall(DMPlexGetPartitioner(dm,&part));
479566063dSJacob Faibussowitsch     PetscCall(PetscPartitionerSetFromOptions(part));
48*445019dbSMatthew G. Knepley     PetscCall(DMPlexDistribute(dm, 0, NULL, &ddm));
49c4762a1bSJed Brown   }
50c4762a1bSJed Brown 
51*445019dbSMatthew G. Knepley   PetscCall(DMCreateGlobalVector(ddm, &v));
529566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) v, "V"));
53*445019dbSMatthew G. Knepley   PetscCall(DMGetCoordinates(ddm, &coord));
549566063dSJacob Faibussowitsch   PetscCall(VecCopy(coord, v));
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   if (verbose) {
57c4762a1bSJed Brown     PetscInt size, bs;
58c4762a1bSJed Brown 
599566063dSJacob Faibussowitsch     PetscCall(VecGetSize(v, &size));
609566063dSJacob Faibussowitsch     PetscCall(VecGetBlockSize(v, &bs));
619566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "==== original V in global ordering. size==%d\tblock size=%d\n", size, bs));
629566063dSJacob Faibussowitsch     PetscCall(VecView(v, PETSC_VIEWER_STDOUT_WORLD));
63c4762a1bSJed Brown   }
64c4762a1bSJed Brown 
659566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &nv));
669566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) nv, "NV"));
67*445019dbSMatthew G. Knepley   PetscCall(DMPlexGlobalToNaturalBegin(ddm, v, nv));
68*445019dbSMatthew G. Knepley   PetscCall(DMPlexGlobalToNaturalEnd(ddm, v, nv));
69c4762a1bSJed Brown 
70c4762a1bSJed Brown   if (verbose) {
71c4762a1bSJed Brown     PetscInt size, bs;
72c4762a1bSJed Brown 
739566063dSJacob Faibussowitsch     PetscCall(VecGetSize(nv, &size));
749566063dSJacob Faibussowitsch     PetscCall(VecGetBlockSize(nv, &bs));
759566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "====  V in natural ordering. size==%d\tblock size=%d\n", size, bs));
769566063dSJacob Faibussowitsch     PetscCall(VecView(nv, PETSC_VIEWER_STDOUT_WORLD));
77c4762a1bSJed Brown   }
78c4762a1bSJed Brown 
799566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(v, NULL, "-global_vec_view"));
80c4762a1bSJed Brown 
81*445019dbSMatthew G. Knepley   if (0 && test_read) {
82*445019dbSMatthew G. Knepley     PetscCall(DMCreateGlobalVector(ddm, &rv));
839566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) rv, "V"));
84c4762a1bSJed Brown     /* Test native read */
859566063dSJacob Faibussowitsch     PetscCall(PetscViewerHDF5Open(comm, "V.h5", FILE_MODE_READ, &hdf5Viewer));
869566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(hdf5Viewer, PETSC_VIEWER_NATIVE));
879566063dSJacob Faibussowitsch     PetscCall(VecLoad(rv, hdf5Viewer));
889566063dSJacob Faibussowitsch     PetscCall(PetscViewerPopFormat(hdf5Viewer));
899566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&hdf5Viewer));
90c4762a1bSJed Brown     if (verbose) {
91c4762a1bSJed Brown       PetscInt size, bs;
92c4762a1bSJed Brown 
939566063dSJacob Faibussowitsch       PetscCall(VecGetSize(rv, &size));
949566063dSJacob Faibussowitsch       PetscCall(VecGetBlockSize(rv, &bs));
959566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "==== Vector from file. size==%d\tblock size=%d\n", size, bs));
969566063dSJacob Faibussowitsch       PetscCall(VecView(rv, PETSC_VIEWER_STDOUT_WORLD));
97c4762a1bSJed Brown     }
989566063dSJacob Faibussowitsch     PetscCall(VecEqual(rv, v, &flg));
99c4762a1bSJed Brown     if (flg) {
1009566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "V and RV are equal\n"));
101c4762a1bSJed Brown     } else {
1029566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "V and RV are not equal\n\n"));
1039566063dSJacob Faibussowitsch       PetscCall(VecAXPY(rv, -1.0, v));
1049566063dSJacob Faibussowitsch       PetscCall(VecNorm(rv, NORM_INFINITY, &norm));
1059566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "diff norm is = %g\n", (double) norm));
1069566063dSJacob Faibussowitsch       PetscCall(VecView(rv, PETSC_VIEWER_STDOUT_WORLD));
107c4762a1bSJed Brown     }
108c4762a1bSJed Brown     /* Test raw read */
1099566063dSJacob Faibussowitsch     PetscCall(PetscViewerHDF5Open(comm, "V.h5", FILE_MODE_READ, &hdf5Viewer));
1109566063dSJacob Faibussowitsch     PetscCall(VecLoad(rv, hdf5Viewer));
1119566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&hdf5Viewer));
112c4762a1bSJed Brown     if (verbose) {
113c4762a1bSJed Brown       PetscInt size, bs;
114c4762a1bSJed Brown 
1159566063dSJacob Faibussowitsch       PetscCall(VecGetSize(rv, &size));
1169566063dSJacob Faibussowitsch       PetscCall(VecGetBlockSize(rv, &bs));
1179566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "==== Vector from file. size==%d\tblock size=%d\n", size, bs));
1189566063dSJacob Faibussowitsch       PetscCall(VecView(rv, PETSC_VIEWER_STDOUT_WORLD));
119c4762a1bSJed Brown     }
1209566063dSJacob Faibussowitsch     PetscCall(VecEqual(rv, nv, &flg));
121c4762a1bSJed Brown     if (flg) {
1229566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "NV and RV are equal\n"));
123c4762a1bSJed Brown     } else {
1249566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "NV and RV are not equal\n\n"));
1259566063dSJacob Faibussowitsch       PetscCall(VecAXPY(rv, -1.0, v));
1269566063dSJacob Faibussowitsch       PetscCall(VecNorm(rv, NORM_INFINITY, &norm));
1279566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "diff norm is = %g\n", (double) norm));
1289566063dSJacob Faibussowitsch       PetscCall(VecView(rv, PETSC_VIEWER_STDOUT_WORLD));
129c4762a1bSJed Brown     }
1309566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&rv));
131c4762a1bSJed Brown   }
1329566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&nv));
1339566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&v));
134*445019dbSMatthew G. Knepley   PetscCall(DMDestroy(&ddm));
1359566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
1369566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
137b122ec5aSJacob Faibussowitsch   return 0;
138c4762a1bSJed Brown }
139c4762a1bSJed Brown 
140c4762a1bSJed Brown /*TEST
141*445019dbSMatthew G. Knepley   #build:
142*445019dbSMatthew G. Knepley   #  requires: triangle hdf5
143*445019dbSMatthew G. Knepley   #test:
144*445019dbSMatthew G. Knepley   #  suffix: 0
145*445019dbSMatthew G. Knepley   #  requires: triangle hdf5
146*445019dbSMatthew G. Knepley   #  nsize: 2
147*445019dbSMatthew G. Knepley   #  args: -petscpartitioner_type simple -verbose -globaltonatural_sf_view
148*445019dbSMatthew G. Knepley   #test:
149*445019dbSMatthew G. Knepley   #  suffix: 1
150*445019dbSMatthew G. Knepley   #  requires: triangle hdf5
151*445019dbSMatthew G. Knepley   #  nsize: 2
152*445019dbSMatthew G. Knepley   #  args: -petscpartitioner_type simple -verbose -global_vec_view hdf5:V.h5:native -test_read
153c4762a1bSJed Brown 
154c4762a1bSJed Brown TEST*/
155