xref: /petsc/src/dm/impls/plex/tests/ex9.c (revision be37439ebbbdb2f81c3420c175a94aa72e59929c)
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 
ProcessOptions(AppCtx * options)24d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(AppCtx *options)
25d71ae5a4SJacob Faibussowitsch {
26c4762a1bSJed Brown   PetscInt  len;
27c4762a1bSJed Brown   PetscBool flg;
28c4762a1bSJed Brown 
29c4762a1bSJed Brown   PetscFunctionBegin;
30c4762a1bSJed Brown   options->dim               = 2;
31c4762a1bSJed Brown   options->cellSimplex       = PETSC_TRUE;
32c4762a1bSJed Brown   options->spectral          = PETSC_FALSE;
33c4762a1bSJed Brown   options->interpolate       = PETSC_FALSE;
34c4762a1bSJed Brown   options->refinementLimit   = 0.0;
35c4762a1bSJed Brown   options->numFields         = 0;
36c4762a1bSJed Brown   options->numComponents     = NULL;
37c4762a1bSJed Brown   options->numDof            = NULL;
38c4762a1bSJed Brown   options->reuseArray        = PETSC_FALSE;
39c4762a1bSJed Brown   options->errors            = PETSC_FALSE;
40c4762a1bSJed Brown   options->iterations        = 1;
41c4762a1bSJed Brown   options->maxConeTime       = 0.0;
42c4762a1bSJed Brown   options->maxClosureTime    = 0.0;
43c4762a1bSJed Brown   options->maxVecClosureTime = 0.0;
44c4762a1bSJed Brown   options->printTimes        = PETSC_FALSE;
45c4762a1bSJed Brown 
46d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_SELF, "", "Meshing Problem Options", "DMPLEX");
479566063dSJacob Faibussowitsch   PetscCall(PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex9.c", options->dim, &options->dim, NULL, 1, 3));
489566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-cellSimplex", "Flag for simplices", "ex9.c", options->cellSimplex, &options->cellSimplex, NULL));
499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-spectral", "Flag for spectral element layout", "ex9.c", options->spectral, &options->spectral, NULL));
509566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-interpolate", "Flag for mesh interpolation", "ex9.c", options->interpolate, &options->interpolate, NULL));
519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-refinement_limit", "The maximum volume of a refined cell", "ex9.c", options->refinementLimit, &options->refinementLimit, NULL));
529566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-num_fields", "The number of section fields", "ex9.c", options->numFields, &options->numFields, NULL, 0));
53c4762a1bSJed Brown   if (options->numFields) {
54c4762a1bSJed Brown     len = options->numFields;
559566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(len, &options->numComponents));
569566063dSJacob Faibussowitsch     PetscCall(PetscOptionsIntArray("-num_components", "The number of components per field", "ex9.c", options->numComponents, &len, &flg));
5763a3b9bcSJacob 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);
58c4762a1bSJed Brown   }
59c4762a1bSJed Brown   len = (options->dim + 1) * PetscMax(1, options->numFields);
609566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(len, &options->numDof));
619566063dSJacob Faibussowitsch   PetscCall(PetscOptionsIntArray("-num_dof", "The dof signature for the section", "ex9.c", options->numDof, &len, &flg));
621dca8a05SBarry 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));
63c4762a1bSJed Brown 
64c4762a1bSJed Brown   /* We are specifying the scalar dof, so augment it for multiple components */
65c4762a1bSJed Brown   {
66c4762a1bSJed Brown     PetscInt f, d;
67c4762a1bSJed Brown 
68c4762a1bSJed Brown     for (f = 0; f < options->numFields; ++f) {
69c4762a1bSJed Brown       for (d = 0; d <= options->dim; ++d) options->numDof[f * (options->dim + 1) + d] *= options->numComponents[f];
70c4762a1bSJed Brown     }
71c4762a1bSJed Brown   }
72c4762a1bSJed Brown 
739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-reuse_array", "Pass in user allocated array to VecGetClosure()", "ex9.c", options->reuseArray, &options->reuseArray, NULL));
749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-errors", "Treat failures as errors", "ex9.c", options->errors, &options->errors, NULL));
759566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-iterations", "The number of iterations for a query", "ex9.c", options->iterations, &options->iterations, NULL, 0));
769566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-max_cone_time", "The maximum time per run for DMPlexGetCone()", "ex9.c", options->maxConeTime, &options->maxConeTime, NULL));
779566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-max_closure_time", "The maximum time per run for DMPlexGetTransitiveClosure()", "ex9.c", options->maxClosureTime, &options->maxClosureTime, NULL));
789566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-max_vec_closure_time", "The maximum time per run for DMPlexVecGetClosure()", "ex9.c", options->maxVecClosureTime, &options->maxVecClosureTime, NULL));
799566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-print_times", "Print total times, do not check limits", "ex9.c", options->printTimes, &options->printTimes, NULL));
80d0609cedSBarry Smith   PetscOptionsEnd();
813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
82c4762a1bSJed Brown }
83c4762a1bSJed Brown 
CreateSimplex_2D(MPI_Comm comm,DM * newdm)84d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateSimplex_2D(MPI_Comm comm, DM *newdm)
85d71ae5a4SJacob Faibussowitsch {
86c4762a1bSJed Brown   DM          dm;
87c4762a1bSJed Brown   PetscInt    numPoints[2]        = {4, 2};
88c4762a1bSJed Brown   PetscInt    coneSize[6]         = {3, 3, 0, 0, 0, 0};
89c4762a1bSJed Brown   PetscInt    cones[6]            = {2, 3, 4, 5, 4, 3};
90c4762a1bSJed Brown   PetscInt    coneOrientations[6] = {0, 0, 0, 0, 0, 0};
91c4762a1bSJed Brown   PetscScalar vertexCoords[8]     = {-0.5, 0.5, 0.0, 0.0, 0.0, 1.0, 0.5, 0.5};
92c4762a1bSJed Brown   PetscInt    markerPoints[8]     = {2, 1, 3, 1, 4, 1, 5, 1};
93c4762a1bSJed Brown   PetscInt    dim = 2, depth = 1, p;
94c4762a1bSJed Brown 
95c4762a1bSJed Brown   PetscFunctionBegin;
969566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, &dm));
979566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)dm, "triangular"));
989566063dSJacob Faibussowitsch   PetscCall(DMSetType(dm, DMPLEX));
999566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(dm, dim));
1009566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
10148a46eb9SPierre Jolivet   for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1]));
102c4762a1bSJed Brown   *newdm = dm;
1033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
104c4762a1bSJed Brown }
105c4762a1bSJed Brown 
CreateSimplex_3D(MPI_Comm comm,DM * newdm)106d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateSimplex_3D(MPI_Comm comm, DM *newdm)
107d71ae5a4SJacob Faibussowitsch {
108c4762a1bSJed Brown   DM          dm;
109c4762a1bSJed Brown   PetscInt    numPoints[2]        = {5, 2};
110c4762a1bSJed Brown   PetscInt    coneSize[23]        = {4, 4, 0, 0, 0, 0, 0};
111c4762a1bSJed Brown   PetscInt    cones[8]            = {2, 4, 3, 5, 3, 4, 6, 5};
112c4762a1bSJed Brown   PetscInt    coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
113c4762a1bSJed 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};
114c4762a1bSJed Brown   PetscInt    markerPoints[10]    = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1};
115c4762a1bSJed Brown   PetscInt    dim = 3, depth = 1, p;
116c4762a1bSJed Brown 
117c4762a1bSJed Brown   PetscFunctionBegin;
1189566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, &dm));
1199566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)dm, "tetrahedral"));
1209566063dSJacob Faibussowitsch   PetscCall(DMSetType(dm, DMPLEX));
1219566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(dm, dim));
1229566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
12348a46eb9SPierre Jolivet   for (p = 0; p < 5; ++p) PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1]));
124c4762a1bSJed Brown   *newdm = dm;
1253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
126c4762a1bSJed Brown }
127c4762a1bSJed Brown 
CreateQuad_2D(MPI_Comm comm,DM * newdm)128d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateQuad_2D(MPI_Comm comm, DM *newdm)
129d71ae5a4SJacob Faibussowitsch {
130c4762a1bSJed Brown   DM          dm;
131c4762a1bSJed Brown   PetscInt    numPoints[2]        = {6, 2};
132c4762a1bSJed Brown   PetscInt    coneSize[8]         = {4, 4, 0, 0, 0, 0, 0, 0};
133c4762a1bSJed Brown   PetscInt    cones[8]            = {2, 3, 4, 5, 3, 6, 7, 4};
134c4762a1bSJed Brown   PetscInt    coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
135c4762a1bSJed 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};
136c4762a1bSJed Brown   PetscInt    markerPoints[12]    = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1, 7, 1};
137c4762a1bSJed Brown   PetscInt    dim = 2, depth = 1, p;
138c4762a1bSJed Brown 
139c4762a1bSJed Brown   PetscFunctionBegin;
1409566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, &dm));
1419566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)dm, "quadrilateral"));
1429566063dSJacob Faibussowitsch   PetscCall(DMSetType(dm, DMPLEX));
1439566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(dm, dim));
1449566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
14548a46eb9SPierre Jolivet   for (p = 0; p < 6; ++p) PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1]));
146c4762a1bSJed Brown   *newdm = dm;
1473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
148c4762a1bSJed Brown }
149c4762a1bSJed Brown 
CreateHex_3D(MPI_Comm comm,DM * newdm)150d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateHex_3D(MPI_Comm comm, DM *newdm)
151d71ae5a4SJacob Faibussowitsch {
152c4762a1bSJed Brown   DM          dm;
153c4762a1bSJed Brown   PetscInt    numPoints[2]         = {12, 2};
154c4762a1bSJed Brown   PetscInt    coneSize[14]         = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
155c4762a1bSJed Brown   PetscInt    cones[16]            = {2, 5, 4, 3, 6, 7, 8, 9, 3, 4, 11, 10, 7, 12, 13, 8};
156c4762a1bSJed Brown   PetscInt    coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
1579371c9d4SSatish 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};
158c4762a1bSJed 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};
159c4762a1bSJed Brown   PetscInt    dim = 3, depth = 1, p;
160c4762a1bSJed Brown 
161c4762a1bSJed Brown   PetscFunctionBegin;
1629566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, &dm));
1639566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)dm, "hexahedral"));
1649566063dSJacob Faibussowitsch   PetscCall(DMSetType(dm, DMPLEX));
1659566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(dm, dim));
1669566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
16748a46eb9SPierre Jolivet   for (p = 0; p < 12; ++p) PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1]));
168c4762a1bSJed Brown   *newdm = dm;
1693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
170c4762a1bSJed Brown }
171c4762a1bSJed Brown 
CreateMesh(MPI_Comm comm,AppCtx * user,DM * newdm)172d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *newdm)
173d71ae5a4SJacob Faibussowitsch {
174c4762a1bSJed Brown   PetscInt  dim         = user->dim;
175c4762a1bSJed Brown   PetscBool cellSimplex = user->cellSimplex;
176c4762a1bSJed Brown 
177c4762a1bSJed Brown   PetscFunctionBegin;
178c4762a1bSJed Brown   switch (dim) {
179c4762a1bSJed Brown   case 2:
180c4762a1bSJed Brown     if (cellSimplex) {
1819566063dSJacob Faibussowitsch       PetscCall(CreateSimplex_2D(comm, newdm));
182c4762a1bSJed Brown     } else {
1839566063dSJacob Faibussowitsch       PetscCall(CreateQuad_2D(comm, newdm));
184c4762a1bSJed Brown     }
185c4762a1bSJed Brown     break;
186c4762a1bSJed Brown   case 3:
187c4762a1bSJed Brown     if (cellSimplex) {
1889566063dSJacob Faibussowitsch       PetscCall(CreateSimplex_3D(comm, newdm));
189c4762a1bSJed Brown     } else {
1909566063dSJacob Faibussowitsch       PetscCall(CreateHex_3D(comm, newdm));
191c4762a1bSJed Brown     }
192c4762a1bSJed Brown     break;
193d71ae5a4SJacob Faibussowitsch   default:
194d71ae5a4SJacob Faibussowitsch     SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %" PetscInt_FMT, dim);
195c4762a1bSJed Brown   }
196c4762a1bSJed Brown   if (user->refinementLimit > 0.0) {
197c4762a1bSJed Brown     DM          rdm;
198c4762a1bSJed Brown     const char *name;
199c4762a1bSJed Brown 
2009566063dSJacob Faibussowitsch     PetscCall(DMPlexSetRefinementUniform(*newdm, PETSC_FALSE));
2019566063dSJacob Faibussowitsch     PetscCall(DMPlexSetRefinementLimit(*newdm, user->refinementLimit));
2029566063dSJacob Faibussowitsch     PetscCall(DMRefine(*newdm, PETSC_COMM_SELF, &rdm));
2039566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName((PetscObject)*newdm, &name));
2049566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)rdm, name));
2059566063dSJacob Faibussowitsch     PetscCall(DMDestroy(newdm));
206c4762a1bSJed Brown     *newdm = rdm;
207c4762a1bSJed Brown   }
208c4762a1bSJed Brown   if (user->interpolate) {
209c4762a1bSJed Brown     DM idm;
210c4762a1bSJed Brown 
2119566063dSJacob Faibussowitsch     PetscCall(DMPlexInterpolate(*newdm, &idm));
2129566063dSJacob Faibussowitsch     PetscCall(DMDestroy(newdm));
213c4762a1bSJed Brown     *newdm = idm;
214c4762a1bSJed Brown   }
2159566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*newdm));
2163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
217c4762a1bSJed Brown }
218c4762a1bSJed Brown 
TestCone(DM dm,AppCtx * user)219d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestCone(DM dm, AppCtx *user)
220d71ae5a4SJacob Faibussowitsch {
221c4762a1bSJed Brown   PetscInt           numRuns, cStart, cEnd, c, i;
222c4762a1bSJed Brown   PetscReal          maxTimePerRun = user->maxConeTime;
223c4762a1bSJed Brown   PetscLogStage      stage;
224c4762a1bSJed Brown   PetscLogEvent      event;
225c4762a1bSJed Brown   PetscEventPerfInfo eventInfo;
226c4762a1bSJed Brown   MPI_Comm           comm;
227c4762a1bSJed Brown   PetscMPIInt        rank;
228c4762a1bSJed Brown 
229c4762a1bSJed Brown   PetscFunctionBegin;
2309566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
2319566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
2329566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("DMPlex Cone Test", &stage));
2339566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("Cone", PETSC_OBJECT_CLASSID, &event));
2349566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(stage));
2359566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
2369566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(event, 0, 0, 0, 0));
237c4762a1bSJed Brown   for (i = 0; i < user->iterations; ++i) {
238c4762a1bSJed Brown     for (c = cStart; c < cEnd; ++c) {
239c4762a1bSJed Brown       const PetscInt *cone;
240c4762a1bSJed Brown 
2419566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCone(dm, c, &cone));
242c4762a1bSJed Brown     }
243c4762a1bSJed Brown   }
2449566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(event, 0, 0, 0, 0));
2459566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePop());
246c4762a1bSJed Brown 
2479566063dSJacob Faibussowitsch   PetscCall(PetscLogEventGetPerfInfo(stage, event, &eventInfo));
248c4762a1bSJed Brown   numRuns = (cEnd - cStart) * user->iterations;
24963a3b9bcSJacob Faibussowitsch   PetscCheck(eventInfo.count == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event calls %d should be 1", eventInfo.count);
25063a3b9bcSJacob Faibussowitsch   PetscCheck((PetscInt)eventInfo.flops == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event flops %" PetscInt_FMT " should be 0", (PetscInt)eventInfo.flops);
251c4762a1bSJed Brown   if (user->printTimes) {
25263a3b9bcSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(comm, "[%d] Cones: %" PetscInt_FMT " Total time: %.3es Average time per cone: %.3es\n", rank, numRuns, eventInfo.time, eventInfo.time / numRuns));
2539566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT));
254c4762a1bSJed Brown   } else if (eventInfo.time > maxTimePerRun * numRuns) {
255*300f1712SStefano Zampini     PetscCall(PetscSynchronizedPrintf(comm, "[%d] Cones: %" PetscInt_FMT " Average time per cone: %gs standard: %gs\n", rank, numRuns, eventInfo.time / numRuns, (double)maxTimePerRun));
2569566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT));
257*300f1712SStefano Zampini     PetscCheck(!user->errors, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Average time for cone %g > standard %g", eventInfo.time / numRuns, (double)maxTimePerRun);
258c4762a1bSJed Brown   }
2593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
260c4762a1bSJed Brown }
261c4762a1bSJed Brown 
TestTransitiveClosure(DM dm,AppCtx * user)262d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestTransitiveClosure(DM dm, AppCtx *user)
263d71ae5a4SJacob Faibussowitsch {
264c4762a1bSJed Brown   PetscInt           numRuns, cStart, cEnd, c, i;
265c4762a1bSJed Brown   PetscReal          maxTimePerRun = user->maxClosureTime;
266c4762a1bSJed Brown   PetscLogStage      stage;
267c4762a1bSJed Brown   PetscLogEvent      event;
268c4762a1bSJed Brown   PetscEventPerfInfo eventInfo;
269c4762a1bSJed Brown   MPI_Comm           comm;
270c4762a1bSJed Brown   PetscMPIInt        rank;
271c4762a1bSJed Brown 
272c4762a1bSJed Brown   PetscFunctionBegin;
2739566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
2749566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
2759566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("DMPlex Transitive Closure Test", &stage));
2769566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("TransitiveClosure", PETSC_OBJECT_CLASSID, &event));
2779566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(stage));
2789566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
2799566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(event, 0, 0, 0, 0));
280c4762a1bSJed Brown   for (i = 0; i < user->iterations; ++i) {
281c4762a1bSJed Brown     for (c = cStart; c < cEnd; ++c) {
282c4762a1bSJed Brown       PetscInt *closure = NULL;
283c4762a1bSJed Brown       PetscInt  closureSize;
284c4762a1bSJed Brown 
2859566063dSJacob Faibussowitsch       PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
2869566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
287c4762a1bSJed Brown     }
288c4762a1bSJed Brown   }
2899566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(event, 0, 0, 0, 0));
2909566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePop());
291c4762a1bSJed Brown 
2929566063dSJacob Faibussowitsch   PetscCall(PetscLogEventGetPerfInfo(stage, event, &eventInfo));
293c4762a1bSJed Brown   numRuns = (cEnd - cStart) * user->iterations;
29463a3b9bcSJacob Faibussowitsch   PetscCheck(eventInfo.count == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event calls %d should be 1", eventInfo.count);
29563a3b9bcSJacob Faibussowitsch   PetscCheck((PetscInt)eventInfo.flops == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event flops %" PetscInt_FMT " should be 0", (PetscInt)eventInfo.flops);
296c4762a1bSJed Brown   if (user->printTimes) {
29763a3b9bcSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(comm, "[%d] Closures: %" PetscInt_FMT " Total time: %.3es Average time per cone: %.3es\n", rank, numRuns, eventInfo.time, eventInfo.time / numRuns));
2989566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT));
299c4762a1bSJed Brown   } else if (eventInfo.time > maxTimePerRun * numRuns) {
300*300f1712SStefano Zampini     PetscCall(PetscSynchronizedPrintf(comm, "[%d] Closures: %" PetscInt_FMT " Average time per cone: %gs standard: %gs\n", rank, numRuns, eventInfo.time / numRuns, (double)maxTimePerRun));
3019566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT));
302*300f1712SStefano Zampini     PetscCheck(!user->errors, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Average time for closure %g > standard %g", eventInfo.time / numRuns, (double)maxTimePerRun);
303c4762a1bSJed Brown   }
3043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
305c4762a1bSJed Brown }
306c4762a1bSJed Brown 
TestVecClosure(DM dm,PetscBool useIndex,PetscBool useSpectral,AppCtx * user)307d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestVecClosure(DM dm, PetscBool useIndex, PetscBool useSpectral, AppCtx *user)
308d71ae5a4SJacob Faibussowitsch {
309c4762a1bSJed Brown   PetscSection       s;
310c4762a1bSJed Brown   Vec                v;
311c4762a1bSJed Brown   PetscInt           numRuns, cStart, cEnd, c, i;
312c4762a1bSJed Brown   PetscScalar        tmpArray[64];
313c4762a1bSJed Brown   PetscScalar       *userArray     = user->reuseArray ? tmpArray : NULL;
314c4762a1bSJed Brown   PetscReal          maxTimePerRun = user->maxVecClosureTime;
315c4762a1bSJed Brown   PetscLogStage      stage;
316c4762a1bSJed Brown   PetscLogEvent      event;
317c4762a1bSJed Brown   PetscEventPerfInfo eventInfo;
318c4762a1bSJed Brown   MPI_Comm           comm;
319c4762a1bSJed Brown   PetscMPIInt        rank;
320c4762a1bSJed Brown 
321c4762a1bSJed Brown   PetscFunctionBegin;
3229566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
3239566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
324c4762a1bSJed Brown   if (useIndex) {
325c4762a1bSJed Brown     if (useSpectral) {
3269566063dSJacob Faibussowitsch       PetscCall(PetscLogStageRegister("DMPlex Vector Closure with Index Test", &stage));
3279566063dSJacob Faibussowitsch       PetscCall(PetscLogEventRegister("VecClosureInd", PETSC_OBJECT_CLASSID, &event));
328c4762a1bSJed Brown     } else {
3299566063dSJacob Faibussowitsch       PetscCall(PetscLogStageRegister("DMPlex Vector Spectral Closure with Index Test", &stage));
3309566063dSJacob Faibussowitsch       PetscCall(PetscLogEventRegister("VecClosureSpecInd", PETSC_OBJECT_CLASSID, &event));
331c4762a1bSJed Brown     }
332c4762a1bSJed Brown   } else {
333c4762a1bSJed Brown     if (useSpectral) {
3349566063dSJacob Faibussowitsch       PetscCall(PetscLogStageRegister("DMPlex Vector Spectral Closure Test", &stage));
3359566063dSJacob Faibussowitsch       PetscCall(PetscLogEventRegister("VecClosureSpec", PETSC_OBJECT_CLASSID, &event));
336c4762a1bSJed Brown     } else {
3379566063dSJacob Faibussowitsch       PetscCall(PetscLogStageRegister("DMPlex Vector Closure Test", &stage));
3389566063dSJacob Faibussowitsch       PetscCall(PetscLogEventRegister("VecClosure", PETSC_OBJECT_CLASSID, &event));
339c4762a1bSJed Brown     }
340c4762a1bSJed Brown   }
3419566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(stage));
3429566063dSJacob Faibussowitsch   PetscCall(DMSetNumFields(dm, user->numFields));
3439566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateSection(dm, NULL, user->numComponents, user->numDof, 0, NULL, NULL, NULL, NULL, &s));
3449566063dSJacob Faibussowitsch   PetscCall(DMSetLocalSection(dm, s));
3459566063dSJacob Faibussowitsch   if (useIndex) PetscCall(DMPlexCreateClosureIndex(dm, s));
3469566063dSJacob Faibussowitsch   if (useSpectral) PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, s));
3479566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&s));
3489566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
3499566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &v));
3509566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(event, 0, 0, 0, 0));
351c4762a1bSJed Brown   for (i = 0; i < user->iterations; ++i) {
352c4762a1bSJed Brown     for (c = cStart; c < cEnd; ++c) {
353c4762a1bSJed Brown       PetscScalar *closure     = userArray;
354c4762a1bSJed Brown       PetscInt     closureSize = 64;
355c4762a1bSJed Brown 
3569566063dSJacob Faibussowitsch       PetscCall(DMPlexVecGetClosure(dm, s, v, c, &closureSize, &closure));
3579566063dSJacob Faibussowitsch       if (!user->reuseArray) PetscCall(DMPlexVecRestoreClosure(dm, s, v, c, &closureSize, &closure));
358c4762a1bSJed Brown     }
359c4762a1bSJed Brown   }
3609566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(event, 0, 0, 0, 0));
3619566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &v));
3629566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePop());
363c4762a1bSJed Brown 
3649566063dSJacob Faibussowitsch   PetscCall(PetscLogEventGetPerfInfo(stage, event, &eventInfo));
365c4762a1bSJed Brown   numRuns = (cEnd - cStart) * user->iterations;
36663a3b9bcSJacob Faibussowitsch   PetscCheck(eventInfo.count == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event calls %d should be 1", eventInfo.count);
36763a3b9bcSJacob Faibussowitsch   PetscCheck((PetscInt)eventInfo.flops == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event flops %" PetscInt_FMT " should be 0", (PetscInt)eventInfo.flops);
368c4762a1bSJed Brown   if (user->printTimes || eventInfo.time > maxTimePerRun * numRuns) {
369c4762a1bSJed Brown     const char *title          = "VecClosures";
370c4762a1bSJed Brown     const char *titleIndex     = "VecClosures with Index";
371c4762a1bSJed Brown     const char *titleSpec      = "VecClosures Spectral";
372c4762a1bSJed Brown     const char *titleSpecIndex = "VecClosures Spectral with Index";
373c4762a1bSJed Brown 
374c4762a1bSJed Brown     if (user->printTimes) {
3759371c9d4SSatish 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,
3769371c9d4SSatish Balay                                         eventInfo.time, eventInfo.time / numRuns));
3779566063dSJacob Faibussowitsch       PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT));
378c4762a1bSJed Brown     } else {
379*300f1712SStefano Zampini       PetscCall(PetscSynchronizedPrintf(comm, "[%d] %s: %" PetscInt_FMT " Average time per vector closure: %gs standard: %gs\n", rank, useIndex ? (useSpectral ? titleSpecIndex : titleIndex) : (useSpectral ? titleSpec : title), numRuns, eventInfo.time / numRuns, (double)maxTimePerRun));
3809566063dSJacob Faibussowitsch       PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT));
381*300f1712SStefano Zampini       PetscCheck(!user->errors, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Average time for vector closure %g > standard %g", eventInfo.time / numRuns, (double)maxTimePerRun);
382c4762a1bSJed Brown     }
383c4762a1bSJed Brown   }
3843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
385c4762a1bSJed Brown }
386c4762a1bSJed Brown 
CleanupContext(AppCtx * user)387d71ae5a4SJacob Faibussowitsch static PetscErrorCode CleanupContext(AppCtx *user)
388d71ae5a4SJacob Faibussowitsch {
389c4762a1bSJed Brown   PetscFunctionBegin;
3909566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->numComponents));
3919566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->numDof));
3923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
393c4762a1bSJed Brown }
394c4762a1bSJed Brown 
main(int argc,char ** argv)395d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
396d71ae5a4SJacob Faibussowitsch {
397c4762a1bSJed Brown   DM     dm;
398c4762a1bSJed Brown   AppCtx user;
399c4762a1bSJed Brown 
400327415f7SBarry Smith   PetscFunctionBeginUser;
4019566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
4029566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(&user));
4039566063dSJacob Faibussowitsch   PetscCall(PetscLogDefaultBegin());
4049566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
4059566063dSJacob Faibussowitsch   PetscCall(TestCone(dm, &user));
4069566063dSJacob Faibussowitsch   PetscCall(TestTransitiveClosure(dm, &user));
4079566063dSJacob Faibussowitsch   PetscCall(TestVecClosure(dm, PETSC_FALSE, PETSC_FALSE, &user));
4089566063dSJacob Faibussowitsch   PetscCall(TestVecClosure(dm, PETSC_TRUE, PETSC_FALSE, &user));
409c4762a1bSJed Brown   if (!user.cellSimplex && user.spectral) {
4109566063dSJacob Faibussowitsch     PetscCall(TestVecClosure(dm, PETSC_FALSE, PETSC_TRUE, &user));
4119566063dSJacob Faibussowitsch     PetscCall(TestVecClosure(dm, PETSC_TRUE, PETSC_TRUE, &user));
412c4762a1bSJed Brown   }
4139566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
4149566063dSJacob Faibussowitsch   PetscCall(CleanupContext(&user));
4159566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
416b122ec5aSJacob Faibussowitsch   return 0;
417c4762a1bSJed Brown }
418c4762a1bSJed Brown 
419c4762a1bSJed Brown /*TEST
420c4762a1bSJed Brown 
421c4762a1bSJed Brown   build:
422dfd57a17SPierre Jolivet     requires: defined(PETSC_USE_LOG)
423c4762a1bSJed Brown 
424c4762a1bSJed Brown   # 2D Simplex P_1 scalar tests
425c4762a1bSJed Brown   testset:
426c4762a1bSJed Brown     args: -num_dof 1,0,0 -iterations 2 -print_times
427c4762a1bSJed Brown     test:
428c4762a1bSJed Brown       suffix: correctness_0
429c4762a1bSJed Brown     test:
430c4762a1bSJed Brown       suffix: correctness_1
431c4762a1bSJed Brown       args: -interpolate -dm_refine 2
432c4762a1bSJed Brown     test:
433c4762a1bSJed Brown       suffix: correctness_2
434c4762a1bSJed Brown       requires: triangle
43556245cf9SMatthew G. Knepley       args: -interpolate -dm_refine 5
436c4762a1bSJed Brown   test:
437c4762a1bSJed Brown     suffix: 0
438c4762a1bSJed Brown     TODO: Only for performance testing
439c4762a1bSJed 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
440c4762a1bSJed Brown   test:
441c4762a1bSJed Brown     suffix: 1
442c4762a1bSJed Brown     requires: triangle
443c4762a1bSJed Brown     TODO: Only for performance testing
444c4762a1bSJed 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
445c4762a1bSJed Brown   test:
446c4762a1bSJed Brown     suffix: 2
447c4762a1bSJed Brown     TODO: Only for performance testing
448c4762a1bSJed 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
449c4762a1bSJed Brown   test:
450c4762a1bSJed Brown     suffix: 3
451c4762a1bSJed Brown     requires: triangle
452c4762a1bSJed Brown     TODO: Only for performance testing
453c4762a1bSJed 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
454c4762a1bSJed Brown   test:
455c4762a1bSJed Brown     suffix: 4
456c4762a1bSJed Brown     TODO: Only for performance testing
457c4762a1bSJed 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
458c4762a1bSJed Brown   test:
459c4762a1bSJed Brown     suffix: 5
460c4762a1bSJed Brown     requires: triangle
461c4762a1bSJed Brown     TODO: Only for performance testing
462c4762a1bSJed 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
463c4762a1bSJed Brown   test:
464c4762a1bSJed Brown     suffix: 6
465c4762a1bSJed Brown     TODO: Only for performance testing
466c4762a1bSJed 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
467c4762a1bSJed Brown   test:
468c4762a1bSJed Brown     suffix: 7
469c4762a1bSJed Brown     requires: triangle
470c4762a1bSJed Brown     TODO: Only for performance testing
471c4762a1bSJed 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
472c4762a1bSJed Brown 
473c4762a1bSJed Brown   # 2D Simplex P_1 vector tests
474c4762a1bSJed Brown   # 2D Simplex P_2 scalar tests
475c4762a1bSJed Brown   # 2D Simplex P_2 vector tests
476c4762a1bSJed Brown   # 2D Simplex P_2/P_1 vector/scalar tests
477c4762a1bSJed Brown   # 2D Quad P_1 scalar tests
478c4762a1bSJed Brown   # 2D Quad P_1 vector tests
479c4762a1bSJed Brown   # 2D Quad P_2 scalar tests
480c4762a1bSJed Brown   # 2D Quad P_2 vector tests
481c4762a1bSJed Brown   # 3D Simplex P_1 scalar tests
482c4762a1bSJed Brown   # 3D Simplex P_1 vector tests
483c4762a1bSJed Brown   # 3D Simplex P_2 scalar tests
484c4762a1bSJed Brown   # 3D Simplex P_2 vector tests
485c4762a1bSJed Brown   # 3D Hex P_1 scalar tests
486c4762a1bSJed Brown   # 3D Hex P_1 vector tests
487c4762a1bSJed Brown   # 3D Hex P_2 scalar tests
488c4762a1bSJed Brown   # 3D Hex P_2 vector tests
489c4762a1bSJed Brown 
490c4762a1bSJed Brown TEST*/
491