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 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(); 81c4762a1bSJed Brown PetscFunctionReturn(0); 82c4762a1bSJed Brown } 83c4762a1bSJed Brown 84c4762a1bSJed Brown static PetscErrorCode CreateSimplex_2D(MPI_Comm comm, DM *newdm) 85c4762a1bSJed Brown { 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)); 101c4762a1bSJed Brown for (p = 0; p < 4; ++p) { 1029566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1])); 103c4762a1bSJed Brown } 104c4762a1bSJed Brown *newdm = dm; 105c4762a1bSJed Brown PetscFunctionReturn(0); 106c4762a1bSJed Brown } 107c4762a1bSJed Brown 108c4762a1bSJed Brown static PetscErrorCode CreateSimplex_3D(MPI_Comm comm, DM *newdm) 109c4762a1bSJed Brown { 110c4762a1bSJed Brown DM dm; 111c4762a1bSJed Brown PetscInt numPoints[2] = {5, 2}; 112c4762a1bSJed Brown PetscInt coneSize[23] = {4, 4, 0, 0, 0, 0, 0}; 113c4762a1bSJed Brown PetscInt cones[8] = {2, 4, 3, 5, 3, 4, 6, 5}; 114c4762a1bSJed Brown PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 115c4762a1bSJed 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}; 116c4762a1bSJed Brown PetscInt markerPoints[10] = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1}; 117c4762a1bSJed Brown PetscInt dim = 3, depth = 1, p; 118c4762a1bSJed Brown 119c4762a1bSJed Brown PetscFunctionBegin; 1209566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, &dm)); 1219566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) dm, "tetrahedral")); 1229566063dSJacob Faibussowitsch PetscCall(DMSetType(dm, DMPLEX)); 1239566063dSJacob Faibussowitsch PetscCall(DMSetDimension(dm, dim)); 1249566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 125c4762a1bSJed Brown for (p = 0; p < 5; ++p) { 1269566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1])); 127c4762a1bSJed Brown } 128c4762a1bSJed Brown *newdm = dm; 129c4762a1bSJed Brown PetscFunctionReturn(0); 130c4762a1bSJed Brown } 131c4762a1bSJed Brown 132c4762a1bSJed Brown static PetscErrorCode CreateQuad_2D(MPI_Comm comm, DM *newdm) 133c4762a1bSJed Brown { 134c4762a1bSJed Brown DM dm; 135c4762a1bSJed Brown PetscInt numPoints[2] = {6, 2}; 136c4762a1bSJed Brown PetscInt coneSize[8] = {4, 4, 0, 0, 0, 0, 0, 0}; 137c4762a1bSJed Brown PetscInt cones[8] = {2, 3, 4, 5, 3, 6, 7, 4}; 138c4762a1bSJed Brown PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 139c4762a1bSJed 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}; 140c4762a1bSJed Brown PetscInt markerPoints[12] = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1, 7, 1}; 141c4762a1bSJed Brown PetscInt dim = 2, depth = 1, p; 142c4762a1bSJed Brown 143c4762a1bSJed Brown PetscFunctionBegin; 1449566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, &dm)); 1459566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) dm, "quadrilateral")); 1469566063dSJacob Faibussowitsch PetscCall(DMSetType(dm, DMPLEX)); 1479566063dSJacob Faibussowitsch PetscCall(DMSetDimension(dm, dim)); 1489566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 149c4762a1bSJed Brown for (p = 0; p < 6; ++p) { 1509566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1])); 151c4762a1bSJed Brown } 152c4762a1bSJed Brown *newdm = dm; 153c4762a1bSJed Brown PetscFunctionReturn(0); 154c4762a1bSJed Brown } 155c4762a1bSJed Brown 156c4762a1bSJed Brown static PetscErrorCode CreateHex_3D(MPI_Comm comm, DM *newdm) 157c4762a1bSJed Brown { 158c4762a1bSJed Brown DM dm; 159c4762a1bSJed Brown PetscInt numPoints[2] = {12, 2}; 160c4762a1bSJed Brown PetscInt coneSize[14] = {8, 8, 0,0,0,0,0,0,0,0,0,0,0,0}; 161c4762a1bSJed Brown PetscInt cones[16] = {2,5,4,3,6,7,8,9, 3,4,11,10,7,12,13,8}; 162c4762a1bSJed Brown PetscInt coneOrientations[16] = {0,0,0,0,0,0,0,0, 0,0, 0, 0,0, 0, 0,0}; 163c4762a1bSJed 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, 164c4762a1bSJed 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, 165c4762a1bSJed 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}; 166c4762a1bSJed 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}; 167c4762a1bSJed Brown PetscInt dim = 3, depth = 1, p; 168c4762a1bSJed Brown 169c4762a1bSJed Brown PetscFunctionBegin; 1709566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, &dm)); 1719566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) dm, "hexahedral")); 1729566063dSJacob Faibussowitsch PetscCall(DMSetType(dm, DMPLEX)); 1739566063dSJacob Faibussowitsch PetscCall(DMSetDimension(dm, dim)); 1749566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 175c4762a1bSJed Brown for (p = 0; p < 12; ++p) { 1769566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1])); 177c4762a1bSJed Brown } 178c4762a1bSJed Brown *newdm = dm; 179c4762a1bSJed Brown PetscFunctionReturn(0); 180c4762a1bSJed Brown } 181c4762a1bSJed Brown 182c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *newdm) 183c4762a1bSJed Brown { 184c4762a1bSJed Brown PetscInt dim = user->dim; 185c4762a1bSJed Brown PetscBool cellSimplex = user->cellSimplex; 186c4762a1bSJed Brown 187c4762a1bSJed Brown PetscFunctionBegin; 188c4762a1bSJed Brown switch (dim) { 189c4762a1bSJed Brown case 2: 190c4762a1bSJed Brown if (cellSimplex) { 1919566063dSJacob Faibussowitsch PetscCall(CreateSimplex_2D(comm, newdm)); 192c4762a1bSJed Brown } else { 1939566063dSJacob Faibussowitsch PetscCall(CreateQuad_2D(comm, newdm)); 194c4762a1bSJed Brown } 195c4762a1bSJed Brown break; 196c4762a1bSJed Brown case 3: 197c4762a1bSJed Brown if (cellSimplex) { 1989566063dSJacob Faibussowitsch PetscCall(CreateSimplex_3D(comm, newdm)); 199c4762a1bSJed Brown } else { 2009566063dSJacob Faibussowitsch PetscCall(CreateHex_3D(comm, newdm)); 201c4762a1bSJed Brown } 202c4762a1bSJed Brown break; 203c4762a1bSJed Brown default: 20463a3b9bcSJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %" PetscInt_FMT, dim); 205c4762a1bSJed Brown } 206c4762a1bSJed Brown if (user->refinementLimit > 0.0) { 207c4762a1bSJed Brown DM rdm; 208c4762a1bSJed Brown const char *name; 209c4762a1bSJed Brown 2109566063dSJacob Faibussowitsch PetscCall(DMPlexSetRefinementUniform(*newdm, PETSC_FALSE)); 2119566063dSJacob Faibussowitsch PetscCall(DMPlexSetRefinementLimit(*newdm, user->refinementLimit)); 2129566063dSJacob Faibussowitsch PetscCall(DMRefine(*newdm, PETSC_COMM_SELF, &rdm)); 2139566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject) *newdm, &name)); 2149566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) rdm, name)); 2159566063dSJacob Faibussowitsch PetscCall(DMDestroy(newdm)); 216c4762a1bSJed Brown *newdm = rdm; 217c4762a1bSJed Brown } 218c4762a1bSJed Brown if (user->interpolate) { 219c4762a1bSJed Brown DM idm; 220c4762a1bSJed Brown 2219566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*newdm, &idm)); 2229566063dSJacob Faibussowitsch PetscCall(DMDestroy(newdm)); 223c4762a1bSJed Brown *newdm = idm; 224c4762a1bSJed Brown } 2259566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*newdm)); 226c4762a1bSJed Brown PetscFunctionReturn(0); 227c4762a1bSJed Brown } 228c4762a1bSJed Brown 229c4762a1bSJed Brown static PetscErrorCode TestCone(DM dm, AppCtx *user) 230c4762a1bSJed Brown { 231c4762a1bSJed Brown PetscInt numRuns, cStart, cEnd, c, i; 232c4762a1bSJed Brown PetscReal maxTimePerRun = user->maxConeTime; 233c4762a1bSJed Brown PetscLogStage stage; 234c4762a1bSJed Brown PetscLogEvent event; 235c4762a1bSJed Brown PetscEventPerfInfo eventInfo; 236c4762a1bSJed Brown MPI_Comm comm; 237c4762a1bSJed Brown PetscMPIInt rank; 238c4762a1bSJed Brown 239c4762a1bSJed Brown PetscFunctionBegin; 2409566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 2419566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 2429566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("DMPlex Cone Test", &stage)); 2439566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("Cone", PETSC_OBJECT_CLASSID, &event)); 2449566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage)); 2459566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 2469566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(event,0,0,0,0)); 247c4762a1bSJed Brown for (i = 0; i < user->iterations; ++i) { 248c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 249c4762a1bSJed Brown const PetscInt *cone; 250c4762a1bSJed Brown 2519566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, c, &cone)); 252c4762a1bSJed Brown } 253c4762a1bSJed Brown } 2549566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(event,0,0,0,0)); 2559566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 256c4762a1bSJed Brown 2579566063dSJacob Faibussowitsch PetscCall(PetscLogEventGetPerfInfo(stage, event, &eventInfo)); 258c4762a1bSJed Brown numRuns = (cEnd-cStart) * user->iterations; 25963a3b9bcSJacob Faibussowitsch PetscCheck(eventInfo.count == 1,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event calls %d should be 1", eventInfo.count); 26063a3b9bcSJacob Faibussowitsch PetscCheck((PetscInt) eventInfo.flops == 0,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event flops %" PetscInt_FMT " should be 0", (PetscInt) eventInfo.flops); 261c4762a1bSJed Brown if (user->printTimes) { 26263a3b9bcSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(comm, "[%d] Cones: %" PetscInt_FMT " Total time: %.3es Average time per cone: %.3es\n", rank, numRuns, eventInfo.time, eventInfo.time/numRuns)); 2639566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT)); 264c4762a1bSJed Brown } else if (eventInfo.time > maxTimePerRun * numRuns) { 26563a3b9bcSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(comm, "[%d] Cones: %" PetscInt_FMT " Average time per cone: %gs standard: %gs\n", rank, numRuns, (double)(eventInfo.time/numRuns), (double)maxTimePerRun)); 2669566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT)); 26763a3b9bcSJacob Faibussowitsch PetscCheck(!user->errors,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Average time for cone %g > standard %g", (double)(eventInfo.time/numRuns), (double)maxTimePerRun); 268c4762a1bSJed Brown } 269c4762a1bSJed Brown PetscFunctionReturn(0); 270c4762a1bSJed Brown } 271c4762a1bSJed Brown 272c4762a1bSJed Brown static PetscErrorCode TestTransitiveClosure(DM dm, AppCtx *user) 273c4762a1bSJed Brown { 274c4762a1bSJed Brown PetscInt numRuns, cStart, cEnd, c, i; 275c4762a1bSJed Brown PetscReal maxTimePerRun = user->maxClosureTime; 276c4762a1bSJed Brown PetscLogStage stage; 277c4762a1bSJed Brown PetscLogEvent event; 278c4762a1bSJed Brown PetscEventPerfInfo eventInfo; 279c4762a1bSJed Brown MPI_Comm comm; 280c4762a1bSJed Brown PetscMPIInt rank; 281c4762a1bSJed Brown 282c4762a1bSJed Brown PetscFunctionBegin; 2839566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 2849566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 2859566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("DMPlex Transitive Closure Test", &stage)); 2869566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("TransitiveClosure", PETSC_OBJECT_CLASSID, &event)); 2879566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage)); 2889566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 2899566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(event,0,0,0,0)); 290c4762a1bSJed Brown for (i = 0; i < user->iterations; ++i) { 291c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 292c4762a1bSJed Brown PetscInt *closure = NULL; 293c4762a1bSJed Brown PetscInt closureSize; 294c4762a1bSJed Brown 2959566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 2969566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 297c4762a1bSJed Brown } 298c4762a1bSJed Brown } 2999566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(event,0,0,0,0)); 3009566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 301c4762a1bSJed Brown 3029566063dSJacob Faibussowitsch PetscCall(PetscLogEventGetPerfInfo(stage, event, &eventInfo)); 303c4762a1bSJed Brown numRuns = (cEnd-cStart) * user->iterations; 30463a3b9bcSJacob Faibussowitsch PetscCheck(eventInfo.count == 1,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event calls %d should be 1", eventInfo.count); 30563a3b9bcSJacob Faibussowitsch PetscCheck((PetscInt) eventInfo.flops == 0,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event flops %" PetscInt_FMT " should be 0", (PetscInt) eventInfo.flops); 306c4762a1bSJed Brown if (user->printTimes) { 30763a3b9bcSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(comm, "[%d] Closures: %" PetscInt_FMT " Total time: %.3es Average time per cone: %.3es\n", rank, numRuns, eventInfo.time, eventInfo.time/numRuns)); 3089566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT)); 309c4762a1bSJed Brown } else if (eventInfo.time > maxTimePerRun * numRuns) { 31063a3b9bcSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(comm, "[%d] Closures: %" PetscInt_FMT " Average time per cone: %gs standard: %gs\n", rank, numRuns, (double)(eventInfo.time/numRuns), (double)maxTimePerRun)); 3119566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT)); 31263a3b9bcSJacob Faibussowitsch PetscCheck(!user->errors,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Average time for closure %g > standard %g", (double)(eventInfo.time/numRuns), (double)maxTimePerRun); 313c4762a1bSJed Brown } 314c4762a1bSJed Brown PetscFunctionReturn(0); 315c4762a1bSJed Brown } 316c4762a1bSJed Brown 317c4762a1bSJed Brown static PetscErrorCode TestVecClosure(DM dm, PetscBool useIndex, PetscBool useSpectral, AppCtx *user) 318c4762a1bSJed Brown { 319c4762a1bSJed Brown PetscSection s; 320c4762a1bSJed Brown Vec v; 321c4762a1bSJed Brown PetscInt numRuns, cStart, cEnd, c, i; 322c4762a1bSJed Brown PetscScalar tmpArray[64]; 323c4762a1bSJed Brown PetscScalar *userArray = user->reuseArray ? tmpArray : NULL; 324c4762a1bSJed Brown PetscReal maxTimePerRun = user->maxVecClosureTime; 325c4762a1bSJed Brown PetscLogStage stage; 326c4762a1bSJed Brown PetscLogEvent event; 327c4762a1bSJed Brown PetscEventPerfInfo eventInfo; 328c4762a1bSJed Brown MPI_Comm comm; 329c4762a1bSJed Brown PetscMPIInt rank; 330c4762a1bSJed Brown 331c4762a1bSJed Brown PetscFunctionBegin; 3329566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 3339566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 334c4762a1bSJed Brown if (useIndex) { 335c4762a1bSJed Brown if (useSpectral) { 3369566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("DMPlex Vector Closure with Index Test", &stage)); 3379566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("VecClosureInd", PETSC_OBJECT_CLASSID, &event)); 338c4762a1bSJed Brown } else { 3399566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("DMPlex Vector Spectral Closure with Index Test", &stage)); 3409566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("VecClosureSpecInd", PETSC_OBJECT_CLASSID, &event)); 341c4762a1bSJed Brown } 342c4762a1bSJed Brown } else { 343c4762a1bSJed Brown if (useSpectral) { 3449566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("DMPlex Vector Spectral Closure Test", &stage)); 3459566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("VecClosureSpec", PETSC_OBJECT_CLASSID, &event)); 346c4762a1bSJed Brown } else { 3479566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("DMPlex Vector Closure Test", &stage)); 3489566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("VecClosure", PETSC_OBJECT_CLASSID, &event)); 349c4762a1bSJed Brown } 350c4762a1bSJed Brown } 3519566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage)); 3529566063dSJacob Faibussowitsch PetscCall(DMSetNumFields(dm, user->numFields)); 3539566063dSJacob Faibussowitsch PetscCall(DMPlexCreateSection(dm, NULL, user->numComponents, user->numDof, 0, NULL, NULL, NULL, NULL, &s)); 3549566063dSJacob Faibussowitsch PetscCall(DMSetLocalSection(dm, s)); 3559566063dSJacob Faibussowitsch if (useIndex) PetscCall(DMPlexCreateClosureIndex(dm, s)); 3569566063dSJacob Faibussowitsch if (useSpectral) PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, s)); 3579566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&s)); 3589566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 3599566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &v)); 3609566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(event,0,0,0,0)); 361c4762a1bSJed Brown for (i = 0; i < user->iterations; ++i) { 362c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 363c4762a1bSJed Brown PetscScalar *closure = userArray; 364c4762a1bSJed Brown PetscInt closureSize = 64; 365c4762a1bSJed Brown 3669566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dm, s, v, c, &closureSize, &closure)); 3679566063dSJacob Faibussowitsch if (!user->reuseArray) PetscCall(DMPlexVecRestoreClosure(dm, s, v, c, &closureSize, &closure)); 368c4762a1bSJed Brown } 369c4762a1bSJed Brown } 3709566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(event,0,0,0,0)); 3719566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &v)); 3729566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 373c4762a1bSJed Brown 3749566063dSJacob Faibussowitsch PetscCall(PetscLogEventGetPerfInfo(stage, event, &eventInfo)); 375c4762a1bSJed Brown numRuns = (cEnd-cStart) * user->iterations; 37663a3b9bcSJacob Faibussowitsch PetscCheck(eventInfo.count == 1,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event calls %d should be 1", eventInfo.count); 37763a3b9bcSJacob Faibussowitsch PetscCheck((PetscInt) eventInfo.flops == 0,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event flops %" PetscInt_FMT " should be 0", (PetscInt) eventInfo.flops); 378c4762a1bSJed Brown if (user->printTimes || eventInfo.time > maxTimePerRun * numRuns) { 379c4762a1bSJed Brown const char *title = "VecClosures"; 380c4762a1bSJed Brown const char *titleIndex = "VecClosures with Index"; 381c4762a1bSJed Brown const char *titleSpec = "VecClosures Spectral"; 382c4762a1bSJed Brown const char *titleSpecIndex = "VecClosures Spectral with Index"; 383c4762a1bSJed Brown 384c4762a1bSJed Brown if (user->printTimes) { 38563a3b9bcSJacob Faibussowitsch 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, eventInfo.time, eventInfo.time/numRuns)); 3869566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT)); 387c4762a1bSJed Brown } else { 38863a3b9bcSJacob Faibussowitsch 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, (double)(eventInfo.time/numRuns), (double)maxTimePerRun)); 3899566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT)); 39063a3b9bcSJacob Faibussowitsch PetscCheck(!user->errors,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Average time for vector closure %g > standard %g", (double)(eventInfo.time/numRuns), (double)maxTimePerRun); 391c4762a1bSJed Brown } 392c4762a1bSJed Brown } 393c4762a1bSJed Brown PetscFunctionReturn(0); 394c4762a1bSJed Brown } 395c4762a1bSJed Brown 396c4762a1bSJed Brown static PetscErrorCode CleanupContext(AppCtx *user) 397c4762a1bSJed Brown { 398c4762a1bSJed Brown PetscFunctionBegin; 3999566063dSJacob Faibussowitsch PetscCall(PetscFree(user->numComponents)); 4009566063dSJacob Faibussowitsch PetscCall(PetscFree(user->numDof)); 401c4762a1bSJed Brown PetscFunctionReturn(0); 402c4762a1bSJed Brown } 403c4762a1bSJed Brown 404c4762a1bSJed Brown int main(int argc, char **argv) 405c4762a1bSJed Brown { 406c4762a1bSJed Brown DM dm; 407c4762a1bSJed Brown AppCtx user; 408c4762a1bSJed Brown 409*327415f7SBarry Smith PetscFunctionBeginUser; 4109566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL,help)); 4119566063dSJacob Faibussowitsch PetscCall(ProcessOptions(&user)); 4129566063dSJacob Faibussowitsch PetscCall(PetscLogDefaultBegin()); 4139566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 4149566063dSJacob Faibussowitsch PetscCall(TestCone(dm, &user)); 4159566063dSJacob Faibussowitsch PetscCall(TestTransitiveClosure(dm, &user)); 4169566063dSJacob Faibussowitsch PetscCall(TestVecClosure(dm, PETSC_FALSE, PETSC_FALSE, &user)); 4179566063dSJacob Faibussowitsch PetscCall(TestVecClosure(dm, PETSC_TRUE, PETSC_FALSE, &user)); 418c4762a1bSJed Brown if (!user.cellSimplex && user.spectral) { 4199566063dSJacob Faibussowitsch PetscCall(TestVecClosure(dm, PETSC_FALSE, PETSC_TRUE, &user)); 4209566063dSJacob Faibussowitsch PetscCall(TestVecClosure(dm, PETSC_TRUE, PETSC_TRUE, &user)); 421c4762a1bSJed Brown } 4229566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 4239566063dSJacob Faibussowitsch PetscCall(CleanupContext(&user)); 4249566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 425b122ec5aSJacob Faibussowitsch return 0; 426c4762a1bSJed Brown } 427c4762a1bSJed Brown 428c4762a1bSJed Brown /*TEST 429c4762a1bSJed Brown 430c4762a1bSJed Brown build: 431dfd57a17SPierre Jolivet requires: defined(PETSC_USE_LOG) 432c4762a1bSJed Brown 433c4762a1bSJed Brown # 2D Simplex P_1 scalar tests 434c4762a1bSJed Brown testset: 435c4762a1bSJed Brown args: -num_dof 1,0,0 -iterations 2 -print_times 436c4762a1bSJed Brown test: 437c4762a1bSJed Brown suffix: correctness_0 438c4762a1bSJed Brown test: 439c4762a1bSJed Brown suffix: correctness_1 440c4762a1bSJed Brown args: -interpolate -dm_refine 2 441c4762a1bSJed Brown test: 442c4762a1bSJed Brown suffix: correctness_2 443c4762a1bSJed Brown requires: triangle 444c4762a1bSJed Brown args: -interpolate -refinement_limit 1.0e-5 445c4762a1bSJed Brown test: 446c4762a1bSJed Brown suffix: 0 447c4762a1bSJed Brown TODO: Only for performance testing 448c4762a1bSJed 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 449c4762a1bSJed Brown test: 450c4762a1bSJed Brown suffix: 1 451c4762a1bSJed Brown requires: triangle 452c4762a1bSJed Brown TODO: Only for performance testing 453c4762a1bSJed 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 454c4762a1bSJed Brown test: 455c4762a1bSJed Brown suffix: 2 456c4762a1bSJed Brown TODO: Only for performance testing 457c4762a1bSJed 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 458c4762a1bSJed Brown test: 459c4762a1bSJed Brown suffix: 3 460c4762a1bSJed Brown requires: triangle 461c4762a1bSJed Brown TODO: Only for performance testing 462c4762a1bSJed 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 463c4762a1bSJed Brown test: 464c4762a1bSJed Brown suffix: 4 465c4762a1bSJed Brown TODO: Only for performance testing 466c4762a1bSJed 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 467c4762a1bSJed Brown test: 468c4762a1bSJed Brown suffix: 5 469c4762a1bSJed Brown requires: triangle 470c4762a1bSJed Brown TODO: Only for performance testing 471c4762a1bSJed 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 472c4762a1bSJed Brown test: 473c4762a1bSJed Brown suffix: 6 474c4762a1bSJed Brown TODO: Only for performance testing 475c4762a1bSJed 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 476c4762a1bSJed Brown test: 477c4762a1bSJed Brown suffix: 7 478c4762a1bSJed Brown requires: triangle 479c4762a1bSJed Brown TODO: Only for performance testing 480c4762a1bSJed 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 481c4762a1bSJed Brown 482c4762a1bSJed Brown # 2D Simplex P_1 vector tests 483c4762a1bSJed Brown # 2D Simplex P_2 scalar tests 484c4762a1bSJed Brown # 2D Simplex P_2 vector tests 485c4762a1bSJed Brown # 2D Simplex P_2/P_1 vector/scalar tests 486c4762a1bSJed Brown # 2D Quad P_1 scalar tests 487c4762a1bSJed Brown # 2D Quad P_1 vector tests 488c4762a1bSJed Brown # 2D Quad P_2 scalar tests 489c4762a1bSJed Brown # 2D Quad P_2 vector tests 490c4762a1bSJed Brown # 3D Simplex P_1 scalar tests 491c4762a1bSJed Brown # 3D Simplex P_1 vector tests 492c4762a1bSJed Brown # 3D Simplex P_2 scalar tests 493c4762a1bSJed Brown # 3D Simplex P_2 vector tests 494c4762a1bSJed Brown # 3D Hex P_1 scalar tests 495c4762a1bSJed Brown # 3D Hex P_1 vector tests 496c4762a1bSJed Brown # 3D Hex P_2 scalar tests 497c4762a1bSJed Brown # 3D Hex P_2 vector tests 498c4762a1bSJed Brown 499c4762a1bSJed Brown TEST*/ 500