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