1c4762a1bSJed Brown static char help[] = "Tests for cell geometry\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdmplex.h> 4c4762a1bSJed Brown 54f99dae5SMatthew G. Knepley typedef enum {RUN_REFERENCE, RUN_HEX_CURVED, RUN_FILE, RUN_DISPLAY} RunType; 6c4762a1bSJed Brown 7c4762a1bSJed Brown typedef struct { 8c4762a1bSJed Brown DM dm; 9c4762a1bSJed Brown RunType runType; /* Type of mesh to use */ 10c4762a1bSJed Brown PetscBool transform; /* Use random coordinate transformations */ 11c4762a1bSJed Brown /* Data for input meshes */ 12c4762a1bSJed Brown PetscReal *v0, *J, *invJ, *detJ; /* FEM data */ 13c4762a1bSJed Brown PetscReal *centroid, *normal, *vol; /* FVM data */ 14c4762a1bSJed Brown } AppCtx; 15c4762a1bSJed Brown 164f99dae5SMatthew G. Knepley static PetscErrorCode ReadMesh(MPI_Comm comm, AppCtx *user, DM *dm) 17c4762a1bSJed Brown { 18c4762a1bSJed Brown PetscFunctionBegin; 199566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 209566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 219566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 229566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(*dm, user)); 239566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) *dm, "Input Mesh")); 249566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 25c4762a1bSJed Brown PetscFunctionReturn(0); 26c4762a1bSJed Brown } 27c4762a1bSJed Brown 28c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 29c4762a1bSJed Brown { 304f99dae5SMatthew G. Knepley const char *runTypes[4] = {"reference", "hex_curved", "file", "display"}; 31c4762a1bSJed Brown PetscInt run; 32c4762a1bSJed Brown PetscErrorCode ierr; 33c4762a1bSJed Brown 34c4762a1bSJed Brown PetscFunctionBeginUser; 35c4762a1bSJed Brown options->runType = RUN_REFERENCE; 36c4762a1bSJed Brown options->transform = PETSC_FALSE; 37c4762a1bSJed Brown 389566063dSJacob Faibussowitsch ierr = PetscOptionsBegin(comm, "", "Geometry Test Options", "DMPLEX");PetscCall(ierr); 39c4762a1bSJed Brown run = options->runType; 409566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-run_type", "The run type", "ex8.c", runTypes, 3, runTypes[options->runType], &run, NULL)); 41c4762a1bSJed Brown options->runType = (RunType) run; 429566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-transform", "Use random transforms", "ex8.c", options->transform, &options->transform, NULL)); 43c4762a1bSJed Brown 44c4762a1bSJed Brown if (options->runType == RUN_FILE) { 45c4762a1bSJed Brown PetscInt dim, cStart, cEnd, numCells, n; 469bf2564aSMatt McGurn PetscBool flg, feFlg; 47c4762a1bSJed Brown 489566063dSJacob Faibussowitsch PetscCall(ReadMesh(PETSC_COMM_WORLD, options, &options->dm)); 499566063dSJacob Faibussowitsch PetscCall(DMGetDimension(options->dm, &dim)); 509566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(options->dm, 0, &cStart, &cEnd)); 51c4762a1bSJed Brown numCells = cEnd-cStart; 529566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(numCells*dim, &options->v0, numCells*dim*dim, &options->J, numCells*dim*dim, &options->invJ, numCells, &options->detJ)); 539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numCells*dim, &options->centroid)); 549566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numCells*dim, &options->normal)); 559566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numCells, &options->vol)); 56c4762a1bSJed Brown n = numCells*dim; 579566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-v0", "Input v0 for each cell", "ex8.c", options->v0, &n, &feFlg)); 58*700a43edSPierre Jolivet PetscCheck(!feFlg || n == numCells*dim, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid size of v0 %D should be %D", n, numCells*dim); 59c4762a1bSJed Brown n = numCells*dim*dim; 609566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-J", "Input Jacobian for each cell", "ex8.c", options->J, &n, &flg)); 61*700a43edSPierre Jolivet PetscCheck(!flg || n == numCells*dim*dim, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid size of J %D should be %D", n, numCells*dim*dim); 62c4762a1bSJed Brown n = numCells*dim*dim; 639566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-invJ", "Input inverse Jacobian for each cell", "ex8.c", options->invJ, &n, &flg)); 64*700a43edSPierre Jolivet PetscCheck(!flg || n == numCells*dim*dim, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid size of invJ %D should be %D", n, numCells*dim*dim); 65c4762a1bSJed Brown n = numCells; 669566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-detJ", "Input Jacobian determinant for each cell", "ex8.c", options->detJ, &n, &flg)); 67*700a43edSPierre Jolivet PetscCheck(!flg || n == numCells, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid size of detJ %D should be %D", n, numCells); 68c4762a1bSJed Brown n = numCells*dim; 694f99dae5SMatthew G. Knepley if (!feFlg) { 709566063dSJacob Faibussowitsch PetscCall(PetscFree4(options->v0, options->J, options->invJ, options->detJ)); 714f99dae5SMatthew G. Knepley options->v0 = options->J = options->invJ = options->detJ = NULL; 724f99dae5SMatthew G. Knepley } 739566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-centroid", "Input centroid for each cell", "ex8.c", options->centroid, &n, &flg)); 74*700a43edSPierre Jolivet PetscCheck(!flg || n == numCells*dim, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid size of centroid %D should be %D", n, numCells*dim); 759bf2564aSMatt McGurn if (!flg) { 769566063dSJacob Faibussowitsch PetscCall(PetscFree(options->centroid)); 779bf2564aSMatt McGurn options->centroid = NULL; 789bf2564aSMatt McGurn } 79c4762a1bSJed Brown n = numCells*dim; 809566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-normal", "Input normal for each cell", "ex8.c", options->normal, &n, &flg)); 81*700a43edSPierre Jolivet PetscCheck(!flg || n == numCells*dim, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid size of normal %D should be %D", n, numCells*dim); 829bf2564aSMatt McGurn if (!flg) { 839566063dSJacob Faibussowitsch PetscCall(PetscFree(options->normal)); 849bf2564aSMatt McGurn options->normal = NULL; 859bf2564aSMatt McGurn } 86c4762a1bSJed Brown n = numCells; 879566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-vol", "Input volume for each cell", "ex8.c", options->vol, &n, &flg)); 88*700a43edSPierre Jolivet PetscCheck(!flg || n == numCells, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid size of vol %D should be %D", n, numCells); 899bf2564aSMatt McGurn if (!flg) { 909566063dSJacob Faibussowitsch PetscCall(PetscFree(options->vol)); 919bf2564aSMatt McGurn options->vol = NULL; 924f99dae5SMatthew G. Knepley } 93c4762a1bSJed Brown } else if (options->runType == RUN_DISPLAY) { 949566063dSJacob Faibussowitsch PetscCall(ReadMesh(PETSC_COMM_WORLD, options, &options->dm)); 95c4762a1bSJed Brown } 969566063dSJacob Faibussowitsch ierr = PetscOptionsEnd();PetscCall(ierr); 97c4762a1bSJed Brown 989566063dSJacob Faibussowitsch if (options->transform) PetscCall(PetscPrintf(comm, "Using random transforms\n")); 99c4762a1bSJed Brown PetscFunctionReturn(0); 100c4762a1bSJed Brown } 101c4762a1bSJed Brown 102c4762a1bSJed Brown static PetscErrorCode ChangeCoordinates(DM dm, PetscInt spaceDim, PetscScalar vertexCoords[]) 103c4762a1bSJed Brown { 104c4762a1bSJed Brown PetscSection coordSection; 105c4762a1bSJed Brown Vec coordinates; 106c4762a1bSJed Brown PetscScalar *coords; 107c4762a1bSJed Brown PetscInt vStart, vEnd, v, d, coordSize; 108c4762a1bSJed Brown 109c4762a1bSJed Brown PetscFunctionBegin; 1109566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 1119566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm, &coordSection)); 1129566063dSJacob Faibussowitsch PetscCall(PetscSectionSetNumFields(coordSection, 1)); 1139566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldComponents(coordSection, 0, spaceDim)); 1149566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(coordSection, vStart, vEnd)); 115c4762a1bSJed Brown for (v = vStart; v < vEnd; ++v) { 1169566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(coordSection, v, spaceDim)); 1179566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, spaceDim)); 118c4762a1bSJed Brown } 1199566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(coordSection)); 1209566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize)); 1219566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates)); 1229566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) coordinates, "coordinates")); 1239566063dSJacob Faibussowitsch PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE)); 1249566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(coordinates)); 1259566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordinates, &coords)); 126c4762a1bSJed Brown for (v = vStart; v < vEnd; ++v) { 127c4762a1bSJed Brown PetscInt off; 128c4762a1bSJed Brown 1299566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(coordSection, v, &off)); 130c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) { 131c4762a1bSJed Brown coords[off+d] = vertexCoords[(v-vStart)*spaceDim+d]; 132c4762a1bSJed Brown } 133c4762a1bSJed Brown } 1349566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coordinates, &coords)); 1359566063dSJacob Faibussowitsch PetscCall(DMSetCoordinateDim(dm, spaceDim)); 1369566063dSJacob Faibussowitsch PetscCall(DMSetCoordinatesLocal(dm, coordinates)); 1379566063dSJacob Faibussowitsch PetscCall(VecDestroy(&coordinates)); 1389566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 139c4762a1bSJed Brown PetscFunctionReturn(0); 140c4762a1bSJed Brown } 141c4762a1bSJed Brown 142c4762a1bSJed Brown #define RelativeError(a,b) PetscAbs(a-b)/(1.0+PetscMax(PetscAbs(a),PetscAbs(b))) 143c4762a1bSJed Brown 144c4762a1bSJed Brown static PetscErrorCode CheckFEMGeometry(DM dm, PetscInt cell, PetscInt spaceDim, PetscReal v0Ex[], PetscReal JEx[], PetscReal invJEx[], PetscReal detJEx) 145c4762a1bSJed Brown { 146c4762a1bSJed Brown PetscReal v0[3], J[9], invJ[9], detJ; 147c4762a1bSJed Brown PetscInt d, i, j; 148c4762a1bSJed Brown 149c4762a1bSJed Brown PetscFunctionBegin; 1509566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFEM(dm, cell, NULL, v0, J, invJ, &detJ)); 151c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) { 152c4762a1bSJed Brown if (v0[d] != v0Ex[d]) { 153c4762a1bSJed Brown switch (spaceDim) { 15498921bdaSJacob Faibussowitsch case 2: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid v0 (%g, %g) != (%g, %g)", (double)v0[0], (double)v0[1], (double)v0Ex[0], (double)v0Ex[1]); 15598921bdaSJacob Faibussowitsch case 3: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid v0 (%g, %g, %g) != (%g, %g, %g)", (double)v0[0], (double)v0[1], (double)v0[2], (double)v0Ex[0], (double)v0Ex[1], (double)v0Ex[2]); 15698921bdaSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid space dimension %D", spaceDim); 157c4762a1bSJed Brown } 158c4762a1bSJed Brown } 159c4762a1bSJed Brown } 160c4762a1bSJed Brown for (i = 0; i < spaceDim; ++i) { 161c4762a1bSJed Brown for (j = 0; j < spaceDim; ++j) { 162*700a43edSPierre Jolivet PetscCheck(RelativeError(J[i*spaceDim+j], JEx[i*spaceDim+j]) < 10*PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid J[%D,%D]: %g != %g", i, j, (double)J[i*spaceDim+j], (double)JEx[i*spaceDim+j]); 163*700a43edSPierre Jolivet PetscCheck(RelativeError(invJ[i*spaceDim+j], invJEx[i*spaceDim+j]) < 10*PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid invJ[%D,%D]: %g != %g", i, j, (double)invJ[i*spaceDim+j], (double)invJEx[i*spaceDim+j]); 164c4762a1bSJed Brown } 165c4762a1bSJed Brown } 166*700a43edSPierre Jolivet PetscCheck(RelativeError(detJ, detJEx) < 10*PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid |J| = %g != %g diff %g", (double)detJ, (double)detJEx,(double)(detJ - detJEx)); 167c4762a1bSJed Brown PetscFunctionReturn(0); 168c4762a1bSJed Brown } 169c4762a1bSJed Brown 170c4762a1bSJed Brown static PetscErrorCode CheckFVMGeometry(DM dm, PetscInt cell, PetscInt spaceDim, PetscReal centroidEx[], PetscReal normalEx[], PetscReal volEx) 171c4762a1bSJed Brown { 1724f99dae5SMatthew G. Knepley PetscReal tol = PetscMax(10*PETSC_SMALL, 1e-10); 173c4762a1bSJed Brown PetscReal centroid[3], normal[3], vol; 174c4762a1bSJed Brown PetscInt d; 175c4762a1bSJed Brown 176c4762a1bSJed Brown PetscFunctionBegin; 1779566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFVM(dm, cell, volEx? &vol : NULL, centroidEx? centroid : NULL, normalEx? normal : NULL)); 178c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) { 1799bf2564aSMatt McGurn if (centroidEx) 180*700a43edSPierre Jolivet PetscCheck(RelativeError(centroid[d], centroidEx[d]) < tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cell %D, Invalid centroid[%D]: %g != %g diff %g", cell, d, (double)centroid[d], (double)centroidEx[d],(double)(centroid[d]-centroidEx[d])); 1819bf2564aSMatt McGurn if (normalEx) 182*700a43edSPierre Jolivet PetscCheck(RelativeError(normal[d], normalEx[d]) < tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cell %D, Invalid normal[%D]: %g != %g", cell, d, (double) normal[d], (double) normalEx[d]); 183c4762a1bSJed Brown } 1849bf2564aSMatt McGurn if (volEx) 185*700a43edSPierre Jolivet PetscCheck(RelativeError(volEx, vol) < tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cell %D, Invalid volume = %g != %g diff %g", cell, (double)vol, (double)volEx,(double)(vol - volEx)); 1864f99dae5SMatthew G. Knepley PetscFunctionReturn(0); 1874f99dae5SMatthew G. Knepley } 1884f99dae5SMatthew G. Knepley 1894f99dae5SMatthew G. Knepley static PetscErrorCode CheckGaussLaw(DM dm, PetscInt cell) 1904f99dae5SMatthew G. Knepley { 1914f99dae5SMatthew G. Knepley DMPolytopeType ct; 1924f99dae5SMatthew G. Knepley PetscReal tol = PetscMax(10*PETSC_SMALL, 1e-10); 1934f99dae5SMatthew G. Knepley PetscReal normal[3], integral[3] = {0., 0., 0.}, area; 1944f99dae5SMatthew G. Knepley const PetscInt *cone, *ornt; 1954f99dae5SMatthew G. Knepley PetscInt coneSize, f, dim, cdim, d; 1964f99dae5SMatthew G. Knepley 1974f99dae5SMatthew G. Knepley PetscFunctionBegin; 1989566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1999566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &cdim)); 2004f99dae5SMatthew G. Knepley if (dim != cdim) PetscFunctionReturn(0); 2019566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, cell, &ct)); 2024f99dae5SMatthew G. Knepley if (ct == DM_POLYTOPE_TRI_PRISM_TENSOR) PetscFunctionReturn(0); 2039566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, cell, &coneSize)); 2049566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, cell, &cone)); 2059566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeOrientation(dm, cell, &ornt)); 2064f99dae5SMatthew G. Knepley for (f = 0; f < coneSize; ++f) { 2079bf2564aSMatt McGurn const PetscInt sgn = dim == 1? (f == 0 ? -1 : 1) : (ornt[f] < 0 ? -1 : 1); 2084f99dae5SMatthew G. Knepley 2099566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFVM(dm, cone[f], &area, NULL, normal)); 2104f99dae5SMatthew G. Knepley for (d = 0; d < cdim; ++d) integral[d] += sgn*area*normal[d]; 2114f99dae5SMatthew G. Knepley } 212*700a43edSPierre Jolivet for (d = 0; d < cdim; ++d) PetscCheck(PetscAbsReal(integral[d]) < tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cell %D Surface integral for component %D: %g != 0. as it should be for a constant field", cell, d, (double) integral[d]); 213c4762a1bSJed Brown PetscFunctionReturn(0); 214c4762a1bSJed Brown } 215c4762a1bSJed Brown 216c4762a1bSJed Brown static PetscErrorCode CheckCell(DM dm, PetscInt cell, PetscBool transform, PetscReal v0Ex[], PetscReal JEx[], PetscReal invJEx[], PetscReal detJEx, PetscReal centroidEx[], PetscReal normalEx[], PetscReal volEx, PetscReal faceCentroidEx[], PetscReal faceNormalEx[], PetscReal faceVolEx[]) 217c4762a1bSJed Brown { 218c4762a1bSJed Brown const PetscInt *cone; 219c4762a1bSJed Brown PetscInt coneSize, c; 220c4762a1bSJed Brown PetscInt dim, depth, cdim; 221c4762a1bSJed Brown 222c4762a1bSJed Brown PetscFunctionBegin; 2239566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepth(dm, &depth)); 2249566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 2259566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &cdim)); 226c4762a1bSJed Brown if (v0Ex) { 2279566063dSJacob Faibussowitsch PetscCall(CheckFEMGeometry(dm, cell, cdim, v0Ex, JEx, invJEx, detJEx)); 228c4762a1bSJed Brown } 229c4762a1bSJed Brown if (dim == depth && centroidEx) { 2309566063dSJacob Faibussowitsch PetscCall(CheckFVMGeometry(dm, cell, cdim, centroidEx, normalEx, volEx)); 2319566063dSJacob Faibussowitsch PetscCall(CheckGaussLaw(dm, cell)); 232c4762a1bSJed Brown if (faceCentroidEx) { 2339566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, cell, &coneSize)); 2349566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, cell, &cone)); 235c4762a1bSJed Brown for (c = 0; c < coneSize; ++c) { 2369566063dSJacob Faibussowitsch PetscCall(CheckFVMGeometry(dm, cone[c], dim, &faceCentroidEx[c*dim], &faceNormalEx[c*dim], faceVolEx[c])); 237c4762a1bSJed Brown } 238c4762a1bSJed Brown } 239c4762a1bSJed Brown } 240c4762a1bSJed Brown if (transform) { 241c4762a1bSJed Brown Vec coordinates; 242c4762a1bSJed Brown PetscSection coordSection; 243c4762a1bSJed Brown PetscScalar *coords = NULL, *origCoords, *newCoords; 244c4762a1bSJed Brown PetscReal *v0ExT, *JExT, *invJExT, detJExT=0, *centroidExT, *normalExT, volExT=0; 245c4762a1bSJed Brown PetscReal *faceCentroidExT, *faceNormalExT, faceVolExT; 246c4762a1bSJed Brown PetscRandom r, ang, ang2; 247c4762a1bSJed Brown PetscInt coordSize, numCorners, t; 248c4762a1bSJed Brown 2499566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 2509566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm, &coordSection)); 2519566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords)); 2529566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(coordSize, &origCoords, coordSize, &newCoords)); 2539566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(cdim, &v0ExT, cdim*cdim, &JExT, cdim*cdim, &invJExT, cdim, ¢roidExT, cdim, &normalExT)); 2549566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(cdim, &faceCentroidExT, cdim, &faceNormalExT)); 255c4762a1bSJed Brown for (c = 0; c < coordSize; ++c) origCoords[c] = coords[c]; 2569566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords)); 257c4762a1bSJed Brown numCorners = coordSize/cdim; 258c4762a1bSJed Brown 2599566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &r)); 2609566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(r)); 2619566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(r, 0.0, 10.0)); 2629566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &ang)); 2639566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(ang)); 2649566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(ang, 0.0, 2*PETSC_PI)); 2659566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &ang2)); 2669566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(ang2)); 2679566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(ang2, 0.0, PETSC_PI)); 268c4762a1bSJed Brown for (t = 0; t < 1; ++t) { 269c4762a1bSJed Brown PetscScalar trans[3]; 270c4762a1bSJed Brown PetscReal R[9], rot[3], rotM[9]; 271c4762a1bSJed Brown PetscReal scale, phi, theta, psi = 0.0, norm; 272c4762a1bSJed Brown PetscInt d, e, f, p; 273c4762a1bSJed Brown 274c4762a1bSJed Brown for (c = 0; c < coordSize; ++c) newCoords[c] = origCoords[c]; 2759566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(r, &scale)); 2769566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(ang, &phi)); 2779566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(ang2, &theta)); 278c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 2799566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(r, &trans[d])); 280c4762a1bSJed Brown } 281c4762a1bSJed Brown switch (cdim) { 282c4762a1bSJed Brown case 2: 283c4762a1bSJed Brown R[0] = PetscCosReal(phi); R[1] = -PetscSinReal(phi); 284c4762a1bSJed Brown R[2] = PetscSinReal(phi); R[3] = PetscCosReal(phi); 285c4762a1bSJed Brown break; 286c4762a1bSJed Brown case 3: 287c4762a1bSJed Brown { 288c4762a1bSJed Brown const PetscReal ct = PetscCosReal(theta), st = PetscSinReal(theta); 289c4762a1bSJed Brown const PetscReal cp = PetscCosReal(phi), sp = PetscSinReal(phi); 290c4762a1bSJed Brown const PetscReal cs = PetscCosReal(psi), ss = PetscSinReal(psi); 291c4762a1bSJed Brown R[0] = ct*cs; R[1] = sp*st*cs - cp*ss; R[2] = sp*ss + cp*st*cs; 292c4762a1bSJed Brown R[3] = ct*ss; R[4] = cp*cs + sp*st*ss; R[5] = cp*st*ss - sp*cs; 293c4762a1bSJed Brown R[6] = -st; R[7] = sp*ct; R[8] = cp*ct; 294c4762a1bSJed Brown break; 295c4762a1bSJed Brown } 29698921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Invalid coordinate dimension %D", cdim); 297c4762a1bSJed Brown } 298c4762a1bSJed Brown if (v0Ex) { 299c4762a1bSJed Brown detJExT = detJEx; 300c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 301c4762a1bSJed Brown v0ExT[d] = v0Ex[d]; 302c4762a1bSJed Brown for (e = 0; e < cdim; ++e) { 303c4762a1bSJed Brown JExT[d*cdim+e] = JEx[d*cdim+e]; 304c4762a1bSJed Brown invJExT[d*cdim+e] = invJEx[d*cdim+e]; 305c4762a1bSJed Brown } 306c4762a1bSJed Brown } 307c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 308c4762a1bSJed Brown v0ExT[d] *= scale; 309c4762a1bSJed Brown v0ExT[d] += PetscRealPart(trans[d]); 310c4762a1bSJed Brown /* Only scale dimensions in the manifold */ 311c4762a1bSJed Brown for (e = 0; e < dim; ++e) { 312c4762a1bSJed Brown JExT[d*cdim+e] *= scale; 313c4762a1bSJed Brown invJExT[d*cdim+e] /= scale; 314c4762a1bSJed Brown } 315c4762a1bSJed Brown if (d < dim) detJExT *= scale; 316c4762a1bSJed Brown } 317c4762a1bSJed Brown /* Do scaling and translation before rotation, so that we can leave out the normal dimension for lower dimensional manifolds */ 318c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 319c4762a1bSJed Brown for (e = 0, rot[d] = 0.0; e < cdim; ++e) { 320c4762a1bSJed Brown rot[d] += R[d*cdim+e] * v0ExT[e]; 321c4762a1bSJed Brown } 322c4762a1bSJed Brown } 323c4762a1bSJed Brown for (d = 0; d < cdim; ++d) v0ExT[d] = rot[d]; 324c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 325c4762a1bSJed Brown for (e = 0; e < cdim; ++e) { 326c4762a1bSJed Brown for (f = 0, rotM[d*cdim+e] = 0.0; f < cdim; ++f) { 327c4762a1bSJed Brown rotM[d*cdim+e] += R[d*cdim+f] * JExT[f*cdim+e]; 328c4762a1bSJed Brown } 329c4762a1bSJed Brown } 330c4762a1bSJed Brown } 331c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 332c4762a1bSJed Brown for (e = 0; e < cdim; ++e) { 333c4762a1bSJed Brown JExT[d*cdim+e] = rotM[d*cdim+e]; 334c4762a1bSJed Brown } 335c4762a1bSJed Brown } 336c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 337c4762a1bSJed Brown for (e = 0; e < cdim; ++e) { 338c4762a1bSJed Brown for (f = 0, rotM[d*cdim+e] = 0.0; f < cdim; ++f) { 339c4762a1bSJed Brown rotM[d*cdim+e] += invJExT[d*cdim+f] * R[e*cdim+f]; 340c4762a1bSJed Brown } 341c4762a1bSJed Brown } 342c4762a1bSJed Brown } 343c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 344c4762a1bSJed Brown for (e = 0; e < cdim; ++e) { 345c4762a1bSJed Brown invJExT[d*cdim+e] = rotM[d*cdim+e]; 346c4762a1bSJed Brown } 347c4762a1bSJed Brown } 348c4762a1bSJed Brown } 349c4762a1bSJed Brown if (centroidEx) { 350c4762a1bSJed Brown volExT = volEx; 351c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 352c4762a1bSJed Brown centroidExT[d] = centroidEx[d]; 353c4762a1bSJed Brown normalExT[d] = normalEx[d]; 354c4762a1bSJed Brown } 355c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 356c4762a1bSJed Brown centroidExT[d] *= scale; 357c4762a1bSJed Brown centroidExT[d] += PetscRealPart(trans[d]); 358c4762a1bSJed Brown normalExT[d] /= scale; 359c4762a1bSJed Brown /* Only scale dimensions in the manifold */ 360c4762a1bSJed Brown if (d < dim) volExT *= scale; 361c4762a1bSJed Brown } 362c4762a1bSJed Brown /* Do scaling and translation before rotation, so that we can leave out the normal dimension for lower dimensional manifolds */ 363c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 364c4762a1bSJed Brown for (e = 0, rot[d] = 0.0; e < cdim; ++e) { 365c4762a1bSJed Brown rot[d] += R[d*cdim+e] * centroidExT[e]; 366c4762a1bSJed Brown } 367c4762a1bSJed Brown } 368c4762a1bSJed Brown for (d = 0; d < cdim; ++d) centroidExT[d] = rot[d]; 369c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 370c4762a1bSJed Brown for (e = 0, rot[d] = 0.0; e < cdim; ++e) { 371c4762a1bSJed Brown rot[d] += R[d*cdim+e] * normalExT[e]; 372c4762a1bSJed Brown } 373c4762a1bSJed Brown } 374c4762a1bSJed Brown for (d = 0; d < cdim; ++d) normalExT[d] = rot[d]; 375c4762a1bSJed Brown for (d = 0, norm = 0.0; d < cdim; ++d) norm += PetscSqr(normalExT[d]); 376c4762a1bSJed Brown norm = PetscSqrtReal(norm); 377*700a43edSPierre Jolivet if (norm != 0.) for (d = 0; d < cdim; ++d) normalExT[d] /= norm; 378c4762a1bSJed Brown } 379c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 380c4762a1bSJed Brown for (p = 0; p < numCorners; ++p) { 381c4762a1bSJed Brown newCoords[p*cdim+d] *= scale; 382c4762a1bSJed Brown newCoords[p*cdim+d] += trans[d]; 383c4762a1bSJed Brown } 384c4762a1bSJed Brown } 385c4762a1bSJed Brown for (p = 0; p < numCorners; ++p) { 386c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 387c4762a1bSJed Brown for (e = 0, rot[d] = 0.0; e < cdim; ++e) { 388c4762a1bSJed Brown rot[d] += R[d*cdim+e] * PetscRealPart(newCoords[p*cdim+e]); 389c4762a1bSJed Brown } 390c4762a1bSJed Brown } 391c4762a1bSJed Brown for (d = 0; d < cdim; ++d) newCoords[p*cdim+d] = rot[d]; 392c4762a1bSJed Brown } 393c4762a1bSJed Brown 3949566063dSJacob Faibussowitsch PetscCall(ChangeCoordinates(dm, cdim, newCoords)); 395c4762a1bSJed Brown if (v0Ex) { 3969566063dSJacob Faibussowitsch PetscCall(CheckFEMGeometry(dm, 0, cdim, v0ExT, JExT, invJExT, detJExT)); 397c4762a1bSJed Brown } 398c4762a1bSJed Brown if (dim == depth && centroidEx) { 3999566063dSJacob Faibussowitsch PetscCall(CheckFVMGeometry(dm, cell, cdim, centroidExT, normalExT, volExT)); 4009566063dSJacob Faibussowitsch PetscCall(CheckGaussLaw(dm, cell)); 401c4762a1bSJed Brown if (faceCentroidEx) { 4029566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, cell, &coneSize)); 4039566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, cell, &cone)); 404c4762a1bSJed Brown for (c = 0; c < coneSize; ++c) { 405c4762a1bSJed Brown PetscInt off = c*cdim; 406c4762a1bSJed Brown 407c4762a1bSJed Brown faceVolExT = faceVolEx[c]; 408c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 409c4762a1bSJed Brown faceCentroidExT[d] = faceCentroidEx[off+d]; 410c4762a1bSJed Brown faceNormalExT[d] = faceNormalEx[off+d]; 411c4762a1bSJed Brown } 412c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 413c4762a1bSJed Brown faceCentroidExT[d] *= scale; 414c4762a1bSJed Brown faceCentroidExT[d] += PetscRealPart(trans[d]); 415c4762a1bSJed Brown faceNormalExT[d] /= scale; 416c4762a1bSJed Brown /* Only scale dimensions in the manifold */ 417c4762a1bSJed Brown if (d < dim-1) { 418c4762a1bSJed Brown faceVolExT *= scale; 419c4762a1bSJed Brown } 420c4762a1bSJed Brown } 421c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 422c4762a1bSJed Brown for (e = 0, rot[d] = 0.0; e < cdim; ++e) { 423c4762a1bSJed Brown rot[d] += R[d*cdim+e] * faceCentroidExT[e]; 424c4762a1bSJed Brown } 425c4762a1bSJed Brown } 426c4762a1bSJed Brown for (d = 0; d < cdim; ++d) faceCentroidExT[d] = rot[d]; 427c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 428c4762a1bSJed Brown for (e = 0, rot[d] = 0.0; e < cdim; ++e) { 429c4762a1bSJed Brown rot[d] += R[d*cdim+e] * faceNormalExT[e]; 430c4762a1bSJed Brown } 431c4762a1bSJed Brown } 432c4762a1bSJed Brown for (d = 0; d < cdim; ++d) faceNormalExT[d] = rot[d]; 433c4762a1bSJed Brown for (d = 0, norm = 0.0; d < cdim; ++d) norm += PetscSqr(faceNormalExT[d]); 434c4762a1bSJed Brown norm = PetscSqrtReal(norm); 435c4762a1bSJed Brown for (d = 0; d < cdim; ++d) faceNormalExT[d] /= norm; 436c4762a1bSJed Brown 4379566063dSJacob Faibussowitsch PetscCall(CheckFVMGeometry(dm, cone[c], cdim, faceCentroidExT, faceNormalExT, faceVolExT)); 438c4762a1bSJed Brown } 439c4762a1bSJed Brown } 440c4762a1bSJed Brown } 441c4762a1bSJed Brown } 4429566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&r)); 4439566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&ang)); 4449566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&ang2)); 4459566063dSJacob Faibussowitsch PetscCall(PetscFree2(origCoords, newCoords)); 4469566063dSJacob Faibussowitsch PetscCall(PetscFree5(v0ExT, JExT, invJExT, centroidExT, normalExT)); 4479566063dSJacob Faibussowitsch PetscCall(PetscFree2(faceCentroidExT, faceNormalExT)); 448c4762a1bSJed Brown } 449c4762a1bSJed Brown PetscFunctionReturn(0); 450c4762a1bSJed Brown } 451c4762a1bSJed Brown 4524f99dae5SMatthew G. Knepley static PetscErrorCode TestTriangle(MPI_Comm comm, PetscBool transform) 453c4762a1bSJed Brown { 454c4762a1bSJed Brown DM dm; 455c4762a1bSJed Brown 456c4762a1bSJed Brown PetscFunctionBegin; 4579566063dSJacob Faibussowitsch PetscCall(DMPlexCreateReferenceCell(comm, DM_POLYTOPE_TRIANGLE, &dm)); 4589566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 459c4762a1bSJed Brown /* Check reference geometry: determinant is scaled by reference volume (2.0) */ 460c4762a1bSJed Brown { 461c4762a1bSJed Brown PetscReal v0Ex[2] = {-1.0, -1.0}; 462c4762a1bSJed Brown PetscReal JEx[4] = {1.0, 0.0, 0.0, 1.0}; 463c4762a1bSJed Brown PetscReal invJEx[4] = {1.0, 0.0, 0.0, 1.0}; 464c4762a1bSJed Brown PetscReal detJEx = 1.0; 465c4762a1bSJed Brown PetscReal centroidEx[2] = {-((PetscReal)1.)/((PetscReal)3.), -((PetscReal)1.)/((PetscReal)3.)}; 466c4762a1bSJed Brown PetscReal normalEx[2] = {0.0, 0.0}; 467c4762a1bSJed Brown PetscReal volEx = 2.0; 468c4762a1bSJed Brown 4699566063dSJacob Faibussowitsch PetscCall(CheckCell(dm, 0, transform, v0Ex, JEx, invJEx, detJEx, centroidEx, normalEx, volEx, NULL, NULL, NULL)); 470c4762a1bSJed Brown } 471c4762a1bSJed Brown /* Move to 3D: Check reference geometry: determinant is scaled by reference volume (2.0) */ 472c4762a1bSJed Brown { 473c4762a1bSJed Brown PetscScalar vertexCoords[9] = {-1.0, -1.0, 0.0, 1.0, -1.0, 0.0, -1.0, 1.0, 0.0}; 474c4762a1bSJed Brown PetscReal v0Ex[3] = {-1.0, -1.0, 0.0}; 475c4762a1bSJed Brown PetscReal JEx[9] = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0}; 476c4762a1bSJed Brown PetscReal invJEx[9] = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0}; 477c4762a1bSJed Brown PetscReal detJEx = 1.0; 478c4762a1bSJed Brown PetscReal centroidEx[3] = {-((PetscReal)1.)/((PetscReal)3.), -((PetscReal)1.)/((PetscReal)3.), 0.0}; 479c4762a1bSJed Brown PetscReal normalEx[3] = {0.0, 0.0, 1.0}; 480c4762a1bSJed Brown PetscReal volEx = 2.0; 481c4762a1bSJed Brown 4829566063dSJacob Faibussowitsch PetscCall(ChangeCoordinates(dm, 3, vertexCoords)); 4839566063dSJacob Faibussowitsch PetscCall(CheckCell(dm, 0, transform, v0Ex, JEx, invJEx, detJEx, centroidEx, normalEx, volEx, NULL, NULL, NULL)); 484c4762a1bSJed Brown } 485c4762a1bSJed Brown /* Cleanup */ 4869566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 487c4762a1bSJed Brown PetscFunctionReturn(0); 488c4762a1bSJed Brown } 489c4762a1bSJed Brown 4904f99dae5SMatthew G. Knepley static PetscErrorCode TestQuadrilateral(MPI_Comm comm, PetscBool transform) 491c4762a1bSJed Brown { 492c4762a1bSJed Brown DM dm; 493c4762a1bSJed Brown 494c4762a1bSJed Brown PetscFunctionBegin; 4959566063dSJacob Faibussowitsch PetscCall(DMPlexCreateReferenceCell(comm, DM_POLYTOPE_QUADRILATERAL, &dm)); 4969566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 497c4762a1bSJed Brown /* Check reference geometry: determinant is scaled by reference volume (2.0) */ 498c4762a1bSJed Brown { 499c4762a1bSJed Brown PetscReal v0Ex[2] = {-1.0, -1.0}; 500c4762a1bSJed Brown PetscReal JEx[4] = {1.0, 0.0, 0.0, 1.0}; 501c4762a1bSJed Brown PetscReal invJEx[4] = {1.0, 0.0, 0.0, 1.0}; 502c4762a1bSJed Brown PetscReal detJEx = 1.0; 503c4762a1bSJed Brown PetscReal centroidEx[2] = {0.0, 0.0}; 504c4762a1bSJed Brown PetscReal normalEx[2] = {0.0, 0.0}; 505c4762a1bSJed Brown PetscReal volEx = 4.0; 506c4762a1bSJed Brown 5079566063dSJacob Faibussowitsch PetscCall(CheckCell(dm, 0, transform, v0Ex, JEx, invJEx, detJEx, centroidEx, normalEx, volEx, NULL, NULL, NULL)); 508c4762a1bSJed Brown } 509c4762a1bSJed Brown /* Move to 3D: Check reference geometry: determinant is scaled by reference volume (4.0) */ 510c4762a1bSJed Brown { 511c4762a1bSJed Brown PetscScalar vertexCoords[12] = {-1.0, -1.0, 0.0, 1.0, -1.0, 0.0, 1.0, 1.0, 0.0, -1.0, 1.0, 0.0}; 512c4762a1bSJed Brown PetscReal v0Ex[3] = {-1.0, -1.0, 0.0}; 513c4762a1bSJed Brown PetscReal JEx[9] = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0}; 514c4762a1bSJed Brown PetscReal invJEx[9] = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0}; 515c4762a1bSJed Brown PetscReal detJEx = 1.0; 516c4762a1bSJed Brown PetscReal centroidEx[3] = {0.0, 0.0, 0.0}; 517c4762a1bSJed Brown PetscReal normalEx[3] = {0.0, 0.0, 1.0}; 518c4762a1bSJed Brown PetscReal volEx = 4.0; 519c4762a1bSJed Brown 5209566063dSJacob Faibussowitsch PetscCall(ChangeCoordinates(dm, 3, vertexCoords)); 5219566063dSJacob Faibussowitsch PetscCall(CheckCell(dm, 0, transform, v0Ex, JEx, invJEx, detJEx, centroidEx, normalEx, volEx, NULL, NULL, NULL)); 522c4762a1bSJed Brown } 523c4762a1bSJed Brown /* Cleanup */ 5249566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 525c4762a1bSJed Brown PetscFunctionReturn(0); 526c4762a1bSJed Brown } 527c4762a1bSJed Brown 5284f99dae5SMatthew G. Knepley static PetscErrorCode TestTetrahedron(MPI_Comm comm, PetscBool transform) 529c4762a1bSJed Brown { 530c4762a1bSJed Brown DM dm; 531c4762a1bSJed Brown 532c4762a1bSJed Brown PetscFunctionBegin; 5339566063dSJacob Faibussowitsch PetscCall(DMPlexCreateReferenceCell(comm, DM_POLYTOPE_TETRAHEDRON, &dm)); 5349566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 535c4762a1bSJed Brown /* Check reference geometry: determinant is scaled by reference volume (4/3) */ 536c4762a1bSJed Brown { 537c4762a1bSJed Brown PetscReal v0Ex[3] = {-1.0, -1.0, -1.0}; 538c4762a1bSJed Brown PetscReal JEx[9] = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0}; 539c4762a1bSJed Brown PetscReal invJEx[9] = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0}; 540c4762a1bSJed Brown PetscReal detJEx = 1.0; 541c4762a1bSJed Brown PetscReal centroidEx[3] = {-0.5, -0.5, -0.5}; 542c4762a1bSJed Brown PetscReal normalEx[3] = {0.0, 0.0, 0.0}; 543c4762a1bSJed Brown PetscReal volEx = (PetscReal)4.0/(PetscReal)3.0; 544c4762a1bSJed Brown 5459566063dSJacob Faibussowitsch PetscCall(CheckCell(dm, 0, transform, v0Ex, JEx, invJEx, detJEx, centroidEx, normalEx, volEx, NULL, NULL, NULL)); 546c4762a1bSJed Brown } 547c4762a1bSJed Brown /* Cleanup */ 5489566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 549c4762a1bSJed Brown PetscFunctionReturn(0); 550c4762a1bSJed Brown } 551c4762a1bSJed Brown 5524f99dae5SMatthew G. Knepley static PetscErrorCode TestHexahedron(MPI_Comm comm, PetscBool transform) 553c4762a1bSJed Brown { 554c4762a1bSJed Brown DM dm; 555c4762a1bSJed Brown 556c4762a1bSJed Brown PetscFunctionBegin; 5579566063dSJacob Faibussowitsch PetscCall(DMPlexCreateReferenceCell(comm, DM_POLYTOPE_HEXAHEDRON, &dm)); 5589566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 559c4762a1bSJed Brown /* Check reference geometry: determinant is scaled by reference volume 8.0 */ 560c4762a1bSJed Brown { 561c4762a1bSJed Brown PetscReal v0Ex[3] = {-1.0, -1.0, -1.0}; 562c4762a1bSJed Brown PetscReal JEx[9] = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0}; 563c4762a1bSJed Brown PetscReal invJEx[9] = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0}; 564c4762a1bSJed Brown PetscReal detJEx = 1.0; 565c4762a1bSJed Brown PetscReal centroidEx[3] = {0.0, 0.0, 0.0}; 566c4762a1bSJed Brown PetscReal normalEx[3] = {0.0, 0.0, 0.0}; 567c4762a1bSJed Brown PetscReal volEx = 8.0; 568c4762a1bSJed Brown 5699566063dSJacob Faibussowitsch PetscCall(CheckCell(dm, 0, transform, v0Ex, JEx, invJEx, detJEx, centroidEx, normalEx, volEx, NULL, NULL, NULL)); 570c4762a1bSJed Brown } 571c4762a1bSJed Brown /* Cleanup */ 5729566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 573c4762a1bSJed Brown PetscFunctionReturn(0); 574c4762a1bSJed Brown } 575c4762a1bSJed Brown 5764f99dae5SMatthew G. Knepley static PetscErrorCode TestHexahedronCurved(MPI_Comm comm) 577c4762a1bSJed Brown { 578c4762a1bSJed Brown DM dm; 5794f99dae5SMatthew G. Knepley PetscScalar coords[24] = {-1.0, -1.0, -1.0, -1.0, 1.0, -1.0, 1.0, 1.0, -1.0, 1.0, -1.0, -1.0, 5804f99dae5SMatthew G. Knepley -1.0, -1.0, 1.1, 1.0, -1.0, 1.0, 1.0, 1.0, 1.1, -1.0, 1.0, 1.0}; 581c4762a1bSJed Brown 582c4762a1bSJed Brown PetscFunctionBegin; 5839566063dSJacob Faibussowitsch PetscCall(DMPlexCreateReferenceCell(comm, DM_POLYTOPE_HEXAHEDRON, &dm)); 5849566063dSJacob Faibussowitsch PetscCall(ChangeCoordinates(dm, 3, coords)); 5859566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 5864f99dae5SMatthew G. Knepley { 5874f99dae5SMatthew G. Knepley PetscReal centroidEx[3] = {0.0, 0.0, 0.016803278688524603}; 5884f99dae5SMatthew G. Knepley PetscReal normalEx[3] = {0.0, 0.0, 0.0}; 5894f99dae5SMatthew G. Knepley PetscReal volEx = 8.1333333333333346; 5904f99dae5SMatthew G. Knepley 5919566063dSJacob Faibussowitsch PetscCall(CheckCell(dm, 0, PETSC_FALSE, NULL, NULL, NULL, 0.0, centroidEx, normalEx, volEx, NULL, NULL, NULL)); 592c4762a1bSJed Brown } 5939566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 5944f99dae5SMatthew G. Knepley PetscFunctionReturn(0); 5954f99dae5SMatthew G. Knepley } 5964f99dae5SMatthew G. Knepley 5974f99dae5SMatthew G. Knepley /* This wedge is a tensor product cell, rather than a normal wedge */ 5984f99dae5SMatthew G. Knepley static PetscErrorCode TestWedge(MPI_Comm comm, PetscBool transform) 5994f99dae5SMatthew G. Knepley { 6004f99dae5SMatthew G. Knepley DM dm; 6014f99dae5SMatthew G. Knepley 6024f99dae5SMatthew G. Knepley PetscFunctionBegin; 6039566063dSJacob Faibussowitsch PetscCall(DMPlexCreateReferenceCell(comm, DM_POLYTOPE_TRI_PRISM_TENSOR, &dm)); 6049566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 605c4762a1bSJed Brown /* Check reference geometry: determinant is scaled by reference volume 4.0 */ 606c4762a1bSJed Brown { 607c4762a1bSJed Brown #if 0 608c4762a1bSJed Brown /* FEM geometry is not functional for wedges */ 609c4762a1bSJed Brown PetscReal v0Ex[3] = {-1.0, -1.0, -1.0}; 610c4762a1bSJed Brown PetscReal JEx[9] = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0}; 611c4762a1bSJed Brown PetscReal invJEx[9] = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0}; 612c4762a1bSJed Brown PetscReal detJEx = 1.0; 613c4762a1bSJed Brown #endif 614c4762a1bSJed Brown 6154f99dae5SMatthew G. Knepley { 616c4762a1bSJed Brown PetscReal centroidEx[3] = {-((PetscReal)1.)/((PetscReal)3.), -((PetscReal)1.)/((PetscReal)3.), 0.0}; 617c4762a1bSJed Brown PetscReal normalEx[3] = {0.0, 0.0, 0.0}; 618c4762a1bSJed Brown PetscReal volEx = 4.0; 619c4762a1bSJed Brown PetscReal faceVolEx[5] = {2.0, 2.0, 4.0, PETSC_SQRT2*4.0, 4.0}; 620c4762a1bSJed Brown PetscReal faceNormalEx[15] = {0.0, 0.0, 1.0, 0.0, 0.0, 1.0, 0.0, -1.0, 0.0, PETSC_SQRT2/2.0, PETSC_SQRT2/2.0, 0.0, -1.0, 0.0, 0.0}; 621c4762a1bSJed Brown PetscReal faceCentroidEx[15] = {-((PetscReal)1.)/((PetscReal)3.), -((PetscReal)1.)/((PetscReal)3.), -1.0, 622c4762a1bSJed Brown -((PetscReal)1.)/((PetscReal)3.), -((PetscReal)1.)/((PetscReal)3.), 1.0, 623c4762a1bSJed Brown 0.0, -1.0, 0.0, 0.0, 0.0, 0.0, -1.0, 0.0, 0.0}; 624c4762a1bSJed Brown 6259566063dSJacob Faibussowitsch PetscCall(CheckCell(dm, 0, transform, NULL, NULL, NULL, 0.0, centroidEx, normalEx, volEx, faceCentroidEx, faceNormalEx, faceVolEx)); 626c4762a1bSJed Brown } 627c4762a1bSJed Brown } 628c4762a1bSJed Brown /* Cleanup */ 6299566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 630c4762a1bSJed Brown PetscFunctionReturn(0); 631c4762a1bSJed Brown } 632c4762a1bSJed Brown 633c4762a1bSJed Brown int main(int argc, char **argv) 634c4762a1bSJed Brown { 635c4762a1bSJed Brown AppCtx user; 636c4762a1bSJed Brown 6379566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL,help)); 6389566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 639c4762a1bSJed Brown if (user.runType == RUN_REFERENCE) { 6409566063dSJacob Faibussowitsch PetscCall(TestTriangle(PETSC_COMM_SELF, user.transform)); 6419566063dSJacob Faibussowitsch PetscCall(TestQuadrilateral(PETSC_COMM_SELF, user.transform)); 6429566063dSJacob Faibussowitsch PetscCall(TestTetrahedron(PETSC_COMM_SELF, user.transform)); 6439566063dSJacob Faibussowitsch PetscCall(TestHexahedron(PETSC_COMM_SELF, user.transform)); 6449566063dSJacob Faibussowitsch PetscCall(TestWedge(PETSC_COMM_SELF, user.transform)); 6454f99dae5SMatthew G. Knepley } else if (user.runType == RUN_HEX_CURVED) { 6469566063dSJacob Faibussowitsch PetscCall(TestHexahedronCurved(PETSC_COMM_SELF)); 647c4762a1bSJed Brown } else if (user.runType == RUN_FILE) { 648c4762a1bSJed Brown PetscInt dim, cStart, cEnd, c; 649c4762a1bSJed Brown 6509566063dSJacob Faibussowitsch PetscCall(DMGetDimension(user.dm, &dim)); 6519566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(user.dm, 0, &cStart, &cEnd)); 652c4762a1bSJed Brown for (c = 0; c < cEnd-cStart; ++c) { 6534f99dae5SMatthew G. Knepley PetscReal *v0 = user.v0 ? &user.v0[c*dim] : NULL; 6544f99dae5SMatthew G. Knepley PetscReal *J = user.J ? &user.J[c*dim*dim] : NULL; 6554f99dae5SMatthew G. Knepley PetscReal *invJ = user.invJ ? &user.invJ[c*dim*dim] : NULL; 6564f99dae5SMatthew G. Knepley PetscReal detJ = user.detJ ? user.detJ[c] : 0.0; 6574f99dae5SMatthew G. Knepley PetscReal *centroid = user.centroid ? &user.centroid[c*dim] : NULL; 6584f99dae5SMatthew G. Knepley PetscReal *normal = user.normal ? &user.normal[c*dim] : NULL; 6594f99dae5SMatthew G. Knepley PetscReal vol = user.vol ? user.vol[c] : 0.0; 6604f99dae5SMatthew G. Knepley 6619566063dSJacob Faibussowitsch PetscCall(CheckCell(user.dm, c+cStart, PETSC_FALSE, v0, J, invJ, detJ, centroid, normal, vol, NULL, NULL, NULL)); 662c4762a1bSJed Brown } 6639566063dSJacob Faibussowitsch PetscCall(PetscFree4(user.v0,user.J,user.invJ,user.detJ)); 6649566063dSJacob Faibussowitsch PetscCall(PetscFree(user.centroid)); 6659566063dSJacob Faibussowitsch PetscCall(PetscFree(user.normal)); 6669566063dSJacob Faibussowitsch PetscCall(PetscFree(user.vol)); 6679566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.dm)); 668c4762a1bSJed Brown } else if (user.runType == RUN_DISPLAY) { 669c4762a1bSJed Brown DM gdm, dmCell; 670c4762a1bSJed Brown Vec cellgeom, facegeom; 671c4762a1bSJed Brown const PetscScalar *cgeom; 672c4762a1bSJed Brown PetscInt dim, d, cStart, cEnd, cEndInterior, c; 673c4762a1bSJed Brown 6749566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(user.dm, &dim)); 6759566063dSJacob Faibussowitsch PetscCall(DMPlexConstructGhostCells(user.dm, NULL, NULL, &gdm)); 676c4762a1bSJed Brown if (gdm) { 6779566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.dm)); 678c4762a1bSJed Brown user.dm = gdm; 679c4762a1bSJed Brown } 6809566063dSJacob Faibussowitsch PetscCall(DMPlexComputeGeometryFVM(user.dm, &cellgeom, &facegeom)); 6819566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(user.dm, 0, &cStart, &cEnd)); 6829566063dSJacob Faibussowitsch PetscCall(DMPlexGetGhostCellStratum(user.dm, &cEndInterior, NULL)); 683c4762a1bSJed Brown if (cEndInterior >= 0) cEnd = cEndInterior; 6849566063dSJacob Faibussowitsch PetscCall(VecGetDM(cellgeom, &dmCell)); 6859566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(cellgeom, &cgeom)); 686c4762a1bSJed Brown for (c = 0; c < cEnd-cStart; ++c) { 687c4762a1bSJed Brown PetscFVCellGeom *cg; 688c4762a1bSJed Brown 6899566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRead(dmCell, c, cgeom, &cg)); 6909566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Cell %4D: Centroid (", c)); 691c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 6929566063dSJacob Faibussowitsch if (d > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", ")); 6939566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "%12.2g", cg->centroid[d])); 694c4762a1bSJed Brown } 6959566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, ") Vol %12.2g\n", cg->volume)); 696c4762a1bSJed Brown } 6979566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(cellgeom, &cgeom)); 6989566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cellgeom)); 6999566063dSJacob Faibussowitsch PetscCall(VecDestroy(&facegeom)); 7009566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.dm)); 701c4762a1bSJed Brown } 7029566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 703b122ec5aSJacob Faibussowitsch return 0; 704c4762a1bSJed Brown } 705c4762a1bSJed Brown 706c4762a1bSJed Brown /*TEST 707c4762a1bSJed Brown 708c4762a1bSJed Brown test: 7094f99dae5SMatthew G. Knepley suffix: 1 710c4762a1bSJed Brown args: -dm_view ascii::ascii_info_detail 711c4762a1bSJed Brown test: 712c4762a1bSJed Brown suffix: 2 7134f99dae5SMatthew G. Knepley args: -run_type hex_curved 714c4762a1bSJed Brown test: 715c4762a1bSJed Brown suffix: 3 7164f99dae5SMatthew G. Knepley args: -transform 717c4762a1bSJed Brown test: 718c4762a1bSJed Brown suffix: 4 719c4762a1bSJed Brown requires: exodusii 7204f99dae5SMatthew G. Knepley args: -run_type file -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/simpleblock-100.exo -dm_view ascii::ascii_info_detail -v0 -1.5,-0.5,0.5,-0.5,-0.5,0.5,0.5,-0.5,0.5 -J 0.0,0.0,0.5,0.0,0.5,0.0,-0.5,0.0,0.0,0.0,0.0,0.5,0.0,0.5,0.0,-0.5,0.0,0.0,0.0,0.0,0.5,0.0,0.5,0.0,-0.5,0.0,0.0 -invJ 0.0,0.0,-2.0,0.0,2.0,0.0,2.0,0.0,0.0,0.0,0.0,-2.0,0.0,2.0,0.0,2.0,0.0,0.0,0.0,0.0,-2.0,0.0,2.0,0.0,2.0,0.0,0.0 -detJ 0.125,0.125,0.125 -centroid -1.0,0.0,0.0,0.0,0.0,0.0,1.0,0.0,0.0 -normal 0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0 -vol 1.0,1.0,1.0 721c4762a1bSJed Brown test: 722c4762a1bSJed Brown suffix: 5 7234f99dae5SMatthew G. Knepley args: -run_type file -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_faces 3,1,1 -dm_plex_box_lower -1.5,-0.5,-0.5 -dm_plex_box_upper 1.5,0.5,0.5 -dm_view ascii::ascii_info_detail -centroid -1.0,0.0,0.0,0.0,0.0,0.0,1.0,0.0,0.0 -normal 0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0 -vol 1.0,1.0,1.0 7249bf2564aSMatt McGurn test: 7259bf2564aSMatt McGurn suffix: 6 7269bf2564aSMatt McGurn args: -run_type file -dm_plex_dim 1 -dm_plex_simplex 0 -dm_plex_box_faces 3 -dm_plex_box_lower -1.5 -dm_plex_box_upper 1.5 -dm_view ascii::ascii_info_detail -centroid -1.0,0.0,1.0 -vol 1.0,1.0,1.0 727c4762a1bSJed Brown TEST*/ 728