xref: /petsc/src/dm/impls/plex/tests/ex15.c (revision d0609ced746bc51b019815ca91d747429db24893)
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 
239566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *) 0, help));
24c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
25c4762a1bSJed Brown 
26*d0609cedSBarry 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));
29*d0609cedSBarry 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     DM               dmDist;
46c4762a1bSJed Brown 
479566063dSJacob Faibussowitsch     PetscCall(DMPlexGetPartitioner(dm,&part));
489566063dSJacob Faibussowitsch     PetscCall(PetscPartitionerSetFromOptions(part));
499566063dSJacob Faibussowitsch     PetscCall(DMPlexDistribute(dm, 0, NULL, &dmDist));
50c4762a1bSJed Brown     if (dmDist) {
519566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&dm));
52c4762a1bSJed Brown       dm   = dmDist;
53c4762a1bSJed Brown     }
54c4762a1bSJed Brown   }
55c4762a1bSJed Brown 
569566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &v));
579566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) v, "V"));
589566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinates(dm, &coord));
599566063dSJacob Faibussowitsch   PetscCall(VecCopy(coord, v));
60c4762a1bSJed Brown 
61c4762a1bSJed Brown   if (verbose) {
62c4762a1bSJed Brown     PetscInt size, bs;
63c4762a1bSJed Brown 
649566063dSJacob Faibussowitsch     PetscCall(VecGetSize(v, &size));
659566063dSJacob Faibussowitsch     PetscCall(VecGetBlockSize(v, &bs));
669566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "==== original V in global ordering. size==%d\tblock size=%d\n", size, bs));
679566063dSJacob Faibussowitsch     PetscCall(VecView(v, PETSC_VIEWER_STDOUT_WORLD));
68c4762a1bSJed Brown   }
69c4762a1bSJed Brown 
709566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &nv));
719566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) nv, "NV"));
729566063dSJacob Faibussowitsch   PetscCall(DMPlexGlobalToNaturalBegin(dm, v, nv));
739566063dSJacob Faibussowitsch   PetscCall(DMPlexGlobalToNaturalEnd(dm, v, nv));
74c4762a1bSJed Brown 
75c4762a1bSJed Brown   if (verbose) {
76c4762a1bSJed Brown     PetscInt size, bs;
77c4762a1bSJed Brown 
789566063dSJacob Faibussowitsch     PetscCall(VecGetSize(nv, &size));
799566063dSJacob Faibussowitsch     PetscCall(VecGetBlockSize(nv, &bs));
809566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "====  V in natural ordering. size==%d\tblock size=%d\n", size, bs));
819566063dSJacob Faibussowitsch     PetscCall(VecView(nv, PETSC_VIEWER_STDOUT_WORLD));
82c4762a1bSJed Brown   }
83c4762a1bSJed Brown 
849566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(v, NULL, "-global_vec_view"));
85c4762a1bSJed Brown 
86c4762a1bSJed Brown   if (test_read) {
879566063dSJacob Faibussowitsch     PetscCall(DMCreateGlobalVector(dm, &rv));
889566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) rv, "V"));
89c4762a1bSJed Brown     /* Test native read */
909566063dSJacob Faibussowitsch     PetscCall(PetscViewerHDF5Open(comm, "V.h5", FILE_MODE_READ, &hdf5Viewer));
919566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(hdf5Viewer, PETSC_VIEWER_NATIVE));
929566063dSJacob Faibussowitsch     PetscCall(VecLoad(rv, hdf5Viewer));
939566063dSJacob Faibussowitsch     PetscCall(PetscViewerPopFormat(hdf5Viewer));
949566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&hdf5Viewer));
95c4762a1bSJed Brown     if (verbose) {
96c4762a1bSJed Brown       PetscInt size, bs;
97c4762a1bSJed Brown 
989566063dSJacob Faibussowitsch       PetscCall(VecGetSize(rv, &size));
999566063dSJacob Faibussowitsch       PetscCall(VecGetBlockSize(rv, &bs));
1009566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "==== Vector from file. size==%d\tblock size=%d\n", size, bs));
1019566063dSJacob Faibussowitsch       PetscCall(VecView(rv, PETSC_VIEWER_STDOUT_WORLD));
102c4762a1bSJed Brown     }
1039566063dSJacob Faibussowitsch     PetscCall(VecEqual(rv, v, &flg));
104c4762a1bSJed Brown     if (flg) {
1059566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "V and RV are equal\n"));
106c4762a1bSJed Brown     } else {
1079566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "V and RV are not equal\n\n"));
1089566063dSJacob Faibussowitsch       PetscCall(VecAXPY(rv, -1.0, v));
1099566063dSJacob Faibussowitsch       PetscCall(VecNorm(rv, NORM_INFINITY, &norm));
1109566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "diff norm is = %g\n", (double) norm));
1119566063dSJacob Faibussowitsch       PetscCall(VecView(rv, PETSC_VIEWER_STDOUT_WORLD));
112c4762a1bSJed Brown     }
113c4762a1bSJed Brown     /* Test raw read */
1149566063dSJacob Faibussowitsch     PetscCall(PetscViewerHDF5Open(comm, "V.h5", FILE_MODE_READ, &hdf5Viewer));
1159566063dSJacob Faibussowitsch     PetscCall(VecLoad(rv, hdf5Viewer));
1169566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&hdf5Viewer));
117c4762a1bSJed Brown     if (verbose) {
118c4762a1bSJed Brown       PetscInt size, bs;
119c4762a1bSJed Brown 
1209566063dSJacob Faibussowitsch       PetscCall(VecGetSize(rv, &size));
1219566063dSJacob Faibussowitsch       PetscCall(VecGetBlockSize(rv, &bs));
1229566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "==== Vector from file. size==%d\tblock size=%d\n", size, bs));
1239566063dSJacob Faibussowitsch       PetscCall(VecView(rv, PETSC_VIEWER_STDOUT_WORLD));
124c4762a1bSJed Brown     }
1259566063dSJacob Faibussowitsch     PetscCall(VecEqual(rv, nv, &flg));
126c4762a1bSJed Brown     if (flg) {
1279566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "NV and RV are equal\n"));
128c4762a1bSJed Brown     } else {
1299566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "NV and RV are not equal\n\n"));
1309566063dSJacob Faibussowitsch       PetscCall(VecAXPY(rv, -1.0, v));
1319566063dSJacob Faibussowitsch       PetscCall(VecNorm(rv, NORM_INFINITY, &norm));
1329566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "diff norm is = %g\n", (double) norm));
1339566063dSJacob Faibussowitsch       PetscCall(VecView(rv, PETSC_VIEWER_STDOUT_WORLD));
134c4762a1bSJed Brown     }
1359566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&rv));
136c4762a1bSJed Brown   }
1379566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&nv));
1389566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&v));
1399566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
1409566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
141b122ec5aSJacob Faibussowitsch   return 0;
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