Actual source code: ex9.c

petsc-3.12.5 2020-03-29
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: } AppCtx;

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

224:     DMPlexInterpolate(*newdm, &idm);
225:     DMDestroy(newdm);
226:     *newdm = idm;
227:   }
228:   DMSetFromOptions(*newdm);
229:   return(0);
230: }

232: static PetscErrorCode TestCone(DM dm, AppCtx *user)
233: {
234:   PetscInt           numRuns, cStart, cEnd, c, i;
235:   PetscReal          maxTimePerRun = user->maxConeTime;
236:   PetscLogStage      stage;
237:   PetscLogEvent      event;
238:   PetscEventPerfInfo eventInfo;
239:   PetscErrorCode     ierr;

242:   PetscLogStageRegister("DMPlex Cone Test", &stage);
243:   PetscLogEventRegister("Cone", PETSC_OBJECT_CLASSID, &event);
244:   PetscLogStagePush(stage);
245:   DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
246:   PetscLogEventBegin(event,0,0,0,0);
247:   for (i = 0; i < user->iterations; ++i) {
248:     for (c = cStart; c < cEnd; ++c) {
249:       const PetscInt *cone;

251:       DMPlexGetCone(dm, c, &cone);
252:     }
253:   }
254:   PetscLogEventEnd(event,0,0,0,0);
255:   PetscLogStagePop();

257:   PetscLogEventGetPerfInfo(stage, event, &eventInfo);
258:   numRuns = (cEnd-cStart) * user->iterations;
259:   if (eventInfo.count != 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event calls %d should be %d", eventInfo.count, 1);
260:   if ((PetscInt) eventInfo.flops != 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event flops %d should be %d", (PetscInt) eventInfo.flops, 0);
261:   if (eventInfo.time > maxTimePerRun * numRuns) {
262:     PetscPrintf(PETSC_COMM_SELF, "Cones: %d Average time per cone: %gs standard: %gs\n", numRuns, eventInfo.time/numRuns, maxTimePerRun);
263:     if (user->errors) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Average time for cone %g > standard %g", eventInfo.time/numRuns, maxTimePerRun);
264:   }
265:   return(0);
266: }

268: static PetscErrorCode TestTransitiveClosure(DM dm, AppCtx *user)
269: {
270:   PetscInt           numRuns, cStart, cEnd, c, i;
271:   PetscReal          maxTimePerRun = user->maxClosureTime;
272:   PetscLogStage      stage;
273:   PetscLogEvent      event;
274:   PetscEventPerfInfo eventInfo;
275:   PetscErrorCode     ierr;

278:   PetscLogStageRegister("DMPlex Transitive Closure Test", &stage);
279:   PetscLogEventRegister("TransitiveClosure", PETSC_OBJECT_CLASSID, &event);
280:   PetscLogStagePush(stage);
281:   DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
282:   PetscLogEventBegin(event,0,0,0,0);
283:   for (i = 0; i < user->iterations; ++i) {
284:     for (c = cStart; c < cEnd; ++c) {
285:       PetscInt *closure = NULL;
286:       PetscInt  closureSize;

288:       DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
289:       DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
290:     }
291:   }
292:   PetscLogEventEnd(event,0,0,0,0);
293:   PetscLogStagePop();

295:   PetscLogEventGetPerfInfo(stage, event, &eventInfo);
296:   numRuns = (cEnd-cStart) * user->iterations;
297:   if (eventInfo.count != 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event calls %d should be %d", eventInfo.count, 1);
298:   if ((PetscInt) eventInfo.flops != 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event flops %d should be %d", (PetscInt) eventInfo.flops, 0);
299:   if (eventInfo.time > maxTimePerRun * numRuns) {
300:     PetscPrintf(PETSC_COMM_SELF, "Closures: %d Average time per cone: %gs standard: %gs\n", numRuns, eventInfo.time/numRuns, maxTimePerRun);
301:     if (user->errors) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Average time for closure %g > standard %g", eventInfo.time/numRuns, maxTimePerRun);
302:   }
303:   return(0);
304: }

306: static PetscErrorCode TestVecClosure(DM dm, PetscBool useIndex, PetscBool useSpectral, AppCtx *user)
307: {
308:   PetscSection       s;
309:   Vec                v;
310:   PetscInt           numRuns, cStart, cEnd, c, i;
311:   PetscScalar        tmpArray[64];
312:   PetscScalar       *userArray     = user->reuseArray ? tmpArray : NULL;
313:   PetscReal          maxTimePerRun = user->maxVecClosureTime;
314:   PetscLogStage      stage;
315:   PetscLogEvent      event;
316:   PetscEventPerfInfo eventInfo;
317:   PetscErrorCode     ierr;

320:   if (useIndex) {
321:     if (useSpectral) {
322:       PetscLogStageRegister("DMPlex Vector Closure with Index Test", &stage);
323:       PetscLogEventRegister("VecClosureInd", PETSC_OBJECT_CLASSID, &event);
324:     } else {
325:       PetscLogStageRegister("DMPlex Vector Spectral Closure with Index Test", &stage);
326:       PetscLogEventRegister("VecClosureSpecInd", PETSC_OBJECT_CLASSID, &event);
327:     }
328:   } else {
329:     if (useSpectral) {
330:       PetscLogStageRegister("DMPlex Vector Spectral Closure Test", &stage);
331:       PetscLogEventRegister("VecClosureSpec", PETSC_OBJECT_CLASSID, &event);
332:     } else {
333:       PetscLogStageRegister("DMPlex Vector Closure Test", &stage);
334:       PetscLogEventRegister("VecClosure", PETSC_OBJECT_CLASSID, &event);
335:     }
336:   }
337:   PetscLogStagePush(stage);
338:   DMSetNumFields(dm, user->numFields);
339:   DMPlexCreateSection(dm, NULL, user->numComponents, user->numDof, 0, NULL, NULL, NULL, NULL, &s);
340:   DMSetLocalSection(dm, s);
341:   if (useIndex) {DMPlexCreateClosureIndex(dm, s);}
342:   if (useSpectral) {DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, s);}
343:   PetscSectionDestroy(&s);
344:   DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
345:   DMGetLocalVector(dm, &v);
346:   PetscLogEventBegin(event,0,0,0,0);
347:   for (i = 0; i < user->iterations; ++i) {
348:     for (c = cStart; c < cEnd; ++c) {
349:       PetscScalar *closure     = userArray;
350:       PetscInt     closureSize = 64;

352:       DMPlexVecGetClosure(dm, s, v, c, &closureSize, &closure);
353:       if (!user->reuseArray) {DMPlexVecRestoreClosure(dm, s, v, c, &closureSize, &closure);}
354:     }
355:   }
356:   PetscLogEventEnd(event,0,0,0,0);
357:   DMRestoreLocalVector(dm, &v);
358:   PetscLogStagePop();

360:   PetscLogEventGetPerfInfo(stage, event, &eventInfo);
361:   numRuns = (cEnd-cStart) * user->iterations;
362:   if (eventInfo.count != 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event calls %d should be %d", eventInfo.count, 1);
363:   if ((PetscInt) eventInfo.flops != 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event flops %d should be %d", (PetscInt) eventInfo.flops, 0);
364:   if (eventInfo.time > maxTimePerRun * numRuns) {
365:     const char *title = "VecClosures";
366:     const char *titleIndex = "VecClosures with Index";
367:     const char *titleSpec = "VecClosures Spectral";
368:     const char *titleSpecIndex = "VecClosures Spectral with Index";

370:     PetscPrintf(PETSC_COMM_SELF, "%s: %d Average time per vector closure: %gs standard: %gs\n", useIndex ? (useSpectral ? titleSpecIndex : titleIndex) : (useSpectral ? titleSpec : title), numRuns, eventInfo.time/numRuns, maxTimePerRun);
371:     if (user->errors) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Average time for vector closure %g > standard %g", eventInfo.time/numRuns, maxTimePerRun);
372:   }
373:   return(0);
374: }

376: static PetscErrorCode CleanupContext(AppCtx *user)
377: {

381:   PetscFree(user->numComponents);
382:   PetscFree(user->numDof);
383:   return(0);
384: }

386: int main(int argc, char **argv)
387: {
388:   DM             dm;
389:   AppCtx         user;

392:   PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
393:   ProcessOptions(&user);
394:   PetscLogDefaultBegin();
395:   CreateMesh(PETSC_COMM_SELF, &user, &dm);
396:   TestCone(dm, &user);
397:   TestTransitiveClosure(dm, &user);
398:   TestVecClosure(dm, PETSC_FALSE, PETSC_FALSE, &user);
399:   TestVecClosure(dm, PETSC_TRUE,  PETSC_FALSE, &user);
400:   if (!user.cellSimplex && user.spectral) {
401:     TestVecClosure(dm, PETSC_FALSE, PETSC_TRUE,  &user);
402:     TestVecClosure(dm, PETSC_TRUE,  PETSC_TRUE,  &user);
403:   }
404:   DMDestroy(&dm);
405:   CleanupContext(&user);
406:   PetscFinalize();
407:   return ierr;
408: }

410: /*TEST

412:   # 2D Simplex P_1 scalar tests
413:   test:
414:     suffix: 0
415:     requires: performance
416:     TODO: missing output file
417:     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
418:   test:
419:     suffix: 1
420:     requires: performance
421:     TODO: missing output file
422:     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
423:   test:
424:     suffix: 2
425:     requires: performance
426:     TODO: missing output file
427:     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
428:   test:
429:     suffix: 3
430:     requires: performance
431:     TODO: missing output file
432:     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
433:   test:
434:     suffix: 4
435:     requires: performance
436:     TODO: missing output file
437:     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
438:   test:
439:     suffix: 5
440:     requires: performance
441:     TODO: missing output file
442:     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
443:   test:
444:     suffix: 6
445:     requires: performance
446:     TODO: missing output file
447:     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
448:   test:
449:     suffix: 7
450:     requires: performance
451:     TODO: missing output file
452:     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

454:   # 2D Simplex P_1 vector tests
455:   # 2D Simplex P_2 scalar tests
456:   # 2D Simplex P_2 vector tests
457:   # 2D Simplex P_2/P_1 vector/scalar tests
458:   # 2D Quad P_1 scalar tests
459:   # 2D Quad P_1 vector tests
460:   # 2D Quad P_2 scalar tests
461:   # 2D Quad P_2 vector tests
462:   # 3D Simplex P_1 scalar tests
463:   # 3D Simplex P_1 vector tests
464:   # 3D Simplex P_2 scalar tests
465:   # 3D Simplex P_2 vector tests
466:   # 3D Hex P_1 scalar tests
467:   # 3D Hex P_1 vector tests
468:   # 3D Hex P_2 scalar tests
469:   # 3D Hex P_2 vector tests

471: TEST*/