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 33c4762a1bSJed Brown PetscFunctionBeginUser; 34c4762a1bSJed Brown options->runType = RUN_REFERENCE; 35c4762a1bSJed Brown options->transform = PETSC_FALSE; 36c4762a1bSJed Brown 37d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Geometry Test Options", "DMPLEX"); 38c4762a1bSJed Brown run = options->runType; 399566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-run_type", "The run type", "ex8.c", runTypes, 3, runTypes[options->runType], &run, NULL)); 40c4762a1bSJed Brown options->runType = (RunType) run; 419566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-transform", "Use random transforms", "ex8.c", options->transform, &options->transform, NULL)); 42c4762a1bSJed Brown 43c4762a1bSJed Brown if (options->runType == RUN_FILE) { 44c4762a1bSJed Brown PetscInt dim, cStart, cEnd, numCells, n; 459bf2564aSMatt McGurn PetscBool flg, feFlg; 46c4762a1bSJed Brown 479566063dSJacob Faibussowitsch PetscCall(ReadMesh(PETSC_COMM_WORLD, options, &options->dm)); 489566063dSJacob Faibussowitsch PetscCall(DMGetDimension(options->dm, &dim)); 499566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(options->dm, 0, &cStart, &cEnd)); 50c4762a1bSJed Brown numCells = cEnd-cStart; 519566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(numCells*dim, &options->v0, numCells*dim*dim, &options->J, numCells*dim*dim, &options->invJ, numCells, &options->detJ)); 529566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numCells*dim, &options->centroid)); 539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numCells*dim, &options->normal)); 549566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numCells, &options->vol)); 55c4762a1bSJed Brown n = numCells*dim; 569566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-v0", "Input v0 for each cell", "ex8.c", options->v0, &n, &feFlg)); 5763a3b9bcSJacob Faibussowitsch PetscCheck(!feFlg || n == numCells*dim, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid size of v0 %" PetscInt_FMT " should be %" PetscInt_FMT, n, numCells*dim); 58c4762a1bSJed Brown n = numCells*dim*dim; 599566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-J", "Input Jacobian for each cell", "ex8.c", options->J, &n, &flg)); 6063a3b9bcSJacob Faibussowitsch PetscCheck(!flg || n == numCells*dim*dim, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid size of J %" PetscInt_FMT " should be %" PetscInt_FMT, n, numCells*dim*dim); 61c4762a1bSJed Brown n = numCells*dim*dim; 629566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-invJ", "Input inverse Jacobian for each cell", "ex8.c", options->invJ, &n, &flg)); 6363a3b9bcSJacob Faibussowitsch PetscCheck(!flg || n == numCells*dim*dim, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid size of invJ %" PetscInt_FMT " should be %" PetscInt_FMT, n, numCells*dim*dim); 64c4762a1bSJed Brown n = numCells; 659566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-detJ", "Input Jacobian determinant for each cell", "ex8.c", options->detJ, &n, &flg)); 6663a3b9bcSJacob Faibussowitsch PetscCheck(!flg || n == numCells, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid size of detJ %" PetscInt_FMT " should be %" PetscInt_FMT, n, numCells); 67c4762a1bSJed Brown n = numCells*dim; 684f99dae5SMatthew G. Knepley if (!feFlg) { 699566063dSJacob Faibussowitsch PetscCall(PetscFree4(options->v0, options->J, options->invJ, options->detJ)); 704f99dae5SMatthew G. Knepley options->v0 = options->J = options->invJ = options->detJ = NULL; 714f99dae5SMatthew G. Knepley } 729566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-centroid", "Input centroid for each cell", "ex8.c", options->centroid, &n, &flg)); 7363a3b9bcSJacob Faibussowitsch PetscCheck(!flg || n == numCells*dim, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid size of centroid %" PetscInt_FMT " should be %" PetscInt_FMT, n, numCells*dim); 749bf2564aSMatt McGurn if (!flg) { 759566063dSJacob Faibussowitsch PetscCall(PetscFree(options->centroid)); 769bf2564aSMatt McGurn options->centroid = NULL; 779bf2564aSMatt McGurn } 78c4762a1bSJed Brown n = numCells*dim; 799566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-normal", "Input normal for each cell", "ex8.c", options->normal, &n, &flg)); 8063a3b9bcSJacob Faibussowitsch PetscCheck(!flg || n == numCells*dim, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid size of normal %" PetscInt_FMT " should be %" PetscInt_FMT, n, numCells*dim); 819bf2564aSMatt McGurn if (!flg) { 829566063dSJacob Faibussowitsch PetscCall(PetscFree(options->normal)); 839bf2564aSMatt McGurn options->normal = NULL; 849bf2564aSMatt McGurn } 85c4762a1bSJed Brown n = numCells; 869566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-vol", "Input volume for each cell", "ex8.c", options->vol, &n, &flg)); 8763a3b9bcSJacob Faibussowitsch PetscCheck(!flg || n == numCells, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid size of vol %" PetscInt_FMT " should be %" PetscInt_FMT, n, numCells); 889bf2564aSMatt McGurn if (!flg) { 899566063dSJacob Faibussowitsch PetscCall(PetscFree(options->vol)); 909bf2564aSMatt McGurn options->vol = NULL; 914f99dae5SMatthew G. Knepley } 92c4762a1bSJed Brown } else if (options->runType == RUN_DISPLAY) { 939566063dSJacob Faibussowitsch PetscCall(ReadMesh(PETSC_COMM_WORLD, options, &options->dm)); 94c4762a1bSJed Brown } 95d0609cedSBarry Smith PetscOptionsEnd(); 96c4762a1bSJed Brown 979566063dSJacob Faibussowitsch if (options->transform) PetscCall(PetscPrintf(comm, "Using random transforms\n")); 98c4762a1bSJed Brown PetscFunctionReturn(0); 99c4762a1bSJed Brown } 100c4762a1bSJed Brown 101c4762a1bSJed Brown static PetscErrorCode ChangeCoordinates(DM dm, PetscInt spaceDim, PetscScalar vertexCoords[]) 102c4762a1bSJed Brown { 103c4762a1bSJed Brown PetscSection coordSection; 104c4762a1bSJed Brown Vec coordinates; 105c4762a1bSJed Brown PetscScalar *coords; 106c4762a1bSJed Brown PetscInt vStart, vEnd, v, d, coordSize; 107c4762a1bSJed Brown 108c4762a1bSJed Brown PetscFunctionBegin; 1099566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 1109566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm, &coordSection)); 1119566063dSJacob Faibussowitsch PetscCall(PetscSectionSetNumFields(coordSection, 1)); 1129566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldComponents(coordSection, 0, spaceDim)); 1139566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(coordSection, vStart, vEnd)); 114c4762a1bSJed Brown for (v = vStart; v < vEnd; ++v) { 1159566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(coordSection, v, spaceDim)); 1169566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, spaceDim)); 117c4762a1bSJed Brown } 1189566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(coordSection)); 1199566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize)); 1209566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates)); 1219566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) coordinates, "coordinates")); 1229566063dSJacob Faibussowitsch PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE)); 1239566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(coordinates)); 1249566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordinates, &coords)); 125c4762a1bSJed Brown for (v = vStart; v < vEnd; ++v) { 126c4762a1bSJed Brown PetscInt off; 127c4762a1bSJed Brown 1289566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(coordSection, v, &off)); 129c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) { 130c4762a1bSJed Brown coords[off+d] = vertexCoords[(v-vStart)*spaceDim+d]; 131c4762a1bSJed Brown } 132c4762a1bSJed Brown } 1339566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coordinates, &coords)); 1349566063dSJacob Faibussowitsch PetscCall(DMSetCoordinateDim(dm, spaceDim)); 1359566063dSJacob Faibussowitsch PetscCall(DMSetCoordinatesLocal(dm, coordinates)); 1369566063dSJacob Faibussowitsch PetscCall(VecDestroy(&coordinates)); 1379566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 138c4762a1bSJed Brown PetscFunctionReturn(0); 139c4762a1bSJed Brown } 140c4762a1bSJed Brown 141c4762a1bSJed Brown #define RelativeError(a,b) PetscAbs(a-b)/(1.0+PetscMax(PetscAbs(a),PetscAbs(b))) 142c4762a1bSJed Brown 143c4762a1bSJed Brown static PetscErrorCode CheckFEMGeometry(DM dm, PetscInt cell, PetscInt spaceDim, PetscReal v0Ex[], PetscReal JEx[], PetscReal invJEx[], PetscReal detJEx) 144c4762a1bSJed Brown { 145c4762a1bSJed Brown PetscReal v0[3], J[9], invJ[9], detJ; 146c4762a1bSJed Brown PetscInt d, i, j; 147c4762a1bSJed Brown 148c4762a1bSJed Brown PetscFunctionBegin; 1499566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFEM(dm, cell, NULL, v0, J, invJ, &detJ)); 150c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) { 151c4762a1bSJed Brown if (v0[d] != v0Ex[d]) { 152c4762a1bSJed Brown switch (spaceDim) { 15398921bdaSJacob 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]); 15498921bdaSJacob 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]); 15563a3b9bcSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid space dimension %" PetscInt_FMT, spaceDim); 156c4762a1bSJed Brown } 157c4762a1bSJed Brown } 158c4762a1bSJed Brown } 159c4762a1bSJed Brown for (i = 0; i < spaceDim; ++i) { 160c4762a1bSJed Brown for (j = 0; j < spaceDim; ++j) { 16163a3b9bcSJacob Faibussowitsch PetscCheck(RelativeError(J[i*spaceDim+j], JEx[i*spaceDim+j]) < 10*PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid J[%" PetscInt_FMT ",%" PetscInt_FMT "]: %g != %g", i, j, (double)J[i*spaceDim+j], (double)JEx[i*spaceDim+j]); 16263a3b9bcSJacob Faibussowitsch PetscCheck(RelativeError(invJ[i*spaceDim+j], invJEx[i*spaceDim+j]) < 10*PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid invJ[%" PetscInt_FMT ",%" PetscInt_FMT "]: %g != %g", i, j, (double)invJ[i*spaceDim+j], (double)invJEx[i*spaceDim+j]); 163c4762a1bSJed Brown } 164c4762a1bSJed Brown } 165700a43edSPierre 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)); 166c4762a1bSJed Brown PetscFunctionReturn(0); 167c4762a1bSJed Brown } 168c4762a1bSJed Brown 169c4762a1bSJed Brown static PetscErrorCode CheckFVMGeometry(DM dm, PetscInt cell, PetscInt spaceDim, PetscReal centroidEx[], PetscReal normalEx[], PetscReal volEx) 170c4762a1bSJed Brown { 1714f99dae5SMatthew G. Knepley PetscReal tol = PetscMax(10*PETSC_SMALL, 1e-10); 172c4762a1bSJed Brown PetscReal centroid[3], normal[3], vol; 173c4762a1bSJed Brown PetscInt d; 174c4762a1bSJed Brown 175c4762a1bSJed Brown PetscFunctionBegin; 1769566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFVM(dm, cell, volEx? &vol : NULL, centroidEx? centroid : NULL, normalEx? normal : NULL)); 177c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) { 1789bf2564aSMatt McGurn if (centroidEx) 17963a3b9bcSJacob Faibussowitsch PetscCheck(RelativeError(centroid[d], centroidEx[d]) < tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cell %" PetscInt_FMT ", Invalid centroid[%" PetscInt_FMT "]: %g != %g diff %g", cell, d, (double)centroid[d], (double)centroidEx[d],(double)(centroid[d]-centroidEx[d])); 1809bf2564aSMatt McGurn if (normalEx) 18163a3b9bcSJacob Faibussowitsch PetscCheck(RelativeError(normal[d], normalEx[d]) < tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cell %" PetscInt_FMT ", Invalid normal[%" PetscInt_FMT "]: %g != %g", cell, d, (double) normal[d], (double) normalEx[d]); 182c4762a1bSJed Brown } 1839bf2564aSMatt McGurn if (volEx) 18463a3b9bcSJacob Faibussowitsch PetscCheck(RelativeError(volEx, vol) < tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cell %" PetscInt_FMT ", Invalid volume = %g != %g diff %g", cell, (double)vol, (double)volEx,(double)(vol - volEx)); 1854f99dae5SMatthew G. Knepley PetscFunctionReturn(0); 1864f99dae5SMatthew G. Knepley } 1874f99dae5SMatthew G. Knepley 1884f99dae5SMatthew G. Knepley static PetscErrorCode CheckGaussLaw(DM dm, PetscInt cell) 1894f99dae5SMatthew G. Knepley { 1904f99dae5SMatthew G. Knepley DMPolytopeType ct; 1914f99dae5SMatthew G. Knepley PetscReal tol = PetscMax(10*PETSC_SMALL, 1e-10); 1924f99dae5SMatthew G. Knepley PetscReal normal[3], integral[3] = {0., 0., 0.}, area; 1934f99dae5SMatthew G. Knepley const PetscInt *cone, *ornt; 1944f99dae5SMatthew G. Knepley PetscInt coneSize, f, dim, cdim, d; 1954f99dae5SMatthew G. Knepley 1964f99dae5SMatthew G. Knepley PetscFunctionBegin; 1979566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1989566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &cdim)); 1994f99dae5SMatthew G. Knepley if (dim != cdim) PetscFunctionReturn(0); 2009566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, cell, &ct)); 2014f99dae5SMatthew G. Knepley if (ct == DM_POLYTOPE_TRI_PRISM_TENSOR) PetscFunctionReturn(0); 2029566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, cell, &coneSize)); 2039566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, cell, &cone)); 2049566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeOrientation(dm, cell, &ornt)); 2054f99dae5SMatthew G. Knepley for (f = 0; f < coneSize; ++f) { 2069bf2564aSMatt McGurn const PetscInt sgn = dim == 1? (f == 0 ? -1 : 1) : (ornt[f] < 0 ? -1 : 1); 2074f99dae5SMatthew G. Knepley 2089566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFVM(dm, cone[f], &area, NULL, normal)); 2094f99dae5SMatthew G. Knepley for (d = 0; d < cdim; ++d) integral[d] += sgn*area*normal[d]; 2104f99dae5SMatthew G. Knepley } 21163a3b9bcSJacob Faibussowitsch for (d = 0; d < cdim; ++d) PetscCheck(PetscAbsReal(integral[d]) < tol, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cell %" PetscInt_FMT " Surface integral for component %" PetscInt_FMT ": %g != 0. as it should be for a constant field", cell, d, (double) integral[d]); 212c4762a1bSJed Brown PetscFunctionReturn(0); 213c4762a1bSJed Brown } 214c4762a1bSJed Brown 215c4762a1bSJed 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[]) 216c4762a1bSJed Brown { 217c4762a1bSJed Brown const PetscInt *cone; 218c4762a1bSJed Brown PetscInt coneSize, c; 219c4762a1bSJed Brown PetscInt dim, depth, cdim; 220c4762a1bSJed Brown 221c4762a1bSJed Brown PetscFunctionBegin; 2229566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepth(dm, &depth)); 2239566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 2249566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &cdim)); 2251baa6e33SBarry Smith if (v0Ex) PetscCall(CheckFEMGeometry(dm, cell, cdim, v0Ex, JEx, invJEx, detJEx)); 226c4762a1bSJed Brown if (dim == depth && centroidEx) { 2279566063dSJacob Faibussowitsch PetscCall(CheckFVMGeometry(dm, cell, cdim, centroidEx, normalEx, volEx)); 2289566063dSJacob Faibussowitsch PetscCall(CheckGaussLaw(dm, cell)); 229c4762a1bSJed Brown if (faceCentroidEx) { 2309566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, cell, &coneSize)); 2319566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, cell, &cone)); 232c4762a1bSJed Brown for (c = 0; c < coneSize; ++c) { 2339566063dSJacob Faibussowitsch PetscCall(CheckFVMGeometry(dm, cone[c], dim, &faceCentroidEx[c*dim], &faceNormalEx[c*dim], faceVolEx[c])); 234c4762a1bSJed Brown } 235c4762a1bSJed Brown } 236c4762a1bSJed Brown } 237c4762a1bSJed Brown if (transform) { 238c4762a1bSJed Brown Vec coordinates; 239c4762a1bSJed Brown PetscSection coordSection; 240c4762a1bSJed Brown PetscScalar *coords = NULL, *origCoords, *newCoords; 241c4762a1bSJed Brown PetscReal *v0ExT, *JExT, *invJExT, detJExT=0, *centroidExT, *normalExT, volExT=0; 242c4762a1bSJed Brown PetscReal *faceCentroidExT, *faceNormalExT, faceVolExT; 243c4762a1bSJed Brown PetscRandom r, ang, ang2; 244c4762a1bSJed Brown PetscInt coordSize, numCorners, t; 245c4762a1bSJed Brown 2469566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 2479566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm, &coordSection)); 2489566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dm, coordSection, coordinates, cell, &coordSize, &coords)); 2499566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(coordSize, &origCoords, coordSize, &newCoords)); 2509566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(cdim, &v0ExT, cdim*cdim, &JExT, cdim*cdim, &invJExT, cdim, ¢roidExT, cdim, &normalExT)); 2519566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(cdim, &faceCentroidExT, cdim, &faceNormalExT)); 252c4762a1bSJed Brown for (c = 0; c < coordSize; ++c) origCoords[c] = coords[c]; 2539566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dm, coordSection, coordinates, cell, &coordSize, &coords)); 254c4762a1bSJed Brown numCorners = coordSize/cdim; 255c4762a1bSJed Brown 2569566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &r)); 2579566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(r)); 2589566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(r, 0.0, 10.0)); 2599566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &ang)); 2609566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(ang)); 2619566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(ang, 0.0, 2*PETSC_PI)); 2629566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &ang2)); 2639566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(ang2)); 2649566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(ang2, 0.0, PETSC_PI)); 265c4762a1bSJed Brown for (t = 0; t < 1; ++t) { 266c4762a1bSJed Brown PetscScalar trans[3]; 267c4762a1bSJed Brown PetscReal R[9], rot[3], rotM[9]; 268c4762a1bSJed Brown PetscReal scale, phi, theta, psi = 0.0, norm; 269c4762a1bSJed Brown PetscInt d, e, f, p; 270c4762a1bSJed Brown 271c4762a1bSJed Brown for (c = 0; c < coordSize; ++c) newCoords[c] = origCoords[c]; 2729566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(r, &scale)); 2739566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(ang, &phi)); 2749566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(ang2, &theta)); 275c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 2769566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(r, &trans[d])); 277c4762a1bSJed Brown } 278c4762a1bSJed Brown switch (cdim) { 279c4762a1bSJed Brown case 2: 280c4762a1bSJed Brown R[0] = PetscCosReal(phi); R[1] = -PetscSinReal(phi); 281c4762a1bSJed Brown R[2] = PetscSinReal(phi); R[3] = PetscCosReal(phi); 282c4762a1bSJed Brown break; 283c4762a1bSJed Brown case 3: 284c4762a1bSJed Brown { 285c4762a1bSJed Brown const PetscReal ct = PetscCosReal(theta), st = PetscSinReal(theta); 286c4762a1bSJed Brown const PetscReal cp = PetscCosReal(phi), sp = PetscSinReal(phi); 287c4762a1bSJed Brown const PetscReal cs = PetscCosReal(psi), ss = PetscSinReal(psi); 288c4762a1bSJed Brown R[0] = ct*cs; R[1] = sp*st*cs - cp*ss; R[2] = sp*ss + cp*st*cs; 289c4762a1bSJed Brown R[3] = ct*ss; R[4] = cp*cs + sp*st*ss; R[5] = cp*st*ss - sp*cs; 290c4762a1bSJed Brown R[6] = -st; R[7] = sp*ct; R[8] = cp*ct; 291c4762a1bSJed Brown break; 292c4762a1bSJed Brown } 29363a3b9bcSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Invalid coordinate dimension %" PetscInt_FMT, cdim); 294c4762a1bSJed Brown } 295c4762a1bSJed Brown if (v0Ex) { 296c4762a1bSJed Brown detJExT = detJEx; 297c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 298c4762a1bSJed Brown v0ExT[d] = v0Ex[d]; 299c4762a1bSJed Brown for (e = 0; e < cdim; ++e) { 300c4762a1bSJed Brown JExT[d*cdim+e] = JEx[d*cdim+e]; 301c4762a1bSJed Brown invJExT[d*cdim+e] = invJEx[d*cdim+e]; 302c4762a1bSJed Brown } 303c4762a1bSJed Brown } 304c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 305c4762a1bSJed Brown v0ExT[d] *= scale; 306c4762a1bSJed Brown v0ExT[d] += PetscRealPart(trans[d]); 307c4762a1bSJed Brown /* Only scale dimensions in the manifold */ 308c4762a1bSJed Brown for (e = 0; e < dim; ++e) { 309c4762a1bSJed Brown JExT[d*cdim+e] *= scale; 310c4762a1bSJed Brown invJExT[d*cdim+e] /= scale; 311c4762a1bSJed Brown } 312c4762a1bSJed Brown if (d < dim) detJExT *= scale; 313c4762a1bSJed Brown } 314c4762a1bSJed Brown /* Do scaling and translation before rotation, so that we can leave out the normal dimension for lower dimensional manifolds */ 315c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 316c4762a1bSJed Brown for (e = 0, rot[d] = 0.0; e < cdim; ++e) { 317c4762a1bSJed Brown rot[d] += R[d*cdim+e] * v0ExT[e]; 318c4762a1bSJed Brown } 319c4762a1bSJed Brown } 320c4762a1bSJed Brown for (d = 0; d < cdim; ++d) v0ExT[d] = rot[d]; 321c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 322c4762a1bSJed Brown for (e = 0; e < cdim; ++e) { 323c4762a1bSJed Brown for (f = 0, rotM[d*cdim+e] = 0.0; f < cdim; ++f) { 324c4762a1bSJed Brown rotM[d*cdim+e] += R[d*cdim+f] * JExT[f*cdim+e]; 325c4762a1bSJed Brown } 326c4762a1bSJed Brown } 327c4762a1bSJed Brown } 328c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 329c4762a1bSJed Brown for (e = 0; e < cdim; ++e) { 330c4762a1bSJed Brown JExT[d*cdim+e] = rotM[d*cdim+e]; 331c4762a1bSJed Brown } 332c4762a1bSJed Brown } 333c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 334c4762a1bSJed Brown for (e = 0; e < cdim; ++e) { 335c4762a1bSJed Brown for (f = 0, rotM[d*cdim+e] = 0.0; f < cdim; ++f) { 336c4762a1bSJed Brown rotM[d*cdim+e] += invJExT[d*cdim+f] * R[e*cdim+f]; 337c4762a1bSJed Brown } 338c4762a1bSJed Brown } 339c4762a1bSJed Brown } 340c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 341c4762a1bSJed Brown for (e = 0; e < cdim; ++e) { 342c4762a1bSJed Brown invJExT[d*cdim+e] = rotM[d*cdim+e]; 343c4762a1bSJed Brown } 344c4762a1bSJed Brown } 345c4762a1bSJed Brown } 346c4762a1bSJed Brown if (centroidEx) { 347c4762a1bSJed Brown volExT = volEx; 348c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 349c4762a1bSJed Brown centroidExT[d] = centroidEx[d]; 350c4762a1bSJed Brown normalExT[d] = normalEx[d]; 351c4762a1bSJed Brown } 352c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 353c4762a1bSJed Brown centroidExT[d] *= scale; 354c4762a1bSJed Brown centroidExT[d] += PetscRealPart(trans[d]); 355c4762a1bSJed Brown normalExT[d] /= scale; 356c4762a1bSJed Brown /* Only scale dimensions in the manifold */ 357c4762a1bSJed Brown if (d < dim) volExT *= scale; 358c4762a1bSJed Brown } 359c4762a1bSJed Brown /* Do scaling and translation before rotation, so that we can leave out the normal dimension for lower dimensional manifolds */ 360c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 361c4762a1bSJed Brown for (e = 0, rot[d] = 0.0; e < cdim; ++e) { 362c4762a1bSJed Brown rot[d] += R[d*cdim+e] * centroidExT[e]; 363c4762a1bSJed Brown } 364c4762a1bSJed Brown } 365c4762a1bSJed Brown for (d = 0; d < cdim; ++d) centroidExT[d] = rot[d]; 366c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 367c4762a1bSJed Brown for (e = 0, rot[d] = 0.0; e < cdim; ++e) { 368c4762a1bSJed Brown rot[d] += R[d*cdim+e] * normalExT[e]; 369c4762a1bSJed Brown } 370c4762a1bSJed Brown } 371c4762a1bSJed Brown for (d = 0; d < cdim; ++d) normalExT[d] = rot[d]; 372c4762a1bSJed Brown for (d = 0, norm = 0.0; d < cdim; ++d) norm += PetscSqr(normalExT[d]); 373c4762a1bSJed Brown norm = PetscSqrtReal(norm); 374700a43edSPierre Jolivet if (norm != 0.) for (d = 0; d < cdim; ++d) normalExT[d] /= norm; 375c4762a1bSJed Brown } 376c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 377c4762a1bSJed Brown for (p = 0; p < numCorners; ++p) { 378c4762a1bSJed Brown newCoords[p*cdim+d] *= scale; 379c4762a1bSJed Brown newCoords[p*cdim+d] += trans[d]; 380c4762a1bSJed Brown } 381c4762a1bSJed Brown } 382c4762a1bSJed Brown for (p = 0; p < numCorners; ++p) { 383c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 384c4762a1bSJed Brown for (e = 0, rot[d] = 0.0; e < cdim; ++e) { 385c4762a1bSJed Brown rot[d] += R[d*cdim+e] * PetscRealPart(newCoords[p*cdim+e]); 386c4762a1bSJed Brown } 387c4762a1bSJed Brown } 388c4762a1bSJed Brown for (d = 0; d < cdim; ++d) newCoords[p*cdim+d] = rot[d]; 389c4762a1bSJed Brown } 390c4762a1bSJed Brown 3919566063dSJacob Faibussowitsch PetscCall(ChangeCoordinates(dm, cdim, newCoords)); 3921baa6e33SBarry Smith if (v0Ex) PetscCall(CheckFEMGeometry(dm, 0, cdim, v0ExT, JExT, invJExT, detJExT)); 393c4762a1bSJed Brown if (dim == depth && centroidEx) { 3949566063dSJacob Faibussowitsch PetscCall(CheckFVMGeometry(dm, cell, cdim, centroidExT, normalExT, volExT)); 3959566063dSJacob Faibussowitsch PetscCall(CheckGaussLaw(dm, cell)); 396c4762a1bSJed Brown if (faceCentroidEx) { 3979566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, cell, &coneSize)); 3989566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, cell, &cone)); 399c4762a1bSJed Brown for (c = 0; c < coneSize; ++c) { 400c4762a1bSJed Brown PetscInt off = c*cdim; 401c4762a1bSJed Brown 402c4762a1bSJed Brown faceVolExT = faceVolEx[c]; 403c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 404c4762a1bSJed Brown faceCentroidExT[d] = faceCentroidEx[off+d]; 405c4762a1bSJed Brown faceNormalExT[d] = faceNormalEx[off+d]; 406c4762a1bSJed Brown } 407c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 408c4762a1bSJed Brown faceCentroidExT[d] *= scale; 409c4762a1bSJed Brown faceCentroidExT[d] += PetscRealPart(trans[d]); 410c4762a1bSJed Brown faceNormalExT[d] /= scale; 411c4762a1bSJed Brown /* Only scale dimensions in the manifold */ 412c4762a1bSJed Brown if (d < dim-1) { 413c4762a1bSJed Brown faceVolExT *= scale; 414c4762a1bSJed Brown } 415c4762a1bSJed Brown } 416c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 417c4762a1bSJed Brown for (e = 0, rot[d] = 0.0; e < cdim; ++e) { 418c4762a1bSJed Brown rot[d] += R[d*cdim+e] * faceCentroidExT[e]; 419c4762a1bSJed Brown } 420c4762a1bSJed Brown } 421c4762a1bSJed Brown for (d = 0; d < cdim; ++d) faceCentroidExT[d] = rot[d]; 422c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 423c4762a1bSJed Brown for (e = 0, rot[d] = 0.0; e < cdim; ++e) { 424c4762a1bSJed Brown rot[d] += R[d*cdim+e] * faceNormalExT[e]; 425c4762a1bSJed Brown } 426c4762a1bSJed Brown } 427c4762a1bSJed Brown for (d = 0; d < cdim; ++d) faceNormalExT[d] = rot[d]; 428c4762a1bSJed Brown for (d = 0, norm = 0.0; d < cdim; ++d) norm += PetscSqr(faceNormalExT[d]); 429c4762a1bSJed Brown norm = PetscSqrtReal(norm); 430c4762a1bSJed Brown for (d = 0; d < cdim; ++d) faceNormalExT[d] /= norm; 431c4762a1bSJed Brown 4329566063dSJacob Faibussowitsch PetscCall(CheckFVMGeometry(dm, cone[c], cdim, faceCentroidExT, faceNormalExT, faceVolExT)); 433c4762a1bSJed Brown } 434c4762a1bSJed Brown } 435c4762a1bSJed Brown } 436c4762a1bSJed Brown } 4379566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&r)); 4389566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&ang)); 4399566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&ang2)); 4409566063dSJacob Faibussowitsch PetscCall(PetscFree2(origCoords, newCoords)); 4419566063dSJacob Faibussowitsch PetscCall(PetscFree5(v0ExT, JExT, invJExT, centroidExT, normalExT)); 4429566063dSJacob Faibussowitsch PetscCall(PetscFree2(faceCentroidExT, faceNormalExT)); 443c4762a1bSJed Brown } 444c4762a1bSJed Brown PetscFunctionReturn(0); 445c4762a1bSJed Brown } 446c4762a1bSJed Brown 4474f99dae5SMatthew G. Knepley static PetscErrorCode TestTriangle(MPI_Comm comm, PetscBool transform) 448c4762a1bSJed Brown { 449c4762a1bSJed Brown DM dm; 450c4762a1bSJed Brown 451c4762a1bSJed Brown PetscFunctionBegin; 4529566063dSJacob Faibussowitsch PetscCall(DMPlexCreateReferenceCell(comm, DM_POLYTOPE_TRIANGLE, &dm)); 4539566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 454c4762a1bSJed Brown /* Check reference geometry: determinant is scaled by reference volume (2.0) */ 455c4762a1bSJed Brown { 456c4762a1bSJed Brown PetscReal v0Ex[2] = {-1.0, -1.0}; 457c4762a1bSJed Brown PetscReal JEx[4] = {1.0, 0.0, 0.0, 1.0}; 458c4762a1bSJed Brown PetscReal invJEx[4] = {1.0, 0.0, 0.0, 1.0}; 459c4762a1bSJed Brown PetscReal detJEx = 1.0; 460c4762a1bSJed Brown PetscReal centroidEx[2] = {-((PetscReal)1.)/((PetscReal)3.), -((PetscReal)1.)/((PetscReal)3.)}; 461c4762a1bSJed Brown PetscReal normalEx[2] = {0.0, 0.0}; 462c4762a1bSJed Brown PetscReal volEx = 2.0; 463c4762a1bSJed Brown 4649566063dSJacob Faibussowitsch PetscCall(CheckCell(dm, 0, transform, v0Ex, JEx, invJEx, detJEx, centroidEx, normalEx, volEx, NULL, NULL, NULL)); 465c4762a1bSJed Brown } 466c4762a1bSJed Brown /* Move to 3D: Check reference geometry: determinant is scaled by reference volume (2.0) */ 467c4762a1bSJed Brown { 468c4762a1bSJed Brown PetscScalar vertexCoords[9] = {-1.0, -1.0, 0.0, 1.0, -1.0, 0.0, -1.0, 1.0, 0.0}; 469c4762a1bSJed Brown PetscReal v0Ex[3] = {-1.0, -1.0, 0.0}; 470c4762a1bSJed Brown PetscReal JEx[9] = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0}; 471c4762a1bSJed Brown PetscReal invJEx[9] = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0}; 472c4762a1bSJed Brown PetscReal detJEx = 1.0; 473c4762a1bSJed Brown PetscReal centroidEx[3] = {-((PetscReal)1.)/((PetscReal)3.), -((PetscReal)1.)/((PetscReal)3.), 0.0}; 474c4762a1bSJed Brown PetscReal normalEx[3] = {0.0, 0.0, 1.0}; 475c4762a1bSJed Brown PetscReal volEx = 2.0; 476c4762a1bSJed Brown 4779566063dSJacob Faibussowitsch PetscCall(ChangeCoordinates(dm, 3, vertexCoords)); 4789566063dSJacob Faibussowitsch PetscCall(CheckCell(dm, 0, transform, v0Ex, JEx, invJEx, detJEx, centroidEx, normalEx, volEx, NULL, NULL, NULL)); 479c4762a1bSJed Brown } 480c4762a1bSJed Brown /* Cleanup */ 4819566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 482c4762a1bSJed Brown PetscFunctionReturn(0); 483c4762a1bSJed Brown } 484c4762a1bSJed Brown 4854f99dae5SMatthew G. Knepley static PetscErrorCode TestQuadrilateral(MPI_Comm comm, PetscBool transform) 486c4762a1bSJed Brown { 487c4762a1bSJed Brown DM dm; 488c4762a1bSJed Brown 489c4762a1bSJed Brown PetscFunctionBegin; 4909566063dSJacob Faibussowitsch PetscCall(DMPlexCreateReferenceCell(comm, DM_POLYTOPE_QUADRILATERAL, &dm)); 4919566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 492c4762a1bSJed Brown /* Check reference geometry: determinant is scaled by reference volume (2.0) */ 493c4762a1bSJed Brown { 494c4762a1bSJed Brown PetscReal v0Ex[2] = {-1.0, -1.0}; 495c4762a1bSJed Brown PetscReal JEx[4] = {1.0, 0.0, 0.0, 1.0}; 496c4762a1bSJed Brown PetscReal invJEx[4] = {1.0, 0.0, 0.0, 1.0}; 497c4762a1bSJed Brown PetscReal detJEx = 1.0; 498c4762a1bSJed Brown PetscReal centroidEx[2] = {0.0, 0.0}; 499c4762a1bSJed Brown PetscReal normalEx[2] = {0.0, 0.0}; 500c4762a1bSJed Brown PetscReal volEx = 4.0; 501c4762a1bSJed Brown 5029566063dSJacob Faibussowitsch PetscCall(CheckCell(dm, 0, transform, v0Ex, JEx, invJEx, detJEx, centroidEx, normalEx, volEx, NULL, NULL, NULL)); 503c4762a1bSJed Brown } 504c4762a1bSJed Brown /* Move to 3D: Check reference geometry: determinant is scaled by reference volume (4.0) */ 505c4762a1bSJed Brown { 506c4762a1bSJed 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}; 507c4762a1bSJed Brown PetscReal v0Ex[3] = {-1.0, -1.0, 0.0}; 508c4762a1bSJed Brown PetscReal JEx[9] = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0}; 509c4762a1bSJed Brown PetscReal invJEx[9] = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0}; 510c4762a1bSJed Brown PetscReal detJEx = 1.0; 511c4762a1bSJed Brown PetscReal centroidEx[3] = {0.0, 0.0, 0.0}; 512c4762a1bSJed Brown PetscReal normalEx[3] = {0.0, 0.0, 1.0}; 513c4762a1bSJed Brown PetscReal volEx = 4.0; 514c4762a1bSJed Brown 5159566063dSJacob Faibussowitsch PetscCall(ChangeCoordinates(dm, 3, vertexCoords)); 5169566063dSJacob Faibussowitsch PetscCall(CheckCell(dm, 0, transform, v0Ex, JEx, invJEx, detJEx, centroidEx, normalEx, volEx, NULL, NULL, NULL)); 517c4762a1bSJed Brown } 518c4762a1bSJed Brown /* Cleanup */ 5199566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 520c4762a1bSJed Brown PetscFunctionReturn(0); 521c4762a1bSJed Brown } 522c4762a1bSJed Brown 5234f99dae5SMatthew G. Knepley static PetscErrorCode TestTetrahedron(MPI_Comm comm, PetscBool transform) 524c4762a1bSJed Brown { 525c4762a1bSJed Brown DM dm; 526c4762a1bSJed Brown 527c4762a1bSJed Brown PetscFunctionBegin; 5289566063dSJacob Faibussowitsch PetscCall(DMPlexCreateReferenceCell(comm, DM_POLYTOPE_TETRAHEDRON, &dm)); 5299566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 530c4762a1bSJed Brown /* Check reference geometry: determinant is scaled by reference volume (4/3) */ 531c4762a1bSJed Brown { 532c4762a1bSJed Brown PetscReal v0Ex[3] = {-1.0, -1.0, -1.0}; 533c4762a1bSJed Brown PetscReal JEx[9] = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0}; 534c4762a1bSJed Brown PetscReal invJEx[9] = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0}; 535c4762a1bSJed Brown PetscReal detJEx = 1.0; 536c4762a1bSJed Brown PetscReal centroidEx[3] = {-0.5, -0.5, -0.5}; 537c4762a1bSJed Brown PetscReal normalEx[3] = {0.0, 0.0, 0.0}; 538c4762a1bSJed Brown PetscReal volEx = (PetscReal)4.0/(PetscReal)3.0; 539c4762a1bSJed Brown 5409566063dSJacob Faibussowitsch PetscCall(CheckCell(dm, 0, transform, v0Ex, JEx, invJEx, detJEx, centroidEx, normalEx, volEx, NULL, NULL, NULL)); 541c4762a1bSJed Brown } 542c4762a1bSJed Brown /* Cleanup */ 5439566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 544c4762a1bSJed Brown PetscFunctionReturn(0); 545c4762a1bSJed Brown } 546c4762a1bSJed Brown 5474f99dae5SMatthew G. Knepley static PetscErrorCode TestHexahedron(MPI_Comm comm, PetscBool transform) 548c4762a1bSJed Brown { 549c4762a1bSJed Brown DM dm; 550c4762a1bSJed Brown 551c4762a1bSJed Brown PetscFunctionBegin; 5529566063dSJacob Faibussowitsch PetscCall(DMPlexCreateReferenceCell(comm, DM_POLYTOPE_HEXAHEDRON, &dm)); 5539566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 554c4762a1bSJed Brown /* Check reference geometry: determinant is scaled by reference volume 8.0 */ 555c4762a1bSJed Brown { 556c4762a1bSJed Brown PetscReal v0Ex[3] = {-1.0, -1.0, -1.0}; 557c4762a1bSJed Brown PetscReal JEx[9] = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0}; 558c4762a1bSJed Brown PetscReal invJEx[9] = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0}; 559c4762a1bSJed Brown PetscReal detJEx = 1.0; 560c4762a1bSJed Brown PetscReal centroidEx[3] = {0.0, 0.0, 0.0}; 561c4762a1bSJed Brown PetscReal normalEx[3] = {0.0, 0.0, 0.0}; 562c4762a1bSJed Brown PetscReal volEx = 8.0; 563c4762a1bSJed Brown 5649566063dSJacob Faibussowitsch PetscCall(CheckCell(dm, 0, transform, v0Ex, JEx, invJEx, detJEx, centroidEx, normalEx, volEx, NULL, NULL, NULL)); 565c4762a1bSJed Brown } 566c4762a1bSJed Brown /* Cleanup */ 5679566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 568c4762a1bSJed Brown PetscFunctionReturn(0); 569c4762a1bSJed Brown } 570c4762a1bSJed Brown 5714f99dae5SMatthew G. Knepley static PetscErrorCode TestHexahedronCurved(MPI_Comm comm) 572c4762a1bSJed Brown { 573c4762a1bSJed Brown DM dm; 5744f99dae5SMatthew 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, 5754f99dae5SMatthew 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}; 576c4762a1bSJed Brown 577c4762a1bSJed Brown PetscFunctionBegin; 5789566063dSJacob Faibussowitsch PetscCall(DMPlexCreateReferenceCell(comm, DM_POLYTOPE_HEXAHEDRON, &dm)); 5799566063dSJacob Faibussowitsch PetscCall(ChangeCoordinates(dm, 3, coords)); 5809566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 5814f99dae5SMatthew G. Knepley { 5824f99dae5SMatthew G. Knepley PetscReal centroidEx[3] = {0.0, 0.0, 0.016803278688524603}; 5834f99dae5SMatthew G. Knepley PetscReal normalEx[3] = {0.0, 0.0, 0.0}; 5844f99dae5SMatthew G. Knepley PetscReal volEx = 8.1333333333333346; 5854f99dae5SMatthew G. Knepley 5869566063dSJacob Faibussowitsch PetscCall(CheckCell(dm, 0, PETSC_FALSE, NULL, NULL, NULL, 0.0, centroidEx, normalEx, volEx, NULL, NULL, NULL)); 587c4762a1bSJed Brown } 5889566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 5894f99dae5SMatthew G. Knepley PetscFunctionReturn(0); 5904f99dae5SMatthew G. Knepley } 5914f99dae5SMatthew G. Knepley 5924f99dae5SMatthew G. Knepley /* This wedge is a tensor product cell, rather than a normal wedge */ 5934f99dae5SMatthew G. Knepley static PetscErrorCode TestWedge(MPI_Comm comm, PetscBool transform) 5944f99dae5SMatthew G. Knepley { 5954f99dae5SMatthew G. Knepley DM dm; 5964f99dae5SMatthew G. Knepley 5974f99dae5SMatthew G. Knepley PetscFunctionBegin; 5989566063dSJacob Faibussowitsch PetscCall(DMPlexCreateReferenceCell(comm, DM_POLYTOPE_TRI_PRISM_TENSOR, &dm)); 5999566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 600c4762a1bSJed Brown /* Check reference geometry: determinant is scaled by reference volume 4.0 */ 601c4762a1bSJed Brown { 602c4762a1bSJed Brown #if 0 603c4762a1bSJed Brown /* FEM geometry is not functional for wedges */ 604c4762a1bSJed Brown PetscReal v0Ex[3] = {-1.0, -1.0, -1.0}; 605c4762a1bSJed Brown PetscReal JEx[9] = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0}; 606c4762a1bSJed Brown PetscReal invJEx[9] = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0}; 607c4762a1bSJed Brown PetscReal detJEx = 1.0; 608c4762a1bSJed Brown #endif 609c4762a1bSJed Brown 6104f99dae5SMatthew G. Knepley { 611c4762a1bSJed Brown PetscReal centroidEx[3] = {-((PetscReal)1.)/((PetscReal)3.), -((PetscReal)1.)/((PetscReal)3.), 0.0}; 612c4762a1bSJed Brown PetscReal normalEx[3] = {0.0, 0.0, 0.0}; 613c4762a1bSJed Brown PetscReal volEx = 4.0; 614c4762a1bSJed Brown PetscReal faceVolEx[5] = {2.0, 2.0, 4.0, PETSC_SQRT2*4.0, 4.0}; 615c4762a1bSJed 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}; 616c4762a1bSJed Brown PetscReal faceCentroidEx[15] = {-((PetscReal)1.)/((PetscReal)3.), -((PetscReal)1.)/((PetscReal)3.), -1.0, 617c4762a1bSJed Brown -((PetscReal)1.)/((PetscReal)3.), -((PetscReal)1.)/((PetscReal)3.), 1.0, 618c4762a1bSJed Brown 0.0, -1.0, 0.0, 0.0, 0.0, 0.0, -1.0, 0.0, 0.0}; 619c4762a1bSJed Brown 6209566063dSJacob Faibussowitsch PetscCall(CheckCell(dm, 0, transform, NULL, NULL, NULL, 0.0, centroidEx, normalEx, volEx, faceCentroidEx, faceNormalEx, faceVolEx)); 621c4762a1bSJed Brown } 622c4762a1bSJed Brown } 623c4762a1bSJed Brown /* Cleanup */ 6249566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 625c4762a1bSJed Brown PetscFunctionReturn(0); 626c4762a1bSJed Brown } 627c4762a1bSJed Brown 628c4762a1bSJed Brown int main(int argc, char **argv) 629c4762a1bSJed Brown { 630c4762a1bSJed Brown AppCtx user; 631c4762a1bSJed Brown 632*327415f7SBarry Smith PetscFunctionBeginUser; 6339566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL,help)); 6349566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 635c4762a1bSJed Brown if (user.runType == RUN_REFERENCE) { 6369566063dSJacob Faibussowitsch PetscCall(TestTriangle(PETSC_COMM_SELF, user.transform)); 6379566063dSJacob Faibussowitsch PetscCall(TestQuadrilateral(PETSC_COMM_SELF, user.transform)); 6389566063dSJacob Faibussowitsch PetscCall(TestTetrahedron(PETSC_COMM_SELF, user.transform)); 6399566063dSJacob Faibussowitsch PetscCall(TestHexahedron(PETSC_COMM_SELF, user.transform)); 6409566063dSJacob Faibussowitsch PetscCall(TestWedge(PETSC_COMM_SELF, user.transform)); 6414f99dae5SMatthew G. Knepley } else if (user.runType == RUN_HEX_CURVED) { 6429566063dSJacob Faibussowitsch PetscCall(TestHexahedronCurved(PETSC_COMM_SELF)); 643c4762a1bSJed Brown } else if (user.runType == RUN_FILE) { 644c4762a1bSJed Brown PetscInt dim, cStart, cEnd, c; 645c4762a1bSJed Brown 6469566063dSJacob Faibussowitsch PetscCall(DMGetDimension(user.dm, &dim)); 6479566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(user.dm, 0, &cStart, &cEnd)); 648c4762a1bSJed Brown for (c = 0; c < cEnd-cStart; ++c) { 6494f99dae5SMatthew G. Knepley PetscReal *v0 = user.v0 ? &user.v0[c*dim] : NULL; 6504f99dae5SMatthew G. Knepley PetscReal *J = user.J ? &user.J[c*dim*dim] : NULL; 6514f99dae5SMatthew G. Knepley PetscReal *invJ = user.invJ ? &user.invJ[c*dim*dim] : NULL; 6524f99dae5SMatthew G. Knepley PetscReal detJ = user.detJ ? user.detJ[c] : 0.0; 6534f99dae5SMatthew G. Knepley PetscReal *centroid = user.centroid ? &user.centroid[c*dim] : NULL; 6544f99dae5SMatthew G. Knepley PetscReal *normal = user.normal ? &user.normal[c*dim] : NULL; 6554f99dae5SMatthew G. Knepley PetscReal vol = user.vol ? user.vol[c] : 0.0; 6564f99dae5SMatthew G. Knepley 6579566063dSJacob Faibussowitsch PetscCall(CheckCell(user.dm, c+cStart, PETSC_FALSE, v0, J, invJ, detJ, centroid, normal, vol, NULL, NULL, NULL)); 658c4762a1bSJed Brown } 6599566063dSJacob Faibussowitsch PetscCall(PetscFree4(user.v0,user.J,user.invJ,user.detJ)); 6609566063dSJacob Faibussowitsch PetscCall(PetscFree(user.centroid)); 6619566063dSJacob Faibussowitsch PetscCall(PetscFree(user.normal)); 6629566063dSJacob Faibussowitsch PetscCall(PetscFree(user.vol)); 6639566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.dm)); 664c4762a1bSJed Brown } else if (user.runType == RUN_DISPLAY) { 665c4762a1bSJed Brown DM gdm, dmCell; 666c4762a1bSJed Brown Vec cellgeom, facegeom; 667c4762a1bSJed Brown const PetscScalar *cgeom; 668c4762a1bSJed Brown PetscInt dim, d, cStart, cEnd, cEndInterior, c; 669c4762a1bSJed Brown 6709566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(user.dm, &dim)); 6719566063dSJacob Faibussowitsch PetscCall(DMPlexConstructGhostCells(user.dm, NULL, NULL, &gdm)); 672c4762a1bSJed Brown if (gdm) { 6739566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.dm)); 674c4762a1bSJed Brown user.dm = gdm; 675c4762a1bSJed Brown } 6769566063dSJacob Faibussowitsch PetscCall(DMPlexComputeGeometryFVM(user.dm, &cellgeom, &facegeom)); 6779566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(user.dm, 0, &cStart, &cEnd)); 6789566063dSJacob Faibussowitsch PetscCall(DMPlexGetGhostCellStratum(user.dm, &cEndInterior, NULL)); 679c4762a1bSJed Brown if (cEndInterior >= 0) cEnd = cEndInterior; 6809566063dSJacob Faibussowitsch PetscCall(VecGetDM(cellgeom, &dmCell)); 6819566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(cellgeom, &cgeom)); 682c4762a1bSJed Brown for (c = 0; c < cEnd-cStart; ++c) { 683c4762a1bSJed Brown PetscFVCellGeom *cg; 684c4762a1bSJed Brown 6859566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRead(dmCell, c, cgeom, &cg)); 68663a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Cell %4" PetscInt_FMT ": Centroid (", c)); 687c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 6889566063dSJacob Faibussowitsch if (d > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", ")); 68963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "%12.2g", (double)cg->centroid[d])); 690c4762a1bSJed Brown } 69163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, ") Vol %12.2g\n", (double)cg->volume)); 692c4762a1bSJed Brown } 6939566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(cellgeom, &cgeom)); 6949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cellgeom)); 6959566063dSJacob Faibussowitsch PetscCall(VecDestroy(&facegeom)); 6969566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.dm)); 697c4762a1bSJed Brown } 6989566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 699b122ec5aSJacob Faibussowitsch return 0; 700c4762a1bSJed Brown } 701c4762a1bSJed Brown 702c4762a1bSJed Brown /*TEST 703c4762a1bSJed Brown 704c4762a1bSJed Brown test: 7054f99dae5SMatthew G. Knepley suffix: 1 706c4762a1bSJed Brown args: -dm_view ascii::ascii_info_detail 707c4762a1bSJed Brown test: 708c4762a1bSJed Brown suffix: 2 7094f99dae5SMatthew G. Knepley args: -run_type hex_curved 710c4762a1bSJed Brown test: 711c4762a1bSJed Brown suffix: 3 7124f99dae5SMatthew G. Knepley args: -transform 713c4762a1bSJed Brown test: 714c4762a1bSJed Brown suffix: 4 715c4762a1bSJed Brown requires: exodusii 7164f99dae5SMatthew 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 717c4762a1bSJed Brown test: 718c4762a1bSJed Brown suffix: 5 7194f99dae5SMatthew 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 7209bf2564aSMatt McGurn test: 7219bf2564aSMatt McGurn suffix: 6 7229bf2564aSMatt 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 723c4762a1bSJed Brown TEST*/ 724