xref: /petsc/src/dm/impls/plex/tests/ex9.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1c4762a1bSJed Brown static char help[] = "Performance tests for DMPlex query operations\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscdmplex.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown typedef struct {
6c4762a1bSJed Brown   PetscInt  dim;             /* The topological mesh dimension */
7c4762a1bSJed Brown   PetscBool cellSimplex;     /* Flag for simplices */
8c4762a1bSJed Brown   PetscBool spectral;        /* Flag for spectral element layout */
9c4762a1bSJed Brown   PetscBool interpolate;     /* Flag for mesh interpolation */
10c4762a1bSJed Brown   PetscReal refinementLimit; /* Maximum volume of a refined cell */
11c4762a1bSJed Brown   PetscInt  numFields;       /* The number of section fields */
12c4762a1bSJed Brown   PetscInt *numComponents;   /* The number of field components */
13c4762a1bSJed Brown   PetscInt *numDof;          /* The dof signature for the section */
14c4762a1bSJed Brown   PetscBool reuseArray;      /* Pass in user allocated array to VecGetClosure() */
15c4762a1bSJed Brown   /* Test data */
16c4762a1bSJed Brown   PetscBool errors;            /* Treat failures as errors */
17c4762a1bSJed Brown   PetscInt  iterations;        /* The number of iterations for a query */
18c4762a1bSJed Brown   PetscReal maxConeTime;       /* Max time per run for DMPlexGetCone() */
19c4762a1bSJed Brown   PetscReal maxClosureTime;    /* Max time per run for DMPlexGetTransitiveClosure() */
20c4762a1bSJed Brown   PetscReal maxVecClosureTime; /* Max time per run for DMPlexVecGetClosure() */
21c4762a1bSJed Brown   PetscBool printTimes;        /* Print total times, do not check limits */
22c4762a1bSJed Brown } AppCtx;
23c4762a1bSJed Brown 
249371c9d4SSatish Balay static PetscErrorCode ProcessOptions(AppCtx *options) {
25c4762a1bSJed Brown   PetscInt  len;
26c4762a1bSJed Brown   PetscBool flg;
27c4762a1bSJed Brown 
28c4762a1bSJed Brown   PetscFunctionBegin;
29c4762a1bSJed Brown   options->dim               = 2;
30c4762a1bSJed Brown   options->cellSimplex       = PETSC_TRUE;
31c4762a1bSJed Brown   options->spectral          = PETSC_FALSE;
32c4762a1bSJed Brown   options->interpolate       = PETSC_FALSE;
33c4762a1bSJed Brown   options->refinementLimit   = 0.0;
34c4762a1bSJed Brown   options->numFields         = 0;
35c4762a1bSJed Brown   options->numComponents     = NULL;
36c4762a1bSJed Brown   options->numDof            = NULL;
37c4762a1bSJed Brown   options->reuseArray        = PETSC_FALSE;
38c4762a1bSJed Brown   options->errors            = PETSC_FALSE;
39c4762a1bSJed Brown   options->iterations        = 1;
40c4762a1bSJed Brown   options->maxConeTime       = 0.0;
41c4762a1bSJed Brown   options->maxClosureTime    = 0.0;
42c4762a1bSJed Brown   options->maxVecClosureTime = 0.0;
43c4762a1bSJed Brown   options->printTimes        = PETSC_FALSE;
44c4762a1bSJed Brown 
45d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_SELF, "", "Meshing Problem Options", "DMPLEX");
469566063dSJacob Faibussowitsch   PetscCall(PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex9.c", options->dim, &options->dim, NULL, 1, 3));
479566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-cellSimplex", "Flag for simplices", "ex9.c", options->cellSimplex, &options->cellSimplex, NULL));
489566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-spectral", "Flag for spectral element layout", "ex9.c", options->spectral, &options->spectral, NULL));
499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-interpolate", "Flag for mesh interpolation", "ex9.c", options->interpolate, &options->interpolate, NULL));
509566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-refinement_limit", "The maximum volume of a refined cell", "ex9.c", options->refinementLimit, &options->refinementLimit, NULL));
519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-num_fields", "The number of section fields", "ex9.c", options->numFields, &options->numFields, NULL, 0));
52c4762a1bSJed Brown   if (options->numFields) {
53c4762a1bSJed Brown     len = options->numFields;
549566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(len, &options->numComponents));
559566063dSJacob Faibussowitsch     PetscCall(PetscOptionsIntArray("-num_components", "The number of components per field", "ex9.c", options->numComponents, &len, &flg));
5663a3b9bcSJacob Faibussowitsch     PetscCheck(!flg || !(len != options->numFields), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Length of components array is %" PetscInt_FMT " should be %" PetscInt_FMT, len, options->numFields);
57c4762a1bSJed Brown   }
58c4762a1bSJed Brown   len = (options->dim + 1) * PetscMax(1, options->numFields);
599566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(len, &options->numDof));
609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsIntArray("-num_dof", "The dof signature for the section", "ex9.c", options->numDof, &len, &flg));
611dca8a05SBarry Smith   PetscCheck(!flg || len == (options->dim + 1) * PetscMax(1, options->numFields), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Length of dof array is %" PetscInt_FMT " should be %" PetscInt_FMT, len, (options->dim + 1) * PetscMax(1, options->numFields));
62c4762a1bSJed Brown 
63c4762a1bSJed Brown   /* We are specifying the scalar dof, so augment it for multiple components */
64c4762a1bSJed Brown   {
65c4762a1bSJed Brown     PetscInt f, d;
66c4762a1bSJed Brown 
67c4762a1bSJed Brown     for (f = 0; f < options->numFields; ++f) {
68c4762a1bSJed Brown       for (d = 0; d <= options->dim; ++d) options->numDof[f * (options->dim + 1) + d] *= options->numComponents[f];
69c4762a1bSJed Brown     }
70c4762a1bSJed Brown   }
71c4762a1bSJed Brown 
729566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-reuse_array", "Pass in user allocated array to VecGetClosure()", "ex9.c", options->reuseArray, &options->reuseArray, NULL));
739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-errors", "Treat failures as errors", "ex9.c", options->errors, &options->errors, NULL));
749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-iterations", "The number of iterations for a query", "ex9.c", options->iterations, &options->iterations, NULL, 0));
759566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-max_cone_time", "The maximum time per run for DMPlexGetCone()", "ex9.c", options->maxConeTime, &options->maxConeTime, NULL));
769566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-max_closure_time", "The maximum time per run for DMPlexGetTransitiveClosure()", "ex9.c", options->maxClosureTime, &options->maxClosureTime, NULL));
779566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-max_vec_closure_time", "The maximum time per run for DMPlexVecGetClosure()", "ex9.c", options->maxVecClosureTime, &options->maxVecClosureTime, NULL));
789566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-print_times", "Print total times, do not check limits", "ex9.c", options->printTimes, &options->printTimes, NULL));
79d0609cedSBarry Smith   PetscOptionsEnd();
80c4762a1bSJed Brown   PetscFunctionReturn(0);
81c4762a1bSJed Brown }
82c4762a1bSJed Brown 
839371c9d4SSatish Balay static PetscErrorCode CreateSimplex_2D(MPI_Comm comm, DM *newdm) {
84c4762a1bSJed Brown   DM          dm;
85c4762a1bSJed Brown   PetscInt    numPoints[2]        = {4, 2};
86c4762a1bSJed Brown   PetscInt    coneSize[6]         = {3, 3, 0, 0, 0, 0};
87c4762a1bSJed Brown   PetscInt    cones[6]            = {2, 3, 4, 5, 4, 3};
88c4762a1bSJed Brown   PetscInt    coneOrientations[6] = {0, 0, 0, 0, 0, 0};
89c4762a1bSJed Brown   PetscScalar vertexCoords[8]     = {-0.5, 0.5, 0.0, 0.0, 0.0, 1.0, 0.5, 0.5};
90c4762a1bSJed Brown   PetscInt    markerPoints[8]     = {2, 1, 3, 1, 4, 1, 5, 1};
91c4762a1bSJed Brown   PetscInt    dim = 2, depth = 1, p;
92c4762a1bSJed Brown 
93c4762a1bSJed Brown   PetscFunctionBegin;
949566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, &dm));
959566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)dm, "triangular"));
969566063dSJacob Faibussowitsch   PetscCall(DMSetType(dm, DMPLEX));
979566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(dm, dim));
989566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
99*48a46eb9SPierre Jolivet   for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1]));
100c4762a1bSJed Brown   *newdm = dm;
101c4762a1bSJed Brown   PetscFunctionReturn(0);
102c4762a1bSJed Brown }
103c4762a1bSJed Brown 
1049371c9d4SSatish Balay static PetscErrorCode CreateSimplex_3D(MPI_Comm comm, DM *newdm) {
105c4762a1bSJed Brown   DM          dm;
106c4762a1bSJed Brown   PetscInt    numPoints[2]        = {5, 2};
107c4762a1bSJed Brown   PetscInt    coneSize[23]        = {4, 4, 0, 0, 0, 0, 0};
108c4762a1bSJed Brown   PetscInt    cones[8]            = {2, 4, 3, 5, 3, 4, 6, 5};
109c4762a1bSJed Brown   PetscInt    coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
110c4762a1bSJed Brown   PetscScalar vertexCoords[15]    = {0.0, 0.0, -0.5, 0.0, -0.5, 0.0, 1.0, 0.0, 0.0, 0.0, 0.5, 0.0, 0.0, 0.0, 0.5};
111c4762a1bSJed Brown   PetscInt    markerPoints[10]    = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1};
112c4762a1bSJed Brown   PetscInt    dim = 3, depth = 1, p;
113c4762a1bSJed Brown 
114c4762a1bSJed Brown   PetscFunctionBegin;
1159566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, &dm));
1169566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)dm, "tetrahedral"));
1179566063dSJacob Faibussowitsch   PetscCall(DMSetType(dm, DMPLEX));
1189566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(dm, dim));
1199566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
120*48a46eb9SPierre Jolivet   for (p = 0; p < 5; ++p) PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1]));
121c4762a1bSJed Brown   *newdm = dm;
122c4762a1bSJed Brown   PetscFunctionReturn(0);
123c4762a1bSJed Brown }
124c4762a1bSJed Brown 
1259371c9d4SSatish Balay static PetscErrorCode CreateQuad_2D(MPI_Comm comm, DM *newdm) {
126c4762a1bSJed Brown   DM          dm;
127c4762a1bSJed Brown   PetscInt    numPoints[2]        = {6, 2};
128c4762a1bSJed Brown   PetscInt    coneSize[8]         = {4, 4, 0, 0, 0, 0, 0, 0};
129c4762a1bSJed Brown   PetscInt    cones[8]            = {2, 3, 4, 5, 3, 6, 7, 4};
130c4762a1bSJed Brown   PetscInt    coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
131c4762a1bSJed Brown   PetscScalar vertexCoords[12]    = {-0.5, 0.0, 0.0, 0.0, 0.0, 1.0, -0.5, 1.0, 0.5, 0.0, 0.5, 1.0};
132c4762a1bSJed Brown   PetscInt    markerPoints[12]    = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1, 7, 1};
133c4762a1bSJed Brown   PetscInt    dim = 2, depth = 1, p;
134c4762a1bSJed Brown 
135c4762a1bSJed Brown   PetscFunctionBegin;
1369566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, &dm));
1379566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)dm, "quadrilateral"));
1389566063dSJacob Faibussowitsch   PetscCall(DMSetType(dm, DMPLEX));
1399566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(dm, dim));
1409566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
141*48a46eb9SPierre Jolivet   for (p = 0; p < 6; ++p) PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1]));
142c4762a1bSJed Brown   *newdm = dm;
143c4762a1bSJed Brown   PetscFunctionReturn(0);
144c4762a1bSJed Brown }
145c4762a1bSJed Brown 
1469371c9d4SSatish Balay static PetscErrorCode CreateHex_3D(MPI_Comm comm, DM *newdm) {
147c4762a1bSJed Brown   DM          dm;
148c4762a1bSJed Brown   PetscInt    numPoints[2]         = {12, 2};
149c4762a1bSJed Brown   PetscInt    coneSize[14]         = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
150c4762a1bSJed Brown   PetscInt    cones[16]            = {2, 5, 4, 3, 6, 7, 8, 9, 3, 4, 11, 10, 7, 12, 13, 8};
151c4762a1bSJed Brown   PetscInt    coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
1529371c9d4SSatish Balay   PetscScalar vertexCoords[36]     = {-0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.0, -0.5, 1.0, 0.0, -0.5, 0.0, 1.0, 0.0, 0.0, 1.0, 0.0, 1.0, 1.0, -0.5, 1.0, 1.0, 0.5, 0.0, 0.0, 0.5, 1.0, 0.0, 0.5, 0.0, 1.0, 0.5, 1.0, 1.0};
153c4762a1bSJed Brown   PetscInt    markerPoints[24]     = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1, 7, 1, 8, 1, 9, 1, 10, 1, 11, 1, 12, 1, 13, 1};
154c4762a1bSJed Brown   PetscInt    dim = 3, depth = 1, p;
155c4762a1bSJed Brown 
156c4762a1bSJed Brown   PetscFunctionBegin;
1579566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, &dm));
1589566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)dm, "hexahedral"));
1599566063dSJacob Faibussowitsch   PetscCall(DMSetType(dm, DMPLEX));
1609566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(dm, dim));
1619566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
162*48a46eb9SPierre Jolivet   for (p = 0; p < 12; ++p) PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1]));
163c4762a1bSJed Brown   *newdm = dm;
164c4762a1bSJed Brown   PetscFunctionReturn(0);
165c4762a1bSJed Brown }
166c4762a1bSJed Brown 
1679371c9d4SSatish Balay static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *newdm) {
168c4762a1bSJed Brown   PetscInt  dim         = user->dim;
169c4762a1bSJed Brown   PetscBool cellSimplex = user->cellSimplex;
170c4762a1bSJed Brown 
171c4762a1bSJed Brown   PetscFunctionBegin;
172c4762a1bSJed Brown   switch (dim) {
173c4762a1bSJed Brown   case 2:
174c4762a1bSJed Brown     if (cellSimplex) {
1759566063dSJacob Faibussowitsch       PetscCall(CreateSimplex_2D(comm, newdm));
176c4762a1bSJed Brown     } else {
1779566063dSJacob Faibussowitsch       PetscCall(CreateQuad_2D(comm, newdm));
178c4762a1bSJed Brown     }
179c4762a1bSJed Brown     break;
180c4762a1bSJed Brown   case 3:
181c4762a1bSJed Brown     if (cellSimplex) {
1829566063dSJacob Faibussowitsch       PetscCall(CreateSimplex_3D(comm, newdm));
183c4762a1bSJed Brown     } else {
1849566063dSJacob Faibussowitsch       PetscCall(CreateHex_3D(comm, newdm));
185c4762a1bSJed Brown     }
186c4762a1bSJed Brown     break;
1879371c9d4SSatish Balay   default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %" PetscInt_FMT, dim);
188c4762a1bSJed Brown   }
189c4762a1bSJed Brown   if (user->refinementLimit > 0.0) {
190c4762a1bSJed Brown     DM          rdm;
191c4762a1bSJed Brown     const char *name;
192c4762a1bSJed Brown 
1939566063dSJacob Faibussowitsch     PetscCall(DMPlexSetRefinementUniform(*newdm, PETSC_FALSE));
1949566063dSJacob Faibussowitsch     PetscCall(DMPlexSetRefinementLimit(*newdm, user->refinementLimit));
1959566063dSJacob Faibussowitsch     PetscCall(DMRefine(*newdm, PETSC_COMM_SELF, &rdm));
1969566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName((PetscObject)*newdm, &name));
1979566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)rdm, name));
1989566063dSJacob Faibussowitsch     PetscCall(DMDestroy(newdm));
199c4762a1bSJed Brown     *newdm = rdm;
200c4762a1bSJed Brown   }
201c4762a1bSJed Brown   if (user->interpolate) {
202c4762a1bSJed Brown     DM idm;
203c4762a1bSJed Brown 
2049566063dSJacob Faibussowitsch     PetscCall(DMPlexInterpolate(*newdm, &idm));
2059566063dSJacob Faibussowitsch     PetscCall(DMDestroy(newdm));
206c4762a1bSJed Brown     *newdm = idm;
207c4762a1bSJed Brown   }
2089566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*newdm));
209c4762a1bSJed Brown   PetscFunctionReturn(0);
210c4762a1bSJed Brown }
211c4762a1bSJed Brown 
2129371c9d4SSatish Balay static PetscErrorCode TestCone(DM dm, AppCtx *user) {
213c4762a1bSJed Brown   PetscInt           numRuns, cStart, cEnd, c, i;
214c4762a1bSJed Brown   PetscReal          maxTimePerRun = user->maxConeTime;
215c4762a1bSJed Brown   PetscLogStage      stage;
216c4762a1bSJed Brown   PetscLogEvent      event;
217c4762a1bSJed Brown   PetscEventPerfInfo eventInfo;
218c4762a1bSJed Brown   MPI_Comm           comm;
219c4762a1bSJed Brown   PetscMPIInt        rank;
220c4762a1bSJed Brown 
221c4762a1bSJed Brown   PetscFunctionBegin;
2229566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
2239566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
2249566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("DMPlex Cone Test", &stage));
2259566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("Cone", PETSC_OBJECT_CLASSID, &event));
2269566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(stage));
2279566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
2289566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(event, 0, 0, 0, 0));
229c4762a1bSJed Brown   for (i = 0; i < user->iterations; ++i) {
230c4762a1bSJed Brown     for (c = cStart; c < cEnd; ++c) {
231c4762a1bSJed Brown       const PetscInt *cone;
232c4762a1bSJed Brown 
2339566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCone(dm, c, &cone));
234c4762a1bSJed Brown     }
235c4762a1bSJed Brown   }
2369566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(event, 0, 0, 0, 0));
2379566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePop());
238c4762a1bSJed Brown 
2399566063dSJacob Faibussowitsch   PetscCall(PetscLogEventGetPerfInfo(stage, event, &eventInfo));
240c4762a1bSJed Brown   numRuns = (cEnd - cStart) * user->iterations;
24163a3b9bcSJacob Faibussowitsch   PetscCheck(eventInfo.count == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event calls %d should be 1", eventInfo.count);
24263a3b9bcSJacob Faibussowitsch   PetscCheck((PetscInt)eventInfo.flops == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event flops %" PetscInt_FMT " should be 0", (PetscInt)eventInfo.flops);
243c4762a1bSJed Brown   if (user->printTimes) {
24463a3b9bcSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(comm, "[%d] Cones: %" PetscInt_FMT " Total time: %.3es Average time per cone: %.3es\n", rank, numRuns, eventInfo.time, eventInfo.time / numRuns));
2459566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT));
246c4762a1bSJed Brown   } else if (eventInfo.time > maxTimePerRun * numRuns) {
24763a3b9bcSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(comm, "[%d] Cones: %" PetscInt_FMT " Average time per cone: %gs standard: %gs\n", rank, numRuns, (double)(eventInfo.time / numRuns), (double)maxTimePerRun));
2489566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT));
24963a3b9bcSJacob Faibussowitsch     PetscCheck(!user->errors, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Average time for cone %g > standard %g", (double)(eventInfo.time / numRuns), (double)maxTimePerRun);
250c4762a1bSJed Brown   }
251c4762a1bSJed Brown   PetscFunctionReturn(0);
252c4762a1bSJed Brown }
253c4762a1bSJed Brown 
2549371c9d4SSatish Balay static PetscErrorCode TestTransitiveClosure(DM dm, AppCtx *user) {
255c4762a1bSJed Brown   PetscInt           numRuns, cStart, cEnd, c, i;
256c4762a1bSJed Brown   PetscReal          maxTimePerRun = user->maxClosureTime;
257c4762a1bSJed Brown   PetscLogStage      stage;
258c4762a1bSJed Brown   PetscLogEvent      event;
259c4762a1bSJed Brown   PetscEventPerfInfo eventInfo;
260c4762a1bSJed Brown   MPI_Comm           comm;
261c4762a1bSJed Brown   PetscMPIInt        rank;
262c4762a1bSJed Brown 
263c4762a1bSJed Brown   PetscFunctionBegin;
2649566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
2659566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
2669566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("DMPlex Transitive Closure Test", &stage));
2679566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("TransitiveClosure", PETSC_OBJECT_CLASSID, &event));
2689566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(stage));
2699566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
2709566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(event, 0, 0, 0, 0));
271c4762a1bSJed Brown   for (i = 0; i < user->iterations; ++i) {
272c4762a1bSJed Brown     for (c = cStart; c < cEnd; ++c) {
273c4762a1bSJed Brown       PetscInt *closure = NULL;
274c4762a1bSJed Brown       PetscInt  closureSize;
275c4762a1bSJed Brown 
2769566063dSJacob Faibussowitsch       PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
2779566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
278c4762a1bSJed Brown     }
279c4762a1bSJed Brown   }
2809566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(event, 0, 0, 0, 0));
2819566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePop());
282c4762a1bSJed Brown 
2839566063dSJacob Faibussowitsch   PetscCall(PetscLogEventGetPerfInfo(stage, event, &eventInfo));
284c4762a1bSJed Brown   numRuns = (cEnd - cStart) * user->iterations;
28563a3b9bcSJacob Faibussowitsch   PetscCheck(eventInfo.count == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event calls %d should be 1", eventInfo.count);
28663a3b9bcSJacob Faibussowitsch   PetscCheck((PetscInt)eventInfo.flops == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event flops %" PetscInt_FMT " should be 0", (PetscInt)eventInfo.flops);
287c4762a1bSJed Brown   if (user->printTimes) {
28863a3b9bcSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(comm, "[%d] Closures: %" PetscInt_FMT " Total time: %.3es Average time per cone: %.3es\n", rank, numRuns, eventInfo.time, eventInfo.time / numRuns));
2899566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT));
290c4762a1bSJed Brown   } else if (eventInfo.time > maxTimePerRun * numRuns) {
29163a3b9bcSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(comm, "[%d] Closures: %" PetscInt_FMT " Average time per cone: %gs standard: %gs\n", rank, numRuns, (double)(eventInfo.time / numRuns), (double)maxTimePerRun));
2929566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT));
29363a3b9bcSJacob Faibussowitsch     PetscCheck(!user->errors, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Average time for closure %g > standard %g", (double)(eventInfo.time / numRuns), (double)maxTimePerRun);
294c4762a1bSJed Brown   }
295c4762a1bSJed Brown   PetscFunctionReturn(0);
296c4762a1bSJed Brown }
297c4762a1bSJed Brown 
2989371c9d4SSatish Balay static PetscErrorCode TestVecClosure(DM dm, PetscBool useIndex, PetscBool useSpectral, AppCtx *user) {
299c4762a1bSJed Brown   PetscSection       s;
300c4762a1bSJed Brown   Vec                v;
301c4762a1bSJed Brown   PetscInt           numRuns, cStart, cEnd, c, i;
302c4762a1bSJed Brown   PetscScalar        tmpArray[64];
303c4762a1bSJed Brown   PetscScalar       *userArray     = user->reuseArray ? tmpArray : NULL;
304c4762a1bSJed Brown   PetscReal          maxTimePerRun = user->maxVecClosureTime;
305c4762a1bSJed Brown   PetscLogStage      stage;
306c4762a1bSJed Brown   PetscLogEvent      event;
307c4762a1bSJed Brown   PetscEventPerfInfo eventInfo;
308c4762a1bSJed Brown   MPI_Comm           comm;
309c4762a1bSJed Brown   PetscMPIInt        rank;
310c4762a1bSJed Brown 
311c4762a1bSJed Brown   PetscFunctionBegin;
3129566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
3139566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
314c4762a1bSJed Brown   if (useIndex) {
315c4762a1bSJed Brown     if (useSpectral) {
3169566063dSJacob Faibussowitsch       PetscCall(PetscLogStageRegister("DMPlex Vector Closure with Index Test", &stage));
3179566063dSJacob Faibussowitsch       PetscCall(PetscLogEventRegister("VecClosureInd", PETSC_OBJECT_CLASSID, &event));
318c4762a1bSJed Brown     } else {
3199566063dSJacob Faibussowitsch       PetscCall(PetscLogStageRegister("DMPlex Vector Spectral Closure with Index Test", &stage));
3209566063dSJacob Faibussowitsch       PetscCall(PetscLogEventRegister("VecClosureSpecInd", PETSC_OBJECT_CLASSID, &event));
321c4762a1bSJed Brown     }
322c4762a1bSJed Brown   } else {
323c4762a1bSJed Brown     if (useSpectral) {
3249566063dSJacob Faibussowitsch       PetscCall(PetscLogStageRegister("DMPlex Vector Spectral Closure Test", &stage));
3259566063dSJacob Faibussowitsch       PetscCall(PetscLogEventRegister("VecClosureSpec", PETSC_OBJECT_CLASSID, &event));
326c4762a1bSJed Brown     } else {
3279566063dSJacob Faibussowitsch       PetscCall(PetscLogStageRegister("DMPlex Vector Closure Test", &stage));
3289566063dSJacob Faibussowitsch       PetscCall(PetscLogEventRegister("VecClosure", PETSC_OBJECT_CLASSID, &event));
329c4762a1bSJed Brown     }
330c4762a1bSJed Brown   }
3319566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(stage));
3329566063dSJacob Faibussowitsch   PetscCall(DMSetNumFields(dm, user->numFields));
3339566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateSection(dm, NULL, user->numComponents, user->numDof, 0, NULL, NULL, NULL, NULL, &s));
3349566063dSJacob Faibussowitsch   PetscCall(DMSetLocalSection(dm, s));
3359566063dSJacob Faibussowitsch   if (useIndex) PetscCall(DMPlexCreateClosureIndex(dm, s));
3369566063dSJacob Faibussowitsch   if (useSpectral) PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, s));
3379566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&s));
3389566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
3399566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &v));
3409566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(event, 0, 0, 0, 0));
341c4762a1bSJed Brown   for (i = 0; i < user->iterations; ++i) {
342c4762a1bSJed Brown     for (c = cStart; c < cEnd; ++c) {
343c4762a1bSJed Brown       PetscScalar *closure     = userArray;
344c4762a1bSJed Brown       PetscInt     closureSize = 64;
345c4762a1bSJed Brown 
3469566063dSJacob Faibussowitsch       PetscCall(DMPlexVecGetClosure(dm, s, v, c, &closureSize, &closure));
3479566063dSJacob Faibussowitsch       if (!user->reuseArray) PetscCall(DMPlexVecRestoreClosure(dm, s, v, c, &closureSize, &closure));
348c4762a1bSJed Brown     }
349c4762a1bSJed Brown   }
3509566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(event, 0, 0, 0, 0));
3519566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &v));
3529566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePop());
353c4762a1bSJed Brown 
3549566063dSJacob Faibussowitsch   PetscCall(PetscLogEventGetPerfInfo(stage, event, &eventInfo));
355c4762a1bSJed Brown   numRuns = (cEnd - cStart) * user->iterations;
35663a3b9bcSJacob Faibussowitsch   PetscCheck(eventInfo.count == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event calls %d should be 1", eventInfo.count);
35763a3b9bcSJacob Faibussowitsch   PetscCheck((PetscInt)eventInfo.flops == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event flops %" PetscInt_FMT " should be 0", (PetscInt)eventInfo.flops);
358c4762a1bSJed Brown   if (user->printTimes || eventInfo.time > maxTimePerRun * numRuns) {
359c4762a1bSJed Brown     const char *title          = "VecClosures";
360c4762a1bSJed Brown     const char *titleIndex     = "VecClosures with Index";
361c4762a1bSJed Brown     const char *titleSpec      = "VecClosures Spectral";
362c4762a1bSJed Brown     const char *titleSpecIndex = "VecClosures Spectral with Index";
363c4762a1bSJed Brown 
364c4762a1bSJed Brown     if (user->printTimes) {
3659371c9d4SSatish Balay       PetscCall(PetscSynchronizedPrintf(comm, "[%d] %s: %" PetscInt_FMT " Total time: %.3es Average time per vector closure: %.3es\n", rank, useIndex ? (useSpectral ? titleSpecIndex : titleIndex) : (useSpectral ? titleSpec : title), numRuns,
3669371c9d4SSatish Balay                                         eventInfo.time, eventInfo.time / numRuns));
3679566063dSJacob Faibussowitsch       PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT));
368c4762a1bSJed Brown     } else {
3699371c9d4SSatish Balay       PetscCall(
3709371c9d4SSatish Balay         PetscSynchronizedPrintf(comm, "[%d] %s: %" PetscInt_FMT " Average time per vector closure: %gs standard: %gs\n", rank, useIndex ? (useSpectral ? titleSpecIndex : titleIndex) : (useSpectral ? titleSpec : title), numRuns, (double)(eventInfo.time / numRuns), (double)maxTimePerRun));
3719566063dSJacob Faibussowitsch       PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT));
37263a3b9bcSJacob Faibussowitsch       PetscCheck(!user->errors, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Average time for vector closure %g > standard %g", (double)(eventInfo.time / numRuns), (double)maxTimePerRun);
373c4762a1bSJed Brown     }
374c4762a1bSJed Brown   }
375c4762a1bSJed Brown   PetscFunctionReturn(0);
376c4762a1bSJed Brown }
377c4762a1bSJed Brown 
3789371c9d4SSatish Balay static PetscErrorCode CleanupContext(AppCtx *user) {
379c4762a1bSJed Brown   PetscFunctionBegin;
3809566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->numComponents));
3819566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->numDof));
382c4762a1bSJed Brown   PetscFunctionReturn(0);
383c4762a1bSJed Brown }
384c4762a1bSJed Brown 
3859371c9d4SSatish Balay int main(int argc, char **argv) {
386c4762a1bSJed Brown   DM     dm;
387c4762a1bSJed Brown   AppCtx user;
388c4762a1bSJed Brown 
389327415f7SBarry Smith   PetscFunctionBeginUser;
3909566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
3919566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(&user));
3929566063dSJacob Faibussowitsch   PetscCall(PetscLogDefaultBegin());
3939566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
3949566063dSJacob Faibussowitsch   PetscCall(TestCone(dm, &user));
3959566063dSJacob Faibussowitsch   PetscCall(TestTransitiveClosure(dm, &user));
3969566063dSJacob Faibussowitsch   PetscCall(TestVecClosure(dm, PETSC_FALSE, PETSC_FALSE, &user));
3979566063dSJacob Faibussowitsch   PetscCall(TestVecClosure(dm, PETSC_TRUE, PETSC_FALSE, &user));
398c4762a1bSJed Brown   if (!user.cellSimplex && user.spectral) {
3999566063dSJacob Faibussowitsch     PetscCall(TestVecClosure(dm, PETSC_FALSE, PETSC_TRUE, &user));
4009566063dSJacob Faibussowitsch     PetscCall(TestVecClosure(dm, PETSC_TRUE, PETSC_TRUE, &user));
401c4762a1bSJed Brown   }
4029566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
4039566063dSJacob Faibussowitsch   PetscCall(CleanupContext(&user));
4049566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
405b122ec5aSJacob Faibussowitsch   return 0;
406c4762a1bSJed Brown }
407c4762a1bSJed Brown 
408c4762a1bSJed Brown /*TEST
409c4762a1bSJed Brown 
410c4762a1bSJed Brown   build:
411dfd57a17SPierre Jolivet     requires: defined(PETSC_USE_LOG)
412c4762a1bSJed Brown 
413c4762a1bSJed Brown   # 2D Simplex P_1 scalar tests
414c4762a1bSJed Brown   testset:
415c4762a1bSJed Brown     args: -num_dof 1,0,0 -iterations 2 -print_times
416c4762a1bSJed Brown     test:
417c4762a1bSJed Brown       suffix: correctness_0
418c4762a1bSJed Brown     test:
419c4762a1bSJed Brown       suffix: correctness_1
420c4762a1bSJed Brown       args: -interpolate -dm_refine 2
421c4762a1bSJed Brown     test:
422c4762a1bSJed Brown       suffix: correctness_2
423c4762a1bSJed Brown       requires: triangle
424c4762a1bSJed Brown       args: -interpolate -refinement_limit 1.0e-5
425c4762a1bSJed Brown   test:
426c4762a1bSJed Brown     suffix: 0
427c4762a1bSJed Brown     TODO: Only for performance testing
428c4762a1bSJed Brown     args: -num_dof 1,0,0 -iterations 10000 -max_cone_time 1.1e-8 -max_closure_time 1.3e-7 -max_vec_closure_time 3.6e-7
429c4762a1bSJed Brown   test:
430c4762a1bSJed Brown     suffix: 1
431c4762a1bSJed Brown     requires: triangle
432c4762a1bSJed Brown     TODO: Only for performance testing
433c4762a1bSJed Brown     args: -refinement_limit 1.0e-5 -num_dof 1,0,0 -iterations 2 -max_cone_time 2.1e-8 -max_closure_time 1.5e-7 -max_vec_closure_time 3.6e-7
434c4762a1bSJed Brown   test:
435c4762a1bSJed Brown     suffix: 2
436c4762a1bSJed Brown     TODO: Only for performance testing
437c4762a1bSJed Brown     args: -num_fields 1 -num_components 1 -num_dof 1,0,0 -iterations 10000 -max_cone_time 1.1e-8 -max_closure_time 1.3e-7 -max_vec_closure_time 4.5e-7
438c4762a1bSJed Brown   test:
439c4762a1bSJed Brown     suffix: 3
440c4762a1bSJed Brown     requires: triangle
441c4762a1bSJed Brown     TODO: Only for performance testing
442c4762a1bSJed Brown     args: -refinement_limit 1.0e-5 -num_fields 1 -num_components 1 -num_dof 1,0,0 -iterations 2 -max_cone_time 2.1e-8 -max_closure_time 1.5e-7 -max_vec_closure_time 4.7e-7
443c4762a1bSJed Brown   test:
444c4762a1bSJed Brown     suffix: 4
445c4762a1bSJed Brown     TODO: Only for performance testing
446c4762a1bSJed Brown     args: -interpolate -num_dof 1,0,0 -iterations 10000 -max_cone_time 1.1e-8 -max_closure_time 6.5e-7 -max_vec_closure_time 1.0e-6
447c4762a1bSJed Brown   test:
448c4762a1bSJed Brown     suffix: 5
449c4762a1bSJed Brown     requires: triangle
450c4762a1bSJed Brown     TODO: Only for performance testing
451c4762a1bSJed Brown     args: -interpolate -refinement_limit 1.0e-4 -num_dof 1,0,0 -iterations 2 -max_cone_time 2.1e-8 -max_closure_time 6.5e-7 -max_vec_closure_time 1.0e-6
452c4762a1bSJed Brown   test:
453c4762a1bSJed Brown     suffix: 6
454c4762a1bSJed Brown     TODO: Only for performance testing
455c4762a1bSJed Brown     args: -interpolate -num_fields 1 -num_components 1 -num_dof 1,0,0 -iterations 10000 -max_cone_time 1.1e-8 -max_closure_time 6.5e-7 -max_vec_closure_time 1.1e-6
456c4762a1bSJed Brown   test:
457c4762a1bSJed Brown     suffix: 7
458c4762a1bSJed Brown     requires: triangle
459c4762a1bSJed Brown     TODO: Only for performance testing
460c4762a1bSJed Brown     args: -interpolate -refinement_limit 1.0e-4 -num_fields 1 -num_components 1 -num_dof 1,0,0 -iterations 2 -max_cone_time 2.1e-8 -max_closure_time 6.5e-7 -max_vec_closure_time 1.2e-6
461c4762a1bSJed Brown 
462c4762a1bSJed Brown   # 2D Simplex P_1 vector tests
463c4762a1bSJed Brown   # 2D Simplex P_2 scalar tests
464c4762a1bSJed Brown   # 2D Simplex P_2 vector tests
465c4762a1bSJed Brown   # 2D Simplex P_2/P_1 vector/scalar tests
466c4762a1bSJed Brown   # 2D Quad P_1 scalar tests
467c4762a1bSJed Brown   # 2D Quad P_1 vector tests
468c4762a1bSJed Brown   # 2D Quad P_2 scalar tests
469c4762a1bSJed Brown   # 2D Quad P_2 vector tests
470c4762a1bSJed Brown   # 3D Simplex P_1 scalar tests
471c4762a1bSJed Brown   # 3D Simplex P_1 vector tests
472c4762a1bSJed Brown   # 3D Simplex P_2 scalar tests
473c4762a1bSJed Brown   # 3D Simplex P_2 vector tests
474c4762a1bSJed Brown   # 3D Hex P_1 scalar tests
475c4762a1bSJed Brown   # 3D Hex P_1 vector tests
476c4762a1bSJed Brown   # 3D Hex P_2 scalar tests
477c4762a1bSJed Brown   # 3D Hex P_2 vector tests
478c4762a1bSJed Brown 
479c4762a1bSJed Brown TEST*/
480