Actual source code: ex9.c

petsc-3.14.6 2021-03-30
Report Typos and Errors
  1: static char help[] = "Performance tests for DMPlex query operations\n\n";

  3: #include <petscdmplex.h>

  5: typedef struct {
  6:   PetscInt  dim;             /* The topological mesh dimension */
  7:   PetscBool cellSimplex;     /* Flag for simplices */
  8:   PetscBool spectral;        /* Flag for spectral element layout */
  9:   PetscBool interpolate;     /* Flag for mesh interpolation */
 10:   PetscReal refinementLimit; /* Maximum volume of a refined cell */
 11:   PetscInt  numFields;       /* The number of section fields */
 12:   PetscInt *numComponents;   /* The number of field components */
 13:   PetscInt *numDof;          /* The dof signature for the section */
 14:   PetscBool reuseArray;      /* Pass in user allocated array to VecGetClosure() */
 15:   /* Test data */
 16:   PetscBool errors;            /* Treat failures as errors */
 17:   PetscInt  iterations;        /* The number of iterations for a query */
 18:   PetscReal maxConeTime;       /* Max time per run for DMPlexGetCone() */
 19:   PetscReal maxClosureTime;    /* Max time per run for DMPlexGetTransitiveClosure() */
 20:   PetscReal maxVecClosureTime; /* Max time per run for DMPlexVecGetClosure() */
 21:   PetscBool printTimes;        /* Print total times, do not check limits */
 22: } AppCtx;

 24: static PetscErrorCode ProcessOptions(AppCtx *options)
 25: {
 26:   PetscInt       len;
 27:   PetscBool      flg;

 31:   options->dim               = 2;
 32:   options->cellSimplex       = PETSC_TRUE;
 33:   options->spectral          = PETSC_FALSE;
 34:   options->interpolate       = PETSC_FALSE;
 35:   options->refinementLimit   = 0.0;
 36:   options->numFields         = 0;
 37:   options->numComponents     = NULL;
 38:   options->numDof            = NULL;
 39:   options->reuseArray        = PETSC_FALSE;
 40:   options->errors            = PETSC_FALSE;
 41:   options->iterations        = 1;
 42:   options->maxConeTime       = 0.0;
 43:   options->maxClosureTime    = 0.0;
 44:   options->maxVecClosureTime = 0.0;
 45:   options->printTimes        = PETSC_FALSE;

 47:   PetscOptionsBegin(PETSC_COMM_SELF, "", "Meshing Problem Options", "DMPLEX");
 48:   PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex9.c", options->dim, &options->dim, NULL,1,3);
 49:   PetscOptionsBool("-cellSimplex", "Flag for simplices", "ex9.c", options->cellSimplex, &options->cellSimplex, NULL);
 50:   PetscOptionsBool("-spectral", "Flag for spectral element layout", "ex9.c", options->spectral, &options->spectral, NULL);
 51:   PetscOptionsBool("-interpolate", "Flag for mesh interpolation", "ex9.c", options->interpolate, &options->interpolate, NULL);
 52:   PetscOptionsReal("-refinement_limit", "The maximum volume of a refined cell", "ex9.c", options->refinementLimit, &options->refinementLimit, NULL);
 53:   PetscOptionsBoundedInt("-num_fields", "The number of section fields", "ex9.c", options->numFields, &options->numFields, NULL, 0);
 54:   if (options->numFields) {
 55:     len  = options->numFields;
 56:     PetscMalloc1(len, &options->numComponents);
 57:     PetscOptionsIntArray("-num_components", "The number of components per field", "ex9.c", options->numComponents, &len, &flg);
 58:     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);
 59:   }
 60:   len  = (options->dim+1) * PetscMax(1, options->numFields);
 61:   PetscMalloc1(len, &options->numDof);
 62:   PetscOptionsIntArray("-num_dof", "The dof signature for the section", "ex9.c", options->numDof, &len, &flg);
 63:   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));

 65:   /* We are specifying the scalar dof, so augment it for multiple components */
 66:   {
 67:     PetscInt f, d;

 69:     for (f = 0; f < options->numFields; ++f) {
 70:       for (d = 0; d <= options->dim; ++d) options->numDof[f*(options->dim+1)+d] *= options->numComponents[f];
 71:     }
 72:   }

 74:   PetscOptionsBool("-reuse_array", "Pass in user allocated array to VecGetClosure()", "ex9.c", options->reuseArray, &options->reuseArray, NULL);
 75:   PetscOptionsBool("-errors", "Treat failures as errors", "ex9.c", options->errors, &options->errors, NULL);
 76:   PetscOptionsBoundedInt("-iterations", "The number of iterations for a query", "ex9.c", options->iterations, &options->iterations, NULL,0);
 77:   PetscOptionsReal("-max_cone_time", "The maximum time per run for DMPlexGetCone()", "ex9.c", options->maxConeTime, &options->maxConeTime, NULL);
 78:   PetscOptionsReal("-max_closure_time", "The maximum time per run for DMPlexGetTransitiveClosure()", "ex9.c", options->maxClosureTime, &options->maxClosureTime, NULL);
 79:   PetscOptionsReal("-max_vec_closure_time", "The maximum time per run for DMPlexVecGetClosure()", "ex9.c", options->maxVecClosureTime, &options->maxVecClosureTime, NULL);
 80:   PetscOptionsBool("-print_times", "Print total times, do not check limits", "ex9.c", options->printTimes, &options->printTimes, NULL);
 81:   PetscOptionsEnd();
 82:   return(0);
 83: }

 85: static PetscErrorCode CreateSimplex_2D(MPI_Comm comm, DM *newdm)
 86: {
 87:   DM             dm;
 88:   PetscInt       numPoints[2]        = {4, 2};
 89:   PetscInt       coneSize[6]         = {3, 3, 0, 0, 0, 0};
 90:   PetscInt       cones[6]            = {2, 3, 4,  5, 4, 3};
 91:   PetscInt       coneOrientations[6] = {0, 0, 0,  0, 0, 0};
 92:   PetscScalar    vertexCoords[8]     = {-0.5, 0.5, 0.0, 0.0, 0.0, 1.0, 0.5, 0.5};
 93:   PetscInt       markerPoints[8]     = {2, 1, 3, 1, 4, 1, 5, 1};
 94:   PetscInt       dim = 2, depth = 1, p;

 98:   DMCreate(comm, &dm);
 99:   PetscObjectSetName((PetscObject) dm, "triangular");
100:   DMSetType(dm, DMPLEX);
101:   DMSetDimension(dm, dim);
102:   DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords);
103:   for (p = 0; p < 4; ++p) {
104:     DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);
105:   }
106:   *newdm = dm;
107:   return(0);
108: }

110: static PetscErrorCode CreateSimplex_3D(MPI_Comm comm, DM *newdm)
111: {
112:   DM             dm;
113:   PetscInt       numPoints[2]        = {5, 2};
114:   PetscInt       coneSize[23]        = {4, 4, 0, 0, 0, 0, 0};
115:   PetscInt       cones[8]            = {2, 4, 3, 5,  3, 4, 6, 5};
116:   PetscInt       coneOrientations[8] = {0, 0, 0, 0,  0, 0, 0, 0};
117:   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};
118:   PetscInt       markerPoints[10]    = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1};
119:   PetscInt       dim = 3, depth = 1, p;

123:   DMCreate(comm, &dm);
124:   PetscObjectSetName((PetscObject) dm, "tetrahedral");
125:   DMSetType(dm, DMPLEX);
126:   DMSetDimension(dm, dim);
127:   DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords);
128:   for (p = 0; p < 5; ++p) {
129:     DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);
130:   }
131:   *newdm = dm;
132:   return(0);
133: }

135: static PetscErrorCode CreateQuad_2D(MPI_Comm comm, DM *newdm)
136: {
137:   DM             dm;
138:   PetscInt       numPoints[2]        = {6, 2};
139:   PetscInt       coneSize[8]         = {4, 4, 0, 0, 0, 0, 0, 0};
140:   PetscInt       cones[8]            = {2, 3, 4, 5,  3, 6, 7, 4};
141:   PetscInt       coneOrientations[8] = {0, 0, 0, 0,  0, 0, 0, 0};
142:   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};
143:   PetscInt       markerPoints[12]    = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1, 7, 1};
144:   PetscInt       dim = 2, depth = 1, p;

148:   DMCreate(comm, &dm);
149:   PetscObjectSetName((PetscObject) dm, "quadrilateral");
150:   DMSetType(dm, DMPLEX);
151:   DMSetDimension(dm, dim);
152:   DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords);
153:   for (p = 0; p < 6; ++p) {
154:     DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);
155:   }
156:   *newdm = dm;
157:   return(0);
158: }

160: static PetscErrorCode CreateHex_3D(MPI_Comm comm, DM *newdm)
161: {
162:   DM             dm;
163:   PetscInt       numPoints[2]         = {12, 2};
164:   PetscInt       coneSize[14]         = {8, 8, 0,0,0,0,0,0,0,0,0,0,0,0};
165:   PetscInt       cones[16]            = {2,5,4,3,6,7,8,9,  3,4,11,10,7,12,13,8};
166:   PetscInt       coneOrientations[16] = {0,0,0,0,0,0,0,0,  0,0, 0, 0,0, 0, 0,0};
167:   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,
168:                                          -0.5,0.0,1.0, 0.0,0.0,1.0, 0.0,1.0,1.0, -0.5,1.0,1.0,
169:                                           0.5,0.0,0.0, 0.5,1.0,0.0, 0.5,0.0,1.0,  0.5,1.0,1.0};
170:   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};
171:   PetscInt       dim = 3, depth = 1, p;

175:   DMCreate(comm, &dm);
176:   PetscObjectSetName((PetscObject) dm, "hexahedral");
177:   DMSetType(dm, DMPLEX);
178:   DMSetDimension(dm, dim);
179:   DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords);
180:   for (p = 0; p < 12; ++p) {
181:     DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);
182:   }
183:   *newdm = dm;
184:   return(0);
185: }

187: static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *newdm)
188: {
189:   PetscInt       dim         = user->dim;
190:   PetscBool      cellSimplex = user->cellSimplex;

194:   switch (dim) {
195:   case 2:
196:     if (cellSimplex) {
197:       CreateSimplex_2D(comm, newdm);
198:     } else {
199:       CreateQuad_2D(comm, newdm);
200:     }
201:     break;
202:   case 3:
203:     if (cellSimplex) {
204:       CreateSimplex_3D(comm, newdm);
205:     } else {
206:       CreateHex_3D(comm, newdm);
207:     }
208:     break;
209:   default:
210:     SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %d", dim);
211:   }
212:   if (user->refinementLimit > 0.0) {
213:     DM rdm;
214:     const char *name;

216:     DMPlexSetRefinementUniform(*newdm, PETSC_FALSE);
217:     DMPlexSetRefinementLimit(*newdm, user->refinementLimit);
218:     DMRefine(*newdm, PETSC_COMM_SELF, &rdm);
219:     PetscObjectGetName((PetscObject) *newdm, &name);
220:     PetscObjectSetName((PetscObject)    rdm,  name);
221:     DMDestroy(newdm);
222:     *newdm = rdm;
223:   }
224:   if (user->interpolate) {
225:     DM idm;

227:     DMPlexInterpolate(*newdm, &idm);
228:     DMDestroy(newdm);
229:     *newdm = idm;
230:   }
231:   DMSetFromOptions(*newdm);
232:   return(0);
233: }

235: static PetscErrorCode TestCone(DM dm, AppCtx *user)
236: {
237:   PetscInt           numRuns, cStart, cEnd, c, i;
238:   PetscReal          maxTimePerRun = user->maxConeTime;
239:   PetscLogStage      stage;
240:   PetscLogEvent      event;
241:   PetscEventPerfInfo eventInfo;
242:   MPI_Comm           comm;
243:   PetscMPIInt        rank;
244:   PetscErrorCode     ierr;

247:   PetscObjectGetComm((PetscObject)dm, &comm);
248:   MPI_Comm_rank(comm, &rank);
249:   PetscLogStageRegister("DMPlex Cone Test", &stage);
250:   PetscLogEventRegister("Cone", PETSC_OBJECT_CLASSID, &event);
251:   PetscLogStagePush(stage);
252:   DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
253:   PetscLogEventBegin(event,0,0,0,0);
254:   for (i = 0; i < user->iterations; ++i) {
255:     for (c = cStart; c < cEnd; ++c) {
256:       const PetscInt *cone;

258:       DMPlexGetCone(dm, c, &cone);
259:     }
260:   }
261:   PetscLogEventEnd(event,0,0,0,0);
262:   PetscLogStagePop();

264:   PetscLogEventGetPerfInfo(stage, event, &eventInfo);
265:   numRuns = (cEnd-cStart) * user->iterations;
266:   if (eventInfo.count != 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event calls %d should be %d", eventInfo.count, 1);
267:   if ((PetscInt) eventInfo.flops != 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event flops %d should be %d", (PetscInt) eventInfo.flops, 0);
268:   if (user->printTimes) {
269:     PetscSynchronizedPrintf(comm, "[%d] Cones: %d Total time: %.3es Average time per cone: %.3es\n", rank, numRuns, eventInfo.time, eventInfo.time/numRuns);
270:     PetscSynchronizedFlush(comm, PETSC_STDOUT);
271:   } else if (eventInfo.time > maxTimePerRun * numRuns) {
272:     PetscSynchronizedPrintf(comm, "[%d] Cones: %d Average time per cone: %gs standard: %gs\n", rank, numRuns, eventInfo.time/numRuns, maxTimePerRun);
273:     PetscSynchronizedFlush(comm, PETSC_STDOUT);
274:     if (user->errors) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Average time for cone %g > standard %g", eventInfo.time/numRuns, maxTimePerRun);
275:   }
276:   return(0);
277: }

279: static PetscErrorCode TestTransitiveClosure(DM dm, AppCtx *user)
280: {
281:   PetscInt           numRuns, cStart, cEnd, c, i;
282:   PetscReal          maxTimePerRun = user->maxClosureTime;
283:   PetscLogStage      stage;
284:   PetscLogEvent      event;
285:   PetscEventPerfInfo eventInfo;
286:   MPI_Comm           comm;
287:   PetscMPIInt        rank;
288:   PetscErrorCode     ierr;

291:   PetscObjectGetComm((PetscObject)dm, &comm);
292:   MPI_Comm_rank(comm, &rank);
293:   PetscLogStageRegister("DMPlex Transitive Closure Test", &stage);
294:   PetscLogEventRegister("TransitiveClosure", PETSC_OBJECT_CLASSID, &event);
295:   PetscLogStagePush(stage);
296:   DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
297:   PetscLogEventBegin(event,0,0,0,0);
298:   for (i = 0; i < user->iterations; ++i) {
299:     for (c = cStart; c < cEnd; ++c) {
300:       PetscInt *closure = NULL;
301:       PetscInt  closureSize;

303:       DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
304:       DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
305:     }
306:   }
307:   PetscLogEventEnd(event,0,0,0,0);
308:   PetscLogStagePop();

310:   PetscLogEventGetPerfInfo(stage, event, &eventInfo);
311:   numRuns = (cEnd-cStart) * user->iterations;
312:   if (eventInfo.count != 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event calls %d should be %d", eventInfo.count, 1);
313:   if ((PetscInt) eventInfo.flops != 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event flops %d should be %d", (PetscInt) eventInfo.flops, 0);
314:   if (user->printTimes) {
315:     PetscSynchronizedPrintf(comm, "[%d] Closures: %d Total time: %.3es Average time per cone: %.3es\n", rank, numRuns, eventInfo.time, eventInfo.time/numRuns);
316:     PetscSynchronizedFlush(comm, PETSC_STDOUT);
317:   } else if (eventInfo.time > maxTimePerRun * numRuns) {
318:     PetscSynchronizedPrintf(comm, "[%d] Closures: %d Average time per cone: %gs standard: %gs\n", rank, numRuns, eventInfo.time/numRuns, maxTimePerRun);
319:     PetscSynchronizedFlush(comm, PETSC_STDOUT);
320:     if (user->errors) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Average time for closure %g > standard %g", eventInfo.time/numRuns, maxTimePerRun);
321:   }
322:   return(0);
323: }

325: static PetscErrorCode TestVecClosure(DM dm, PetscBool useIndex, PetscBool useSpectral, AppCtx *user)
326: {
327:   PetscSection       s;
328:   Vec                v;
329:   PetscInt           numRuns, cStart, cEnd, c, i;
330:   PetscScalar        tmpArray[64];
331:   PetscScalar       *userArray     = user->reuseArray ? tmpArray : NULL;
332:   PetscReal          maxTimePerRun = user->maxVecClosureTime;
333:   PetscLogStage      stage;
334:   PetscLogEvent      event;
335:   PetscEventPerfInfo eventInfo;
336:   MPI_Comm           comm;
337:   PetscMPIInt        rank;
338:   PetscErrorCode     ierr;

341:   PetscObjectGetComm((PetscObject)dm, &comm);
342:   MPI_Comm_rank(comm, &rank);
343:   if (useIndex) {
344:     if (useSpectral) {
345:       PetscLogStageRegister("DMPlex Vector Closure with Index Test", &stage);
346:       PetscLogEventRegister("VecClosureInd", PETSC_OBJECT_CLASSID, &event);
347:     } else {
348:       PetscLogStageRegister("DMPlex Vector Spectral Closure with Index Test", &stage);
349:       PetscLogEventRegister("VecClosureSpecInd", PETSC_OBJECT_CLASSID, &event);
350:     }
351:   } else {
352:     if (useSpectral) {
353:       PetscLogStageRegister("DMPlex Vector Spectral Closure Test", &stage);
354:       PetscLogEventRegister("VecClosureSpec", PETSC_OBJECT_CLASSID, &event);
355:     } else {
356:       PetscLogStageRegister("DMPlex Vector Closure Test", &stage);
357:       PetscLogEventRegister("VecClosure", PETSC_OBJECT_CLASSID, &event);
358:     }
359:   }
360:   PetscLogStagePush(stage);
361:   DMSetNumFields(dm, user->numFields);
362:   DMPlexCreateSection(dm, NULL, user->numComponents, user->numDof, 0, NULL, NULL, NULL, NULL, &s);
363:   DMSetLocalSection(dm, s);
364:   if (useIndex) {DMPlexCreateClosureIndex(dm, s);}
365:   if (useSpectral) {DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, s);}
366:   PetscSectionDestroy(&s);
367:   DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
368:   DMGetLocalVector(dm, &v);
369:   PetscLogEventBegin(event,0,0,0,0);
370:   for (i = 0; i < user->iterations; ++i) {
371:     for (c = cStart; c < cEnd; ++c) {
372:       PetscScalar *closure     = userArray;
373:       PetscInt     closureSize = 64;

375:       DMPlexVecGetClosure(dm, s, v, c, &closureSize, &closure);
376:       if (!user->reuseArray) {DMPlexVecRestoreClosure(dm, s, v, c, &closureSize, &closure);}
377:     }
378:   }
379:   PetscLogEventEnd(event,0,0,0,0);
380:   DMRestoreLocalVector(dm, &v);
381:   PetscLogStagePop();

383:   PetscLogEventGetPerfInfo(stage, event, &eventInfo);
384:   numRuns = (cEnd-cStart) * user->iterations;
385:   if (eventInfo.count != 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event calls %d should be %d", eventInfo.count, 1);
386:   if ((PetscInt) eventInfo.flops != 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event flops %d should be %d", (PetscInt) eventInfo.flops, 0);
387:   if (user->printTimes || eventInfo.time > maxTimePerRun * numRuns) {
388:     const char *title = "VecClosures";
389:     const char *titleIndex = "VecClosures with Index";
390:     const char *titleSpec = "VecClosures Spectral";
391:     const char *titleSpecIndex = "VecClosures Spectral with Index";

393:     if (user->printTimes) {
394:       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);
395:       PetscSynchronizedFlush(comm, PETSC_STDOUT);
396:     } else {
397:       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);
398:       PetscSynchronizedFlush(comm, PETSC_STDOUT);
399:       if (user->errors) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Average time for vector closure %g > standard %g", eventInfo.time/numRuns, maxTimePerRun);
400:     }
401:   }
402:   return(0);
403: }

405: static PetscErrorCode CleanupContext(AppCtx *user)
406: {

410:   PetscFree(user->numComponents);
411:   PetscFree(user->numDof);
412:   return(0);
413: }

415: int main(int argc, char **argv)
416: {
417:   DM             dm;
418:   AppCtx         user;

421:   PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
422:   ProcessOptions(&user);
423:   PetscLogDefaultBegin();
424:   CreateMesh(PETSC_COMM_WORLD, &user, &dm);
425:   TestCone(dm, &user);
426:   TestTransitiveClosure(dm, &user);
427:   TestVecClosure(dm, PETSC_FALSE, PETSC_FALSE, &user);
428:   TestVecClosure(dm, PETSC_TRUE,  PETSC_FALSE, &user);
429:   if (!user.cellSimplex && user.spectral) {
430:     TestVecClosure(dm, PETSC_FALSE, PETSC_TRUE,  &user);
431:     TestVecClosure(dm, PETSC_TRUE,  PETSC_TRUE,  &user);
432:   }
433:   DMDestroy(&dm);
434:   CleanupContext(&user);
435:   PetscFinalize();
436:   return ierr;
437: }

439: /*TEST

441:   build:
442:     requires: define(PETSC_USE_LOG)

444:   # 2D Simplex P_1 scalar tests
445:   testset:
446:     args: -num_dof 1,0,0 -iterations 2 -print_times
447:     test:
448:       suffix: correctness_0
449:     test:
450:       suffix: correctness_1
451:       args: -interpolate -dm_refine 2
452:     test:
453:       suffix: correctness_2
454:       requires: triangle
455:       args: -interpolate -refinement_limit 1.0e-5
456:   test:
457:     suffix: 0
458:     TODO: Only for performance testing
459:     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
460:   test:
461:     suffix: 1
462:     requires: triangle
463:     TODO: Only for performance testing
464:     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
465:   test:
466:     suffix: 2
467:     TODO: Only for performance testing
468:     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
469:   test:
470:     suffix: 3
471:     requires: triangle
472:     TODO: Only for performance testing
473:     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
474:   test:
475:     suffix: 4
476:     TODO: Only for performance testing
477:     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
478:   test:
479:     suffix: 5
480:     requires: triangle
481:     TODO: Only for performance testing
482:     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
483:   test:
484:     suffix: 6
485:     TODO: Only for performance testing
486:     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
487:   test:
488:     suffix: 7
489:     requires: triangle
490:     TODO: Only for performance testing
491:     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

493:   # 2D Simplex P_1 vector tests
494:   # 2D Simplex P_2 scalar tests
495:   # 2D Simplex P_2 vector tests
496:   # 2D Simplex P_2/P_1 vector/scalar tests
497:   # 2D Quad P_1 scalar tests
498:   # 2D Quad P_1 vector tests
499:   # 2D Quad P_2 scalar tests
500:   # 2D Quad P_2 vector tests
501:   # 3D Simplex P_1 scalar tests
502:   # 3D Simplex P_1 vector tests
503:   # 3D Simplex P_2 scalar tests
504:   # 3D Simplex P_2 vector tests
505:   # 3D Hex P_1 scalar tests
506:   # 3D Hex P_1 vector tests
507:   # 3D Hex P_2 scalar tests
508:   # 3D Hex P_2 vector tests

510: TEST*/