xref: /petsc/src/dm/impls/plex/tests/ex9.c (revision dfd57a172ac9fa6c7b5fe6de6ab5df85cefc2996)
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 
24c4762a1bSJed Brown static PetscErrorCode ProcessOptions(AppCtx *options)
25c4762a1bSJed Brown {
26c4762a1bSJed Brown   PetscInt       len;
27c4762a1bSJed Brown   PetscBool      flg;
28c4762a1bSJed Brown   PetscErrorCode ierr;
29c4762a1bSJed Brown 
30c4762a1bSJed Brown   PetscFunctionBegin;
31c4762a1bSJed Brown   options->dim               = 2;
32c4762a1bSJed Brown   options->cellSimplex       = PETSC_TRUE;
33c4762a1bSJed Brown   options->spectral          = PETSC_FALSE;
34c4762a1bSJed Brown   options->interpolate       = PETSC_FALSE;
35c4762a1bSJed Brown   options->refinementLimit   = 0.0;
36c4762a1bSJed Brown   options->numFields         = 0;
37c4762a1bSJed Brown   options->numComponents     = NULL;
38c4762a1bSJed Brown   options->numDof            = NULL;
39c4762a1bSJed Brown   options->reuseArray        = PETSC_FALSE;
40c4762a1bSJed Brown   options->errors            = PETSC_FALSE;
41c4762a1bSJed Brown   options->iterations        = 1;
42c4762a1bSJed Brown   options->maxConeTime       = 0.0;
43c4762a1bSJed Brown   options->maxClosureTime    = 0.0;
44c4762a1bSJed Brown   options->maxVecClosureTime = 0.0;
45c4762a1bSJed Brown   options->printTimes        = PETSC_FALSE;
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   ierr = PetscOptionsBegin(PETSC_COMM_SELF, "", "Meshing Problem Options", "DMPLEX");CHKERRQ(ierr);
48c4762a1bSJed Brown   ierr = PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex9.c", options->dim, &options->dim, NULL,1,3);CHKERRQ(ierr);
49c4762a1bSJed Brown   ierr = PetscOptionsBool("-cellSimplex", "Flag for simplices", "ex9.c", options->cellSimplex, &options->cellSimplex, NULL);CHKERRQ(ierr);
50c4762a1bSJed Brown   ierr = PetscOptionsBool("-spectral", "Flag for spectral element layout", "ex9.c", options->spectral, &options->spectral, NULL);CHKERRQ(ierr);
51c4762a1bSJed Brown   ierr = PetscOptionsBool("-interpolate", "Flag for mesh interpolation", "ex9.c", options->interpolate, &options->interpolate, NULL);CHKERRQ(ierr);
52c4762a1bSJed Brown   ierr = PetscOptionsReal("-refinement_limit", "The maximum volume of a refined cell", "ex9.c", options->refinementLimit, &options->refinementLimit, NULL);CHKERRQ(ierr);
53c4762a1bSJed Brown   ierr = PetscOptionsBoundedInt("-num_fields", "The number of section fields", "ex9.c", options->numFields, &options->numFields, NULL, 0);CHKERRQ(ierr);
54c4762a1bSJed Brown   if (options->numFields) {
55c4762a1bSJed Brown     len  = options->numFields;
56c4762a1bSJed Brown     ierr = PetscMalloc1(len, &options->numComponents);CHKERRQ(ierr);
57c4762a1bSJed Brown     ierr = PetscOptionsIntArray("-num_components", "The number of components per field", "ex9.c", options->numComponents, &len, &flg);CHKERRQ(ierr);
58c4762a1bSJed Brown     if (flg && (len != options->numFields)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Length of components array is %d should be %d", len, options->numFields);
59c4762a1bSJed Brown   }
60c4762a1bSJed Brown   len  = (options->dim+1) * PetscMax(1, options->numFields);
61c4762a1bSJed Brown   ierr = PetscMalloc1(len, &options->numDof);CHKERRQ(ierr);
62c4762a1bSJed Brown   ierr = PetscOptionsIntArray("-num_dof", "The dof signature for the section", "ex9.c", options->numDof, &len, &flg);CHKERRQ(ierr);
63c4762a1bSJed Brown   if (flg && (len != (options->dim+1) * PetscMax(1, options->numFields))) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Length of dof array is %d should be %d", len, (options->dim+1) * PetscMax(1, options->numFields));
64c4762a1bSJed Brown 
65c4762a1bSJed Brown   /* We are specifying the scalar dof, so augment it for multiple components */
66c4762a1bSJed Brown   {
67c4762a1bSJed Brown     PetscInt f, d;
68c4762a1bSJed Brown 
69c4762a1bSJed Brown     for (f = 0; f < options->numFields; ++f) {
70c4762a1bSJed Brown       for (d = 0; d <= options->dim; ++d) options->numDof[f*(options->dim+1)+d] *= options->numComponents[f];
71c4762a1bSJed Brown     }
72c4762a1bSJed Brown   }
73c4762a1bSJed Brown 
74c4762a1bSJed Brown   ierr = PetscOptionsBool("-reuse_array", "Pass in user allocated array to VecGetClosure()", "ex9.c", options->reuseArray, &options->reuseArray, NULL);CHKERRQ(ierr);
75c4762a1bSJed Brown   ierr = PetscOptionsBool("-errors", "Treat failures as errors", "ex9.c", options->errors, &options->errors, NULL);CHKERRQ(ierr);
76c4762a1bSJed Brown   ierr = PetscOptionsBoundedInt("-iterations", "The number of iterations for a query", "ex9.c", options->iterations, &options->iterations, NULL,0);CHKERRQ(ierr);
77c4762a1bSJed Brown   ierr = PetscOptionsReal("-max_cone_time", "The maximum time per run for DMPlexGetCone()", "ex9.c", options->maxConeTime, &options->maxConeTime, NULL);CHKERRQ(ierr);
78c4762a1bSJed Brown   ierr = PetscOptionsReal("-max_closure_time", "The maximum time per run for DMPlexGetTransitiveClosure()", "ex9.c", options->maxClosureTime, &options->maxClosureTime, NULL);CHKERRQ(ierr);
79c4762a1bSJed Brown   ierr = PetscOptionsReal("-max_vec_closure_time", "The maximum time per run for DMPlexVecGetClosure()", "ex9.c", options->maxVecClosureTime, &options->maxVecClosureTime, NULL);CHKERRQ(ierr);
80c4762a1bSJed Brown   ierr = PetscOptionsBool("-print_times", "Print total times, do not check limits", "ex9.c", options->printTimes, &options->printTimes, NULL);CHKERRQ(ierr);
81c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
82c4762a1bSJed Brown   PetscFunctionReturn(0);
83c4762a1bSJed Brown }
84c4762a1bSJed Brown 
85c4762a1bSJed Brown static PetscErrorCode CreateSimplex_2D(MPI_Comm comm, DM *newdm)
86c4762a1bSJed Brown {
87c4762a1bSJed Brown   DM             dm;
88c4762a1bSJed Brown   PetscInt       numPoints[2]        = {4, 2};
89c4762a1bSJed Brown   PetscInt       coneSize[6]         = {3, 3, 0, 0, 0, 0};
90c4762a1bSJed Brown   PetscInt       cones[6]            = {2, 3, 4,  5, 4, 3};
91c4762a1bSJed Brown   PetscInt       coneOrientations[6] = {0, 0, 0,  0, 0, 0};
92c4762a1bSJed Brown   PetscScalar    vertexCoords[8]     = {-0.5, 0.5, 0.0, 0.0, 0.0, 1.0, 0.5, 0.5};
93c4762a1bSJed Brown   PetscInt       markerPoints[8]     = {2, 1, 3, 1, 4, 1, 5, 1};
94c4762a1bSJed Brown   PetscInt       dim = 2, depth = 1, p;
95c4762a1bSJed Brown   PetscErrorCode ierr;
96c4762a1bSJed Brown 
97c4762a1bSJed Brown   PetscFunctionBegin;
98c4762a1bSJed Brown   ierr = DMCreate(comm, &dm);CHKERRQ(ierr);
99c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) dm, "triangular");CHKERRQ(ierr);
100c4762a1bSJed Brown   ierr = DMSetType(dm, DMPLEX);CHKERRQ(ierr);
101c4762a1bSJed Brown   ierr = DMSetDimension(dm, dim);CHKERRQ(ierr);
102c4762a1bSJed Brown   ierr = DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
103c4762a1bSJed Brown   for (p = 0; p < 4; ++p) {
104c4762a1bSJed Brown     ierr = DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);
105c4762a1bSJed Brown   }
106c4762a1bSJed Brown   *newdm = dm;
107c4762a1bSJed Brown   PetscFunctionReturn(0);
108c4762a1bSJed Brown }
109c4762a1bSJed Brown 
110c4762a1bSJed Brown static PetscErrorCode CreateSimplex_3D(MPI_Comm comm, DM *newdm)
111c4762a1bSJed Brown {
112c4762a1bSJed Brown   DM             dm;
113c4762a1bSJed Brown   PetscInt       numPoints[2]        = {5, 2};
114c4762a1bSJed Brown   PetscInt       coneSize[23]        = {4, 4, 0, 0, 0, 0, 0};
115c4762a1bSJed Brown   PetscInt       cones[8]            = {2, 4, 3, 5,  3, 4, 6, 5};
116c4762a1bSJed Brown   PetscInt       coneOrientations[8] = {0, 0, 0, 0,  0, 0, 0, 0};
117c4762a1bSJed 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};
118c4762a1bSJed Brown   PetscInt       markerPoints[10]    = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1};
119c4762a1bSJed Brown   PetscInt       dim = 3, depth = 1, p;
120c4762a1bSJed Brown   PetscErrorCode ierr;
121c4762a1bSJed Brown 
122c4762a1bSJed Brown   PetscFunctionBegin;
123c4762a1bSJed Brown   ierr = DMCreate(comm, &dm);CHKERRQ(ierr);
124c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) dm, "tetrahedral");CHKERRQ(ierr);
125c4762a1bSJed Brown   ierr = DMSetType(dm, DMPLEX);CHKERRQ(ierr);
126c4762a1bSJed Brown   ierr = DMSetDimension(dm, dim);CHKERRQ(ierr);
127c4762a1bSJed Brown   ierr = DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
128c4762a1bSJed Brown   for (p = 0; p < 5; ++p) {
129c4762a1bSJed Brown     ierr = DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);
130c4762a1bSJed Brown   }
131c4762a1bSJed Brown   *newdm = dm;
132c4762a1bSJed Brown   PetscFunctionReturn(0);
133c4762a1bSJed Brown }
134c4762a1bSJed Brown 
135c4762a1bSJed Brown static PetscErrorCode CreateQuad_2D(MPI_Comm comm, DM *newdm)
136c4762a1bSJed Brown {
137c4762a1bSJed Brown   DM             dm;
138c4762a1bSJed Brown   PetscInt       numPoints[2]        = {6, 2};
139c4762a1bSJed Brown   PetscInt       coneSize[8]         = {4, 4, 0, 0, 0, 0, 0, 0};
140c4762a1bSJed Brown   PetscInt       cones[8]            = {2, 3, 4, 5,  3, 6, 7, 4};
141c4762a1bSJed Brown   PetscInt       coneOrientations[8] = {0, 0, 0, 0,  0, 0, 0, 0};
142c4762a1bSJed 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};
143c4762a1bSJed Brown   PetscInt       markerPoints[12]    = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1, 7, 1};
144c4762a1bSJed Brown   PetscInt       dim = 2, depth = 1, p;
145c4762a1bSJed Brown   PetscErrorCode ierr;
146c4762a1bSJed Brown 
147c4762a1bSJed Brown   PetscFunctionBegin;
148c4762a1bSJed Brown   ierr = DMCreate(comm, &dm);CHKERRQ(ierr);
149c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) dm, "quadrilateral");CHKERRQ(ierr);
150c4762a1bSJed Brown   ierr = DMSetType(dm, DMPLEX);CHKERRQ(ierr);
151c4762a1bSJed Brown   ierr = DMSetDimension(dm, dim);CHKERRQ(ierr);
152c4762a1bSJed Brown   ierr = DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
153c4762a1bSJed Brown   for (p = 0; p < 6; ++p) {
154c4762a1bSJed Brown     ierr = DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);
155c4762a1bSJed Brown   }
156c4762a1bSJed Brown   *newdm = dm;
157c4762a1bSJed Brown   PetscFunctionReturn(0);
158c4762a1bSJed Brown }
159c4762a1bSJed Brown 
160c4762a1bSJed Brown static PetscErrorCode CreateHex_3D(MPI_Comm comm, DM *newdm)
161c4762a1bSJed Brown {
162c4762a1bSJed Brown   DM             dm;
163c4762a1bSJed Brown   PetscInt       numPoints[2]         = {12, 2};
164c4762a1bSJed Brown   PetscInt       coneSize[14]         = {8, 8, 0,0,0,0,0,0,0,0,0,0,0,0};
165c4762a1bSJed Brown   PetscInt       cones[16]            = {2,5,4,3,6,7,8,9,  3,4,11,10,7,12,13,8};
166c4762a1bSJed Brown   PetscInt       coneOrientations[16] = {0,0,0,0,0,0,0,0,  0,0, 0, 0,0, 0, 0,0};
167c4762a1bSJed Brown   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,
168c4762a1bSJed Brown                                          -0.5,0.0,1.0, 0.0,0.0,1.0, 0.0,1.0,1.0, -0.5,1.0,1.0,
169c4762a1bSJed Brown                                           0.5,0.0,0.0, 0.5,1.0,0.0, 0.5,0.0,1.0,  0.5,1.0,1.0};
170c4762a1bSJed 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};
171c4762a1bSJed Brown   PetscInt       dim = 3, depth = 1, p;
172c4762a1bSJed Brown   PetscErrorCode ierr;
173c4762a1bSJed Brown 
174c4762a1bSJed Brown   PetscFunctionBegin;
175c4762a1bSJed Brown   ierr = DMCreate(comm, &dm);CHKERRQ(ierr);
176c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) dm, "hexahedral");CHKERRQ(ierr);
177c4762a1bSJed Brown   ierr = DMSetType(dm, DMPLEX);CHKERRQ(ierr);
178c4762a1bSJed Brown   ierr = DMSetDimension(dm, dim);CHKERRQ(ierr);
179c4762a1bSJed Brown   ierr = DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
180c4762a1bSJed Brown   for (p = 0; p < 12; ++p) {
181c4762a1bSJed Brown     ierr = DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);
182c4762a1bSJed Brown   }
183c4762a1bSJed Brown   *newdm = dm;
184c4762a1bSJed Brown   PetscFunctionReturn(0);
185c4762a1bSJed Brown }
186c4762a1bSJed Brown 
187c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *newdm)
188c4762a1bSJed Brown {
189c4762a1bSJed Brown   PetscInt       dim         = user->dim;
190c4762a1bSJed Brown   PetscBool      cellSimplex = user->cellSimplex;
191c4762a1bSJed Brown   PetscErrorCode ierr;
192c4762a1bSJed Brown 
193c4762a1bSJed Brown   PetscFunctionBegin;
194c4762a1bSJed Brown   switch (dim) {
195c4762a1bSJed Brown   case 2:
196c4762a1bSJed Brown     if (cellSimplex) {
197c4762a1bSJed Brown       ierr = CreateSimplex_2D(comm, newdm);CHKERRQ(ierr);
198c4762a1bSJed Brown     } else {
199c4762a1bSJed Brown       ierr = CreateQuad_2D(comm, newdm);CHKERRQ(ierr);
200c4762a1bSJed Brown     }
201c4762a1bSJed Brown     break;
202c4762a1bSJed Brown   case 3:
203c4762a1bSJed Brown     if (cellSimplex) {
204c4762a1bSJed Brown       ierr = CreateSimplex_3D(comm, newdm);CHKERRQ(ierr);
205c4762a1bSJed Brown     } else {
206c4762a1bSJed Brown       ierr = CreateHex_3D(comm, newdm);CHKERRQ(ierr);
207c4762a1bSJed Brown     }
208c4762a1bSJed Brown     break;
209c4762a1bSJed Brown   default:
210c4762a1bSJed Brown     SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %d", dim);
211c4762a1bSJed Brown   }
212c4762a1bSJed Brown   if (user->refinementLimit > 0.0) {
213c4762a1bSJed Brown     DM rdm;
214c4762a1bSJed Brown     const char *name;
215c4762a1bSJed Brown 
216c4762a1bSJed Brown     ierr = DMPlexSetRefinementUniform(*newdm, PETSC_FALSE);CHKERRQ(ierr);
217c4762a1bSJed Brown     ierr = DMPlexSetRefinementLimit(*newdm, user->refinementLimit);CHKERRQ(ierr);
218c4762a1bSJed Brown     ierr = DMRefine(*newdm, PETSC_COMM_SELF, &rdm);CHKERRQ(ierr);
219c4762a1bSJed Brown     ierr = PetscObjectGetName((PetscObject) *newdm, &name);CHKERRQ(ierr);
220c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject)    rdm,  name);CHKERRQ(ierr);
221c4762a1bSJed Brown     ierr = DMDestroy(newdm);CHKERRQ(ierr);
222c4762a1bSJed Brown     *newdm = rdm;
223c4762a1bSJed Brown   }
224c4762a1bSJed Brown   if (user->interpolate) {
225c4762a1bSJed Brown     DM idm;
226c4762a1bSJed Brown 
227c4762a1bSJed Brown     ierr = DMPlexInterpolate(*newdm, &idm);CHKERRQ(ierr);
228c4762a1bSJed Brown     ierr = DMDestroy(newdm);CHKERRQ(ierr);
229c4762a1bSJed Brown     *newdm = idm;
230c4762a1bSJed Brown   }
231c4762a1bSJed Brown   ierr = DMSetFromOptions(*newdm);CHKERRQ(ierr);
232c4762a1bSJed Brown   PetscFunctionReturn(0);
233c4762a1bSJed Brown }
234c4762a1bSJed Brown 
235c4762a1bSJed Brown static PetscErrorCode TestCone(DM dm, AppCtx *user)
236c4762a1bSJed Brown {
237c4762a1bSJed Brown   PetscInt           numRuns, cStart, cEnd, c, i;
238c4762a1bSJed Brown   PetscReal          maxTimePerRun = user->maxConeTime;
239c4762a1bSJed Brown   PetscLogStage      stage;
240c4762a1bSJed Brown   PetscLogEvent      event;
241c4762a1bSJed Brown   PetscEventPerfInfo eventInfo;
242c4762a1bSJed Brown   MPI_Comm           comm;
243c4762a1bSJed Brown   PetscMPIInt        rank;
244c4762a1bSJed Brown   PetscErrorCode     ierr;
245c4762a1bSJed Brown 
246c4762a1bSJed Brown   PetscFunctionBegin;
247c4762a1bSJed Brown   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
248ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
249c4762a1bSJed Brown   ierr = PetscLogStageRegister("DMPlex Cone Test", &stage);CHKERRQ(ierr);
250c4762a1bSJed Brown   ierr = PetscLogEventRegister("Cone", PETSC_OBJECT_CLASSID, &event);CHKERRQ(ierr);
251c4762a1bSJed Brown   ierr = PetscLogStagePush(stage);CHKERRQ(ierr);
252c4762a1bSJed Brown   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
253c4762a1bSJed Brown   ierr = PetscLogEventBegin(event,0,0,0,0);CHKERRQ(ierr);
254c4762a1bSJed Brown   for (i = 0; i < user->iterations; ++i) {
255c4762a1bSJed Brown     for (c = cStart; c < cEnd; ++c) {
256c4762a1bSJed Brown       const PetscInt *cone;
257c4762a1bSJed Brown 
258c4762a1bSJed Brown       ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr);
259c4762a1bSJed Brown     }
260c4762a1bSJed Brown   }
261c4762a1bSJed Brown   ierr = PetscLogEventEnd(event,0,0,0,0);CHKERRQ(ierr);
262c4762a1bSJed Brown   ierr = PetscLogStagePop();CHKERRQ(ierr);
263c4762a1bSJed Brown 
264c4762a1bSJed Brown   ierr = PetscLogEventGetPerfInfo(stage, event, &eventInfo);CHKERRQ(ierr);
265c4762a1bSJed Brown   numRuns = (cEnd-cStart) * user->iterations;
266c4762a1bSJed Brown   if (eventInfo.count != 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event calls %d should be %d", eventInfo.count, 1);
267c4762a1bSJed Brown   if ((PetscInt) eventInfo.flops != 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event flops %d should be %d", (PetscInt) eventInfo.flops, 0);
268c4762a1bSJed Brown   if (user->printTimes) {
269c4762a1bSJed Brown     ierr = PetscSynchronizedPrintf(comm, "[%d] Cones: %d Total time: %.3es Average time per cone: %.3es\n", rank, numRuns, eventInfo.time, eventInfo.time/numRuns);CHKERRQ(ierr);
270c4762a1bSJed Brown     ierr = PetscSynchronizedFlush(comm, PETSC_STDOUT);CHKERRQ(ierr);
271c4762a1bSJed Brown   } else if (eventInfo.time > maxTimePerRun * numRuns) {
272c4762a1bSJed Brown     ierr = PetscSynchronizedPrintf(comm, "[%d] Cones: %d Average time per cone: %gs standard: %gs\n", rank, numRuns, eventInfo.time/numRuns, maxTimePerRun);CHKERRQ(ierr);
273c4762a1bSJed Brown     ierr = PetscSynchronizedFlush(comm, PETSC_STDOUT);CHKERRQ(ierr);
274c4762a1bSJed Brown     if (user->errors) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Average time for cone %g > standard %g", eventInfo.time/numRuns, maxTimePerRun);
275c4762a1bSJed Brown   }
276c4762a1bSJed Brown   PetscFunctionReturn(0);
277c4762a1bSJed Brown }
278c4762a1bSJed Brown 
279c4762a1bSJed Brown static PetscErrorCode TestTransitiveClosure(DM dm, AppCtx *user)
280c4762a1bSJed Brown {
281c4762a1bSJed Brown   PetscInt           numRuns, cStart, cEnd, c, i;
282c4762a1bSJed Brown   PetscReal          maxTimePerRun = user->maxClosureTime;
283c4762a1bSJed Brown   PetscLogStage      stage;
284c4762a1bSJed Brown   PetscLogEvent      event;
285c4762a1bSJed Brown   PetscEventPerfInfo eventInfo;
286c4762a1bSJed Brown   MPI_Comm           comm;
287c4762a1bSJed Brown   PetscMPIInt        rank;
288c4762a1bSJed Brown   PetscErrorCode     ierr;
289c4762a1bSJed Brown 
290c4762a1bSJed Brown   PetscFunctionBegin;
291c4762a1bSJed Brown   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
292ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
293c4762a1bSJed Brown   ierr = PetscLogStageRegister("DMPlex Transitive Closure Test", &stage);CHKERRQ(ierr);
294c4762a1bSJed Brown   ierr = PetscLogEventRegister("TransitiveClosure", PETSC_OBJECT_CLASSID, &event);CHKERRQ(ierr);
295c4762a1bSJed Brown   ierr = PetscLogStagePush(stage);CHKERRQ(ierr);
296c4762a1bSJed Brown   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
297c4762a1bSJed Brown   ierr = PetscLogEventBegin(event,0,0,0,0);CHKERRQ(ierr);
298c4762a1bSJed Brown   for (i = 0; i < user->iterations; ++i) {
299c4762a1bSJed Brown     for (c = cStart; c < cEnd; ++c) {
300c4762a1bSJed Brown       PetscInt *closure = NULL;
301c4762a1bSJed Brown       PetscInt  closureSize;
302c4762a1bSJed Brown 
303c4762a1bSJed Brown       ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
304c4762a1bSJed Brown       ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
305c4762a1bSJed Brown     }
306c4762a1bSJed Brown   }
307c4762a1bSJed Brown   ierr = PetscLogEventEnd(event,0,0,0,0);CHKERRQ(ierr);
308c4762a1bSJed Brown   ierr = PetscLogStagePop();CHKERRQ(ierr);
309c4762a1bSJed Brown 
310c4762a1bSJed Brown   ierr = PetscLogEventGetPerfInfo(stage, event, &eventInfo);CHKERRQ(ierr);
311c4762a1bSJed Brown   numRuns = (cEnd-cStart) * user->iterations;
312c4762a1bSJed Brown   if (eventInfo.count != 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event calls %d should be %d", eventInfo.count, 1);
313c4762a1bSJed Brown   if ((PetscInt) eventInfo.flops != 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event flops %d should be %d", (PetscInt) eventInfo.flops, 0);
314c4762a1bSJed Brown   if (user->printTimes) {
315c4762a1bSJed Brown     ierr = PetscSynchronizedPrintf(comm, "[%d] Closures: %d Total time: %.3es Average time per cone: %.3es\n", rank, numRuns, eventInfo.time, eventInfo.time/numRuns);CHKERRQ(ierr);
316c4762a1bSJed Brown     ierr = PetscSynchronizedFlush(comm, PETSC_STDOUT);CHKERRQ(ierr);
317c4762a1bSJed Brown   } else if (eventInfo.time > maxTimePerRun * numRuns) {
318c4762a1bSJed Brown     ierr = PetscSynchronizedPrintf(comm, "[%d] Closures: %d Average time per cone: %gs standard: %gs\n", rank, numRuns, eventInfo.time/numRuns, maxTimePerRun);CHKERRQ(ierr);
319c4762a1bSJed Brown     ierr = PetscSynchronizedFlush(comm, PETSC_STDOUT);CHKERRQ(ierr);
320c4762a1bSJed Brown     if (user->errors) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Average time for closure %g > standard %g", eventInfo.time/numRuns, maxTimePerRun);
321c4762a1bSJed Brown   }
322c4762a1bSJed Brown   PetscFunctionReturn(0);
323c4762a1bSJed Brown }
324c4762a1bSJed Brown 
325c4762a1bSJed Brown static PetscErrorCode TestVecClosure(DM dm, PetscBool useIndex, PetscBool useSpectral, AppCtx *user)
326c4762a1bSJed Brown {
327c4762a1bSJed Brown   PetscSection       s;
328c4762a1bSJed Brown   Vec                v;
329c4762a1bSJed Brown   PetscInt           numRuns, cStart, cEnd, c, i;
330c4762a1bSJed Brown   PetscScalar        tmpArray[64];
331c4762a1bSJed Brown   PetscScalar       *userArray     = user->reuseArray ? tmpArray : NULL;
332c4762a1bSJed Brown   PetscReal          maxTimePerRun = user->maxVecClosureTime;
333c4762a1bSJed Brown   PetscLogStage      stage;
334c4762a1bSJed Brown   PetscLogEvent      event;
335c4762a1bSJed Brown   PetscEventPerfInfo eventInfo;
336c4762a1bSJed Brown   MPI_Comm           comm;
337c4762a1bSJed Brown   PetscMPIInt        rank;
338c4762a1bSJed Brown   PetscErrorCode     ierr;
339c4762a1bSJed Brown 
340c4762a1bSJed Brown   PetscFunctionBegin;
341c4762a1bSJed Brown   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
342ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
343c4762a1bSJed Brown   if (useIndex) {
344c4762a1bSJed Brown     if (useSpectral) {
345c4762a1bSJed Brown       ierr = PetscLogStageRegister("DMPlex Vector Closure with Index Test", &stage);CHKERRQ(ierr);
346c4762a1bSJed Brown       ierr = PetscLogEventRegister("VecClosureInd", PETSC_OBJECT_CLASSID, &event);CHKERRQ(ierr);
347c4762a1bSJed Brown     } else {
348c4762a1bSJed Brown       ierr = PetscLogStageRegister("DMPlex Vector Spectral Closure with Index Test", &stage);CHKERRQ(ierr);
349c4762a1bSJed Brown       ierr = PetscLogEventRegister("VecClosureSpecInd", PETSC_OBJECT_CLASSID, &event);CHKERRQ(ierr);
350c4762a1bSJed Brown     }
351c4762a1bSJed Brown   } else {
352c4762a1bSJed Brown     if (useSpectral) {
353c4762a1bSJed Brown       ierr = PetscLogStageRegister("DMPlex Vector Spectral Closure Test", &stage);CHKERRQ(ierr);
354c4762a1bSJed Brown       ierr = PetscLogEventRegister("VecClosureSpec", PETSC_OBJECT_CLASSID, &event);CHKERRQ(ierr);
355c4762a1bSJed Brown     } else {
356c4762a1bSJed Brown       ierr = PetscLogStageRegister("DMPlex Vector Closure Test", &stage);CHKERRQ(ierr);
357c4762a1bSJed Brown       ierr = PetscLogEventRegister("VecClosure", PETSC_OBJECT_CLASSID, &event);CHKERRQ(ierr);
358c4762a1bSJed Brown     }
359c4762a1bSJed Brown   }
360c4762a1bSJed Brown   ierr = PetscLogStagePush(stage);CHKERRQ(ierr);
361c4762a1bSJed Brown   ierr = DMSetNumFields(dm, user->numFields);CHKERRQ(ierr);
362c4762a1bSJed Brown   ierr = DMPlexCreateSection(dm, NULL, user->numComponents, user->numDof, 0, NULL, NULL, NULL, NULL, &s);CHKERRQ(ierr);
363c4762a1bSJed Brown   ierr = DMSetLocalSection(dm, s);CHKERRQ(ierr);
364c4762a1bSJed Brown   if (useIndex) {ierr = DMPlexCreateClosureIndex(dm, s);CHKERRQ(ierr);}
365c4762a1bSJed Brown   if (useSpectral) {ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, s);CHKERRQ(ierr);}
366c4762a1bSJed Brown   ierr = PetscSectionDestroy(&s);CHKERRQ(ierr);
367c4762a1bSJed Brown   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
368c4762a1bSJed Brown   ierr = DMGetLocalVector(dm, &v);CHKERRQ(ierr);
369c4762a1bSJed Brown   ierr = PetscLogEventBegin(event,0,0,0,0);CHKERRQ(ierr);
370c4762a1bSJed Brown   for (i = 0; i < user->iterations; ++i) {
371c4762a1bSJed Brown     for (c = cStart; c < cEnd; ++c) {
372c4762a1bSJed Brown       PetscScalar *closure     = userArray;
373c4762a1bSJed Brown       PetscInt     closureSize = 64;
374c4762a1bSJed Brown 
375c4762a1bSJed Brown       ierr = DMPlexVecGetClosure(dm, s, v, c, &closureSize, &closure);CHKERRQ(ierr);
376c4762a1bSJed Brown       if (!user->reuseArray) {ierr = DMPlexVecRestoreClosure(dm, s, v, c, &closureSize, &closure);CHKERRQ(ierr);}
377c4762a1bSJed Brown     }
378c4762a1bSJed Brown   }
379c4762a1bSJed Brown   ierr = PetscLogEventEnd(event,0,0,0,0);CHKERRQ(ierr);
380c4762a1bSJed Brown   ierr = DMRestoreLocalVector(dm, &v);CHKERRQ(ierr);
381c4762a1bSJed Brown   ierr = PetscLogStagePop();CHKERRQ(ierr);
382c4762a1bSJed Brown 
383c4762a1bSJed Brown   ierr = PetscLogEventGetPerfInfo(stage, event, &eventInfo);CHKERRQ(ierr);
384c4762a1bSJed Brown   numRuns = (cEnd-cStart) * user->iterations;
385c4762a1bSJed Brown   if (eventInfo.count != 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event calls %d should be %d", eventInfo.count, 1);
386c4762a1bSJed Brown   if ((PetscInt) eventInfo.flops != 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event flops %d should be %d", (PetscInt) eventInfo.flops, 0);
387c4762a1bSJed Brown   if (user->printTimes || eventInfo.time > maxTimePerRun * numRuns) {
388c4762a1bSJed Brown     const char *title = "VecClosures";
389c4762a1bSJed Brown     const char *titleIndex = "VecClosures with Index";
390c4762a1bSJed Brown     const char *titleSpec = "VecClosures Spectral";
391c4762a1bSJed Brown     const char *titleSpecIndex = "VecClosures Spectral with Index";
392c4762a1bSJed Brown 
393c4762a1bSJed Brown     if (user->printTimes) {
394c4762a1bSJed Brown       ierr = PetscSynchronizedPrintf(comm, "[%d] %s: %d Total time: %.3es Average time per vector closure: %.3es\n", rank, useIndex ? (useSpectral ? titleSpecIndex : titleIndex) : (useSpectral ? titleSpec : title), numRuns, eventInfo.time, eventInfo.time/numRuns);CHKERRQ(ierr);
395c4762a1bSJed Brown       ierr = PetscSynchronizedFlush(comm, PETSC_STDOUT);CHKERRQ(ierr);
396c4762a1bSJed Brown     } else {
397c4762a1bSJed Brown       ierr = PetscSynchronizedPrintf(comm, "[%d] %s: %d Average time per vector closure: %gs standard: %gs\n", rank, useIndex ? (useSpectral ? titleSpecIndex : titleIndex) : (useSpectral ? titleSpec : title), numRuns, eventInfo.time/numRuns, maxTimePerRun);CHKERRQ(ierr);
398c4762a1bSJed Brown       ierr = PetscSynchronizedFlush(comm, PETSC_STDOUT);CHKERRQ(ierr);
399c4762a1bSJed Brown       if (user->errors) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Average time for vector closure %g > standard %g", eventInfo.time/numRuns, maxTimePerRun);
400c4762a1bSJed Brown     }
401c4762a1bSJed Brown   }
402c4762a1bSJed Brown   PetscFunctionReturn(0);
403c4762a1bSJed Brown }
404c4762a1bSJed Brown 
405c4762a1bSJed Brown static PetscErrorCode CleanupContext(AppCtx *user)
406c4762a1bSJed Brown {
407c4762a1bSJed Brown   PetscErrorCode ierr;
408c4762a1bSJed Brown 
409c4762a1bSJed Brown   PetscFunctionBegin;
410c4762a1bSJed Brown   ierr = PetscFree(user->numComponents);CHKERRQ(ierr);
411c4762a1bSJed Brown   ierr = PetscFree(user->numDof);CHKERRQ(ierr);
412c4762a1bSJed Brown   PetscFunctionReturn(0);
413c4762a1bSJed Brown }
414c4762a1bSJed Brown 
415c4762a1bSJed Brown int main(int argc, char **argv)
416c4762a1bSJed Brown {
417c4762a1bSJed Brown   DM             dm;
418c4762a1bSJed Brown   AppCtx         user;
419c4762a1bSJed Brown   PetscErrorCode ierr;
420c4762a1bSJed Brown 
421c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
422c4762a1bSJed Brown   ierr = ProcessOptions(&user);CHKERRQ(ierr);
423c4762a1bSJed Brown   ierr = PetscLogDefaultBegin();CHKERRQ(ierr);
424c4762a1bSJed Brown   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
425c4762a1bSJed Brown   ierr = TestCone(dm, &user);CHKERRQ(ierr);
426c4762a1bSJed Brown   ierr = TestTransitiveClosure(dm, &user);CHKERRQ(ierr);
427c4762a1bSJed Brown   ierr = TestVecClosure(dm, PETSC_FALSE, PETSC_FALSE, &user);CHKERRQ(ierr);
428c4762a1bSJed Brown   ierr = TestVecClosure(dm, PETSC_TRUE,  PETSC_FALSE, &user);CHKERRQ(ierr);
429c4762a1bSJed Brown   if (!user.cellSimplex && user.spectral) {
430c4762a1bSJed Brown     ierr = TestVecClosure(dm, PETSC_FALSE, PETSC_TRUE,  &user);CHKERRQ(ierr);
431c4762a1bSJed Brown     ierr = TestVecClosure(dm, PETSC_TRUE,  PETSC_TRUE,  &user);CHKERRQ(ierr);
432c4762a1bSJed Brown   }
433c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
434c4762a1bSJed Brown   ierr = CleanupContext(&user);CHKERRQ(ierr);
435c4762a1bSJed Brown   ierr = PetscFinalize();
436c4762a1bSJed Brown   return ierr;
437c4762a1bSJed Brown }
438c4762a1bSJed Brown 
439c4762a1bSJed Brown /*TEST
440c4762a1bSJed Brown 
441c4762a1bSJed Brown   build:
442*dfd57a17SPierre Jolivet     requires: defined(PETSC_USE_LOG)
443c4762a1bSJed Brown 
444c4762a1bSJed Brown   # 2D Simplex P_1 scalar tests
445c4762a1bSJed Brown   testset:
446c4762a1bSJed Brown     args: -num_dof 1,0,0 -iterations 2 -print_times
447c4762a1bSJed Brown     test:
448c4762a1bSJed Brown       suffix: correctness_0
449c4762a1bSJed Brown     test:
450c4762a1bSJed Brown       suffix: correctness_1
451c4762a1bSJed Brown       args: -interpolate -dm_refine 2
452c4762a1bSJed Brown     test:
453c4762a1bSJed Brown       suffix: correctness_2
454c4762a1bSJed Brown       requires: triangle
455c4762a1bSJed Brown       args: -interpolate -refinement_limit 1.0e-5
456c4762a1bSJed Brown   test:
457c4762a1bSJed Brown     suffix: 0
458c4762a1bSJed Brown     TODO: Only for performance testing
459c4762a1bSJed 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
460c4762a1bSJed Brown   test:
461c4762a1bSJed Brown     suffix: 1
462c4762a1bSJed Brown     requires: triangle
463c4762a1bSJed Brown     TODO: Only for performance testing
464c4762a1bSJed 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
465c4762a1bSJed Brown   test:
466c4762a1bSJed Brown     suffix: 2
467c4762a1bSJed Brown     TODO: Only for performance testing
468c4762a1bSJed 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
469c4762a1bSJed Brown   test:
470c4762a1bSJed Brown     suffix: 3
471c4762a1bSJed Brown     requires: triangle
472c4762a1bSJed Brown     TODO: Only for performance testing
473c4762a1bSJed 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
474c4762a1bSJed Brown   test:
475c4762a1bSJed Brown     suffix: 4
476c4762a1bSJed Brown     TODO: Only for performance testing
477c4762a1bSJed 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
478c4762a1bSJed Brown   test:
479c4762a1bSJed Brown     suffix: 5
480c4762a1bSJed Brown     requires: triangle
481c4762a1bSJed Brown     TODO: Only for performance testing
482c4762a1bSJed 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
483c4762a1bSJed Brown   test:
484c4762a1bSJed Brown     suffix: 6
485c4762a1bSJed Brown     TODO: Only for performance testing
486c4762a1bSJed 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
487c4762a1bSJed Brown   test:
488c4762a1bSJed Brown     suffix: 7
489c4762a1bSJed Brown     requires: triangle
490c4762a1bSJed Brown     TODO: Only for performance testing
491c4762a1bSJed 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
492c4762a1bSJed Brown 
493c4762a1bSJed Brown   # 2D Simplex P_1 vector tests
494c4762a1bSJed Brown   # 2D Simplex P_2 scalar tests
495c4762a1bSJed Brown   # 2D Simplex P_2 vector tests
496c4762a1bSJed Brown   # 2D Simplex P_2/P_1 vector/scalar tests
497c4762a1bSJed Brown   # 2D Quad P_1 scalar tests
498c4762a1bSJed Brown   # 2D Quad P_1 vector tests
499c4762a1bSJed Brown   # 2D Quad P_2 scalar tests
500c4762a1bSJed Brown   # 2D Quad P_2 vector tests
501c4762a1bSJed Brown   # 3D Simplex P_1 scalar tests
502c4762a1bSJed Brown   # 3D Simplex P_1 vector tests
503c4762a1bSJed Brown   # 3D Simplex P_2 scalar tests
504c4762a1bSJed Brown   # 3D Simplex P_2 vector tests
505c4762a1bSJed Brown   # 3D Hex P_1 scalar tests
506c4762a1bSJed Brown   # 3D Hex P_1 vector tests
507c4762a1bSJed Brown   # 3D Hex P_2 scalar tests
508c4762a1bSJed Brown   # 3D Hex P_2 vector tests
509c4762a1bSJed Brown 
510c4762a1bSJed Brown TEST*/
511