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); 485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex9.c", options->dim, &options->dim, NULL,1,3)); 495f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-cellSimplex", "Flag for simplices", "ex9.c", options->cellSimplex, &options->cellSimplex, NULL)); 505f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-spectral", "Flag for spectral element layout", "ex9.c", options->spectral, &options->spectral, NULL)); 515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-interpolate", "Flag for mesh interpolation", "ex9.c", options->interpolate, &options->interpolate, NULL)); 525f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-refinement_limit", "The maximum volume of a refined cell", "ex9.c", options->refinementLimit, &options->refinementLimit, NULL)); 535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBoundedInt("-num_fields", "The number of section fields", "ex9.c", options->numFields, &options->numFields, NULL, 0)); 54c4762a1bSJed Brown if (options->numFields) { 55c4762a1bSJed Brown len = options->numFields; 565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(len, &options->numComponents)); 575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsIntArray("-num_components", "The number of components per field", "ex9.c", options->numComponents, &len, &flg)); 582c71b3e2SJacob Faibussowitsch PetscCheckFalse(flg && (len != options->numFields),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); 615f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(len, &options->numDof)); 625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsIntArray("-num_dof", "The dof signature for the section", "ex9.c", options->numDof, &len, &flg)); 632c71b3e2SJacob Faibussowitsch PetscCheckFalse(flg && (len != (options->dim+1) * PetscMax(1, options->numFields)),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 745f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-reuse_array", "Pass in user allocated array to VecGetClosure()", "ex9.c", options->reuseArray, &options->reuseArray, NULL)); 755f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-errors", "Treat failures as errors", "ex9.c", options->errors, &options->errors, NULL)); 765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBoundedInt("-iterations", "The number of iterations for a query", "ex9.c", options->iterations, &options->iterations, NULL,0)); 775f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-max_cone_time", "The maximum time per run for DMPlexGetCone()", "ex9.c", options->maxConeTime, &options->maxConeTime, NULL)); 785f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-max_closure_time", "The maximum time per run for DMPlexGetTransitiveClosure()", "ex9.c", options->maxClosureTime, &options->maxClosureTime, NULL)); 795f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-max_vec_closure_time", "The maximum time per run for DMPlexVecGetClosure()", "ex9.c", options->maxVecClosureTime, &options->maxVecClosureTime, NULL)); 805f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-print_times", "Print total times, do not check limits", "ex9.c", options->printTimes, &options->printTimes, NULL)); 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 96c4762a1bSJed Brown PetscFunctionBegin; 975f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(comm, &dm)); 985f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) dm, "triangular")); 995f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(dm, DMPLEX)); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetDimension(dm, dim)); 1015f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 102c4762a1bSJed Brown for (p = 0; p < 4; ++p) { 1035f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1])); 104c4762a1bSJed Brown } 105c4762a1bSJed Brown *newdm = dm; 106c4762a1bSJed Brown PetscFunctionReturn(0); 107c4762a1bSJed Brown } 108c4762a1bSJed Brown 109c4762a1bSJed Brown static PetscErrorCode CreateSimplex_3D(MPI_Comm comm, DM *newdm) 110c4762a1bSJed Brown { 111c4762a1bSJed Brown DM dm; 112c4762a1bSJed Brown PetscInt numPoints[2] = {5, 2}; 113c4762a1bSJed Brown PetscInt coneSize[23] = {4, 4, 0, 0, 0, 0, 0}; 114c4762a1bSJed Brown PetscInt cones[8] = {2, 4, 3, 5, 3, 4, 6, 5}; 115c4762a1bSJed Brown PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 116c4762a1bSJed 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}; 117c4762a1bSJed Brown PetscInt markerPoints[10] = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1}; 118c4762a1bSJed Brown PetscInt dim = 3, depth = 1, p; 119c4762a1bSJed Brown 120c4762a1bSJed Brown PetscFunctionBegin; 1215f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(comm, &dm)); 1225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) dm, "tetrahedral")); 1235f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(dm, DMPLEX)); 1245f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetDimension(dm, dim)); 1255f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 126c4762a1bSJed Brown for (p = 0; p < 5; ++p) { 1275f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1])); 128c4762a1bSJed Brown } 129c4762a1bSJed Brown *newdm = dm; 130c4762a1bSJed Brown PetscFunctionReturn(0); 131c4762a1bSJed Brown } 132c4762a1bSJed Brown 133c4762a1bSJed Brown static PetscErrorCode CreateQuad_2D(MPI_Comm comm, DM *newdm) 134c4762a1bSJed Brown { 135c4762a1bSJed Brown DM dm; 136c4762a1bSJed Brown PetscInt numPoints[2] = {6, 2}; 137c4762a1bSJed Brown PetscInt coneSize[8] = {4, 4, 0, 0, 0, 0, 0, 0}; 138c4762a1bSJed Brown PetscInt cones[8] = {2, 3, 4, 5, 3, 6, 7, 4}; 139c4762a1bSJed Brown PetscInt coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0}; 140c4762a1bSJed 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}; 141c4762a1bSJed Brown PetscInt markerPoints[12] = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1, 7, 1}; 142c4762a1bSJed Brown PetscInt dim = 2, depth = 1, p; 143c4762a1bSJed Brown 144c4762a1bSJed Brown PetscFunctionBegin; 1455f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(comm, &dm)); 1465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) dm, "quadrilateral")); 1475f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(dm, DMPLEX)); 1485f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetDimension(dm, dim)); 1495f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 150c4762a1bSJed Brown for (p = 0; p < 6; ++p) { 1515f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1])); 152c4762a1bSJed Brown } 153c4762a1bSJed Brown *newdm = dm; 154c4762a1bSJed Brown PetscFunctionReturn(0); 155c4762a1bSJed Brown } 156c4762a1bSJed Brown 157c4762a1bSJed Brown static PetscErrorCode CreateHex_3D(MPI_Comm comm, DM *newdm) 158c4762a1bSJed Brown { 159c4762a1bSJed Brown DM dm; 160c4762a1bSJed Brown PetscInt numPoints[2] = {12, 2}; 161c4762a1bSJed Brown PetscInt coneSize[14] = {8, 8, 0,0,0,0,0,0,0,0,0,0,0,0}; 162c4762a1bSJed Brown PetscInt cones[16] = {2,5,4,3,6,7,8,9, 3,4,11,10,7,12,13,8}; 163c4762a1bSJed Brown PetscInt coneOrientations[16] = {0,0,0,0,0,0,0,0, 0,0, 0, 0,0, 0, 0,0}; 164c4762a1bSJed 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, 165c4762a1bSJed 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, 166c4762a1bSJed 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}; 167c4762a1bSJed 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}; 168c4762a1bSJed Brown PetscInt dim = 3, depth = 1, p; 169c4762a1bSJed Brown 170c4762a1bSJed Brown PetscFunctionBegin; 1715f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(comm, &dm)); 1725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) dm, "hexahedral")); 1735f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(dm, DMPLEX)); 1745f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetDimension(dm, dim)); 1755f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords)); 176c4762a1bSJed Brown for (p = 0; p < 12; ++p) { 1775f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1])); 178c4762a1bSJed Brown } 179c4762a1bSJed Brown *newdm = dm; 180c4762a1bSJed Brown PetscFunctionReturn(0); 181c4762a1bSJed Brown } 182c4762a1bSJed Brown 183c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *newdm) 184c4762a1bSJed Brown { 185c4762a1bSJed Brown PetscInt dim = user->dim; 186c4762a1bSJed Brown PetscBool cellSimplex = user->cellSimplex; 187c4762a1bSJed Brown 188c4762a1bSJed Brown PetscFunctionBegin; 189c4762a1bSJed Brown switch (dim) { 190c4762a1bSJed Brown case 2: 191c4762a1bSJed Brown if (cellSimplex) { 1925f80ce2aSJacob Faibussowitsch CHKERRQ(CreateSimplex_2D(comm, newdm)); 193c4762a1bSJed Brown } else { 1945f80ce2aSJacob Faibussowitsch CHKERRQ(CreateQuad_2D(comm, newdm)); 195c4762a1bSJed Brown } 196c4762a1bSJed Brown break; 197c4762a1bSJed Brown case 3: 198c4762a1bSJed Brown if (cellSimplex) { 1995f80ce2aSJacob Faibussowitsch CHKERRQ(CreateSimplex_3D(comm, newdm)); 200c4762a1bSJed Brown } else { 2015f80ce2aSJacob Faibussowitsch CHKERRQ(CreateHex_3D(comm, newdm)); 202c4762a1bSJed Brown } 203c4762a1bSJed Brown break; 204c4762a1bSJed Brown default: 20598921bdaSJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %d", dim); 206c4762a1bSJed Brown } 207c4762a1bSJed Brown if (user->refinementLimit > 0.0) { 208c4762a1bSJed Brown DM rdm; 209c4762a1bSJed Brown const char *name; 210c4762a1bSJed Brown 2115f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexSetRefinementUniform(*newdm, PETSC_FALSE)); 2125f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexSetRefinementLimit(*newdm, user->refinementLimit)); 2135f80ce2aSJacob Faibussowitsch CHKERRQ(DMRefine(*newdm, PETSC_COMM_SELF, &rdm)); 2145f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetName((PetscObject) *newdm, &name)); 2155f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) rdm, name)); 2165f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(newdm)); 217c4762a1bSJed Brown *newdm = rdm; 218c4762a1bSJed Brown } 219c4762a1bSJed Brown if (user->interpolate) { 220c4762a1bSJed Brown DM idm; 221c4762a1bSJed Brown 2225f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexInterpolate(*newdm, &idm)); 2235f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(newdm)); 224c4762a1bSJed Brown *newdm = idm; 225c4762a1bSJed Brown } 2265f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(*newdm)); 227c4762a1bSJed Brown PetscFunctionReturn(0); 228c4762a1bSJed Brown } 229c4762a1bSJed Brown 230c4762a1bSJed Brown static PetscErrorCode TestCone(DM dm, AppCtx *user) 231c4762a1bSJed Brown { 232c4762a1bSJed Brown PetscInt numRuns, cStart, cEnd, c, i; 233c4762a1bSJed Brown PetscReal maxTimePerRun = user->maxConeTime; 234c4762a1bSJed Brown PetscLogStage stage; 235c4762a1bSJed Brown PetscLogEvent event; 236c4762a1bSJed Brown PetscEventPerfInfo eventInfo; 237c4762a1bSJed Brown MPI_Comm comm; 238c4762a1bSJed Brown PetscMPIInt rank; 239c4762a1bSJed Brown 240c4762a1bSJed Brown PetscFunctionBegin; 2415f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject)dm, &comm)); 2425f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(comm, &rank)); 2435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStageRegister("DMPlex Cone Test", &stage)); 2445f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventRegister("Cone", PETSC_OBJECT_CLASSID, &event)); 2455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePush(stage)); 2465f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 2475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventBegin(event,0,0,0,0)); 248c4762a1bSJed Brown for (i = 0; i < user->iterations; ++i) { 249c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 250c4762a1bSJed Brown const PetscInt *cone; 251c4762a1bSJed Brown 2525f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetCone(dm, c, &cone)); 253c4762a1bSJed Brown } 254c4762a1bSJed Brown } 2555f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventEnd(event,0,0,0,0)); 2565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePop()); 257c4762a1bSJed Brown 2585f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventGetPerfInfo(stage, event, &eventInfo)); 259c4762a1bSJed Brown numRuns = (cEnd-cStart) * user->iterations; 2602c71b3e2SJacob Faibussowitsch PetscCheckFalse(eventInfo.count != 1,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event calls %d should be %d", eventInfo.count, 1); 2612c71b3e2SJacob Faibussowitsch PetscCheckFalse((PetscInt) eventInfo.flops != 0,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event flops %d should be %d", (PetscInt) eventInfo.flops, 0); 262c4762a1bSJed Brown if (user->printTimes) { 2635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(comm, "[%d] Cones: %d Total time: %.3es Average time per cone: %.3es\n", rank, numRuns, eventInfo.time, eventInfo.time/numRuns)); 2645f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedFlush(comm, PETSC_STDOUT)); 265c4762a1bSJed Brown } else if (eventInfo.time > maxTimePerRun * numRuns) { 2665f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(comm, "[%d] Cones: %d Average time per cone: %gs standard: %gs\n", rank, numRuns, eventInfo.time/numRuns, maxTimePerRun)); 2675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedFlush(comm, PETSC_STDOUT)); 26828b400f6SJacob Faibussowitsch PetscCheck(!user->errors,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Average time for cone %g > standard %g", eventInfo.time/numRuns, maxTimePerRun); 269c4762a1bSJed Brown } 270c4762a1bSJed Brown PetscFunctionReturn(0); 271c4762a1bSJed Brown } 272c4762a1bSJed Brown 273c4762a1bSJed Brown static PetscErrorCode TestTransitiveClosure(DM dm, AppCtx *user) 274c4762a1bSJed Brown { 275c4762a1bSJed Brown PetscInt numRuns, cStart, cEnd, c, i; 276c4762a1bSJed Brown PetscReal maxTimePerRun = user->maxClosureTime; 277c4762a1bSJed Brown PetscLogStage stage; 278c4762a1bSJed Brown PetscLogEvent event; 279c4762a1bSJed Brown PetscEventPerfInfo eventInfo; 280c4762a1bSJed Brown MPI_Comm comm; 281c4762a1bSJed Brown PetscMPIInt rank; 282c4762a1bSJed Brown 283c4762a1bSJed Brown PetscFunctionBegin; 2845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject)dm, &comm)); 2855f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(comm, &rank)); 2865f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStageRegister("DMPlex Transitive Closure Test", &stage)); 2875f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventRegister("TransitiveClosure", PETSC_OBJECT_CLASSID, &event)); 2885f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePush(stage)); 2895f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 2905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventBegin(event,0,0,0,0)); 291c4762a1bSJed Brown for (i = 0; i < user->iterations; ++i) { 292c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 293c4762a1bSJed Brown PetscInt *closure = NULL; 294c4762a1bSJed Brown PetscInt closureSize; 295c4762a1bSJed Brown 2965f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 2975f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 298c4762a1bSJed Brown } 299c4762a1bSJed Brown } 3005f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventEnd(event,0,0,0,0)); 3015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePop()); 302c4762a1bSJed Brown 3035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventGetPerfInfo(stage, event, &eventInfo)); 304c4762a1bSJed Brown numRuns = (cEnd-cStart) * user->iterations; 3052c71b3e2SJacob Faibussowitsch PetscCheckFalse(eventInfo.count != 1,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event calls %d should be %d", eventInfo.count, 1); 3062c71b3e2SJacob Faibussowitsch PetscCheckFalse((PetscInt) eventInfo.flops != 0,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event flops %d should be %d", (PetscInt) eventInfo.flops, 0); 307c4762a1bSJed Brown if (user->printTimes) { 3085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(comm, "[%d] Closures: %d Total time: %.3es Average time per cone: %.3es\n", rank, numRuns, eventInfo.time, eventInfo.time/numRuns)); 3095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedFlush(comm, PETSC_STDOUT)); 310c4762a1bSJed Brown } else if (eventInfo.time > maxTimePerRun * numRuns) { 3115f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(comm, "[%d] Closures: %d Average time per cone: %gs standard: %gs\n", rank, numRuns, eventInfo.time/numRuns, maxTimePerRun)); 3125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedFlush(comm, PETSC_STDOUT)); 31328b400f6SJacob Faibussowitsch PetscCheck(!user->errors,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Average time for closure %g > standard %g", eventInfo.time/numRuns, maxTimePerRun); 314c4762a1bSJed Brown } 315c4762a1bSJed Brown PetscFunctionReturn(0); 316c4762a1bSJed Brown } 317c4762a1bSJed Brown 318c4762a1bSJed Brown static PetscErrorCode TestVecClosure(DM dm, PetscBool useIndex, PetscBool useSpectral, AppCtx *user) 319c4762a1bSJed Brown { 320c4762a1bSJed Brown PetscSection s; 321c4762a1bSJed Brown Vec v; 322c4762a1bSJed Brown PetscInt numRuns, cStart, cEnd, c, i; 323c4762a1bSJed Brown PetscScalar tmpArray[64]; 324c4762a1bSJed Brown PetscScalar *userArray = user->reuseArray ? tmpArray : NULL; 325c4762a1bSJed Brown PetscReal maxTimePerRun = user->maxVecClosureTime; 326c4762a1bSJed Brown PetscLogStage stage; 327c4762a1bSJed Brown PetscLogEvent event; 328c4762a1bSJed Brown PetscEventPerfInfo eventInfo; 329c4762a1bSJed Brown MPI_Comm comm; 330c4762a1bSJed Brown PetscMPIInt rank; 331c4762a1bSJed Brown 332c4762a1bSJed Brown PetscFunctionBegin; 3335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject)dm, &comm)); 3345f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(comm, &rank)); 335c4762a1bSJed Brown if (useIndex) { 336c4762a1bSJed Brown if (useSpectral) { 3375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStageRegister("DMPlex Vector Closure with Index Test", &stage)); 3385f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventRegister("VecClosureInd", PETSC_OBJECT_CLASSID, &event)); 339c4762a1bSJed Brown } else { 3405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStageRegister("DMPlex Vector Spectral Closure with Index Test", &stage)); 3415f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventRegister("VecClosureSpecInd", PETSC_OBJECT_CLASSID, &event)); 342c4762a1bSJed Brown } 343c4762a1bSJed Brown } else { 344c4762a1bSJed Brown if (useSpectral) { 3455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStageRegister("DMPlex Vector Spectral Closure Test", &stage)); 3465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventRegister("VecClosureSpec", PETSC_OBJECT_CLASSID, &event)); 347c4762a1bSJed Brown } else { 3485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStageRegister("DMPlex Vector Closure Test", &stage)); 3495f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventRegister("VecClosure", PETSC_OBJECT_CLASSID, &event)); 350c4762a1bSJed Brown } 351c4762a1bSJed Brown } 3525f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePush(stage)); 3535f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetNumFields(dm, user->numFields)); 3545f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateSection(dm, NULL, user->numComponents, user->numDof, 0, NULL, NULL, NULL, NULL, &s)); 3555f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetLocalSection(dm, s)); 3565f80ce2aSJacob Faibussowitsch if (useIndex) CHKERRQ(DMPlexCreateClosureIndex(dm, s)); 3575f80ce2aSJacob Faibussowitsch if (useSpectral) CHKERRQ(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, s)); 3585f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionDestroy(&s)); 3595f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 3605f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(dm, &v)); 3615f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventBegin(event,0,0,0,0)); 362c4762a1bSJed Brown for (i = 0; i < user->iterations; ++i) { 363c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 364c4762a1bSJed Brown PetscScalar *closure = userArray; 365c4762a1bSJed Brown PetscInt closureSize = 64; 366c4762a1bSJed Brown 3675f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexVecGetClosure(dm, s, v, c, &closureSize, &closure)); 3685f80ce2aSJacob Faibussowitsch if (!user->reuseArray) CHKERRQ(DMPlexVecRestoreClosure(dm, s, v, c, &closureSize, &closure)); 369c4762a1bSJed Brown } 370c4762a1bSJed Brown } 3715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventEnd(event,0,0,0,0)); 3725f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(dm, &v)); 3735f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStagePop()); 374c4762a1bSJed Brown 3755f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogEventGetPerfInfo(stage, event, &eventInfo)); 376c4762a1bSJed Brown numRuns = (cEnd-cStart) * user->iterations; 3772c71b3e2SJacob Faibussowitsch PetscCheckFalse(eventInfo.count != 1,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event calls %d should be %d", eventInfo.count, 1); 3782c71b3e2SJacob Faibussowitsch PetscCheckFalse((PetscInt) eventInfo.flops != 0,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of event flops %d should be %d", (PetscInt) eventInfo.flops, 0); 379c4762a1bSJed Brown if (user->printTimes || eventInfo.time > maxTimePerRun * numRuns) { 380c4762a1bSJed Brown const char *title = "VecClosures"; 381c4762a1bSJed Brown const char *titleIndex = "VecClosures with Index"; 382c4762a1bSJed Brown const char *titleSpec = "VecClosures Spectral"; 383c4762a1bSJed Brown const char *titleSpecIndex = "VecClosures Spectral with Index"; 384c4762a1bSJed Brown 385c4762a1bSJed Brown if (user->printTimes) { 3865f80ce2aSJacob Faibussowitsch CHKERRQ(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)); 3875f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedFlush(comm, PETSC_STDOUT)); 388c4762a1bSJed Brown } else { 3895f80ce2aSJacob Faibussowitsch CHKERRQ(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)); 3905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedFlush(comm, PETSC_STDOUT)); 39128b400f6SJacob Faibussowitsch PetscCheck(!user->errors,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Average time for vector closure %g > standard %g", eventInfo.time/numRuns, maxTimePerRun); 392c4762a1bSJed Brown } 393c4762a1bSJed Brown } 394c4762a1bSJed Brown PetscFunctionReturn(0); 395c4762a1bSJed Brown } 396c4762a1bSJed Brown 397c4762a1bSJed Brown static PetscErrorCode CleanupContext(AppCtx *user) 398c4762a1bSJed Brown { 399c4762a1bSJed Brown PetscFunctionBegin; 4005f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(user->numComponents)); 4015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(user->numDof)); 402c4762a1bSJed Brown PetscFunctionReturn(0); 403c4762a1bSJed Brown } 404c4762a1bSJed Brown 405c4762a1bSJed Brown int main(int argc, char **argv) 406c4762a1bSJed Brown { 407c4762a1bSJed Brown DM dm; 408c4762a1bSJed Brown AppCtx user; 409c4762a1bSJed Brown 410*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc, &argv, NULL,help)); 4115f80ce2aSJacob Faibussowitsch CHKERRQ(ProcessOptions(&user)); 4125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogDefaultBegin()); 4135f80ce2aSJacob Faibussowitsch CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 4145f80ce2aSJacob Faibussowitsch CHKERRQ(TestCone(dm, &user)); 4155f80ce2aSJacob Faibussowitsch CHKERRQ(TestTransitiveClosure(dm, &user)); 4165f80ce2aSJacob Faibussowitsch CHKERRQ(TestVecClosure(dm, PETSC_FALSE, PETSC_FALSE, &user)); 4175f80ce2aSJacob Faibussowitsch CHKERRQ(TestVecClosure(dm, PETSC_TRUE, PETSC_FALSE, &user)); 418c4762a1bSJed Brown if (!user.cellSimplex && user.spectral) { 4195f80ce2aSJacob Faibussowitsch CHKERRQ(TestVecClosure(dm, PETSC_FALSE, PETSC_TRUE, &user)); 4205f80ce2aSJacob Faibussowitsch CHKERRQ(TestVecClosure(dm, PETSC_TRUE, PETSC_TRUE, &user)); 421c4762a1bSJed Brown } 4225f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 4235f80ce2aSJacob Faibussowitsch CHKERRQ(CleanupContext(&user)); 424*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 425*b122ec5aSJacob 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