1c4762a1bSJed Brown static char help[] = "Tests for cell geometry\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdmplex.h> 4c4762a1bSJed Brown 59371c9d4SSatish Balay typedef enum { 69371c9d4SSatish Balay RUN_REFERENCE, 79371c9d4SSatish Balay RUN_HEX_CURVED, 89371c9d4SSatish Balay RUN_FILE, 99371c9d4SSatish Balay RUN_DISPLAY 109371c9d4SSatish Balay } RunType; 11c4762a1bSJed Brown 12c4762a1bSJed Brown typedef struct { 13c4762a1bSJed Brown DM dm; 14c4762a1bSJed Brown RunType runType; /* Type of mesh to use */ 15c4762a1bSJed Brown PetscBool transform; /* Use random coordinate transformations */ 16c4762a1bSJed Brown /* Data for input meshes */ 17c4762a1bSJed Brown PetscReal *v0, *J, *invJ, *detJ; /* FEM data */ 18c4762a1bSJed Brown PetscReal *centroid, *normal, *vol; /* FVM data */ 19c4762a1bSJed Brown } AppCtx; 20c4762a1bSJed Brown 21d71ae5a4SJacob Faibussowitsch static PetscErrorCode ReadMesh(MPI_Comm comm, AppCtx *user, DM *dm) 22d71ae5a4SJacob Faibussowitsch { 23c4762a1bSJed Brown PetscFunctionBegin; 249566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 259566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 269566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 279566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(*dm, user)); 289566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*dm, "Input Mesh")); 299566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 31c4762a1bSJed Brown } 32c4762a1bSJed Brown 33d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 34d71ae5a4SJacob Faibussowitsch { 354f99dae5SMatthew G. Knepley const char *runTypes[4] = {"reference", "hex_curved", "file", "display"}; 36c4762a1bSJed Brown PetscInt run; 37c4762a1bSJed Brown 38c4762a1bSJed Brown PetscFunctionBeginUser; 39c4762a1bSJed Brown options->runType = RUN_REFERENCE; 40c4762a1bSJed Brown options->transform = PETSC_FALSE; 41c4762a1bSJed Brown 42d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Geometry Test Options", "DMPLEX"); 43c4762a1bSJed Brown run = options->runType; 449566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-run_type", "The run type", "ex8.c", runTypes, 3, runTypes[options->runType], &run, NULL)); 45c4762a1bSJed Brown options->runType = (RunType)run; 469566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-transform", "Use random transforms", "ex8.c", options->transform, &options->transform, NULL)); 47c4762a1bSJed Brown 48c4762a1bSJed Brown if (options->runType == RUN_FILE) { 49c4762a1bSJed Brown PetscInt dim, cStart, cEnd, numCells, n; 509bf2564aSMatt McGurn PetscBool flg, feFlg; 51c4762a1bSJed Brown 529566063dSJacob Faibussowitsch PetscCall(ReadMesh(PETSC_COMM_WORLD, options, &options->dm)); 539566063dSJacob Faibussowitsch PetscCall(DMGetDimension(options->dm, &dim)); 549566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(options->dm, 0, &cStart, &cEnd)); 55c4762a1bSJed Brown numCells = cEnd - cStart; 569566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(numCells * dim, &options->v0, numCells * dim * dim, &options->J, numCells * dim * dim, &options->invJ, numCells, &options->detJ)); 579566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numCells * dim, &options->centroid)); 589566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numCells * dim, &options->normal)); 599566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numCells, &options->vol)); 60c4762a1bSJed Brown n = numCells * dim; 619566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-v0", "Input v0 for each cell", "ex8.c", options->v0, &n, &feFlg)); 6263a3b9bcSJacob 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); 63c4762a1bSJed Brown n = numCells * dim * dim; 649566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-J", "Input Jacobian for each cell", "ex8.c", options->J, &n, &flg)); 6563a3b9bcSJacob 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); 66c4762a1bSJed Brown n = numCells * dim * dim; 679566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-invJ", "Input inverse Jacobian for each cell", "ex8.c", options->invJ, &n, &flg)); 6863a3b9bcSJacob 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); 69c4762a1bSJed Brown n = numCells; 709566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-detJ", "Input Jacobian determinant for each cell", "ex8.c", options->detJ, &n, &flg)); 7163a3b9bcSJacob Faibussowitsch PetscCheck(!flg || n == numCells, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid size of detJ %" PetscInt_FMT " should be %" PetscInt_FMT, n, numCells); 72c4762a1bSJed Brown n = numCells * dim; 734f99dae5SMatthew G. Knepley if (!feFlg) { 749566063dSJacob Faibussowitsch PetscCall(PetscFree4(options->v0, options->J, options->invJ, options->detJ)); 754f99dae5SMatthew G. Knepley options->v0 = options->J = options->invJ = options->detJ = NULL; 764f99dae5SMatthew G. Knepley } 779566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-centroid", "Input centroid for each cell", "ex8.c", options->centroid, &n, &flg)); 7863a3b9bcSJacob 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); 799bf2564aSMatt McGurn if (!flg) { 809566063dSJacob Faibussowitsch PetscCall(PetscFree(options->centroid)); 819bf2564aSMatt McGurn options->centroid = NULL; 829bf2564aSMatt McGurn } 83c4762a1bSJed Brown n = numCells * dim; 849566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-normal", "Input normal for each cell", "ex8.c", options->normal, &n, &flg)); 8563a3b9bcSJacob 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); 869bf2564aSMatt McGurn if (!flg) { 879566063dSJacob Faibussowitsch PetscCall(PetscFree(options->normal)); 889bf2564aSMatt McGurn options->normal = NULL; 899bf2564aSMatt McGurn } 90c4762a1bSJed Brown n = numCells; 919566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-vol", "Input volume for each cell", "ex8.c", options->vol, &n, &flg)); 9263a3b9bcSJacob Faibussowitsch PetscCheck(!flg || n == numCells, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid size of vol %" PetscInt_FMT " should be %" PetscInt_FMT, n, numCells); 939bf2564aSMatt McGurn if (!flg) { 949566063dSJacob Faibussowitsch PetscCall(PetscFree(options->vol)); 959bf2564aSMatt McGurn options->vol = NULL; 964f99dae5SMatthew G. Knepley } 97c4762a1bSJed Brown } else if (options->runType == RUN_DISPLAY) { 989566063dSJacob Faibussowitsch PetscCall(ReadMesh(PETSC_COMM_WORLD, options, &options->dm)); 99c4762a1bSJed Brown } 100d0609cedSBarry Smith PetscOptionsEnd(); 101c4762a1bSJed Brown 1029566063dSJacob Faibussowitsch if (options->transform) PetscCall(PetscPrintf(comm, "Using random transforms\n")); 1033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 104c4762a1bSJed Brown } 105c4762a1bSJed Brown 106d71ae5a4SJacob Faibussowitsch static PetscErrorCode ChangeCoordinates(DM dm, PetscInt spaceDim, PetscScalar vertexCoords[]) 107d71ae5a4SJacob Faibussowitsch { 108c4762a1bSJed Brown PetscSection coordSection; 109c4762a1bSJed Brown Vec coordinates; 110c4762a1bSJed Brown PetscScalar *coords; 111c4762a1bSJed Brown PetscInt vStart, vEnd, v, d, coordSize; 112c4762a1bSJed Brown 113c4762a1bSJed Brown PetscFunctionBegin; 1149566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 1159566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm, &coordSection)); 1169566063dSJacob Faibussowitsch PetscCall(PetscSectionSetNumFields(coordSection, 1)); 1179566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldComponents(coordSection, 0, spaceDim)); 1189566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(coordSection, vStart, vEnd)); 119c4762a1bSJed Brown for (v = vStart; v < vEnd; ++v) { 1209566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(coordSection, v, spaceDim)); 1219566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, spaceDim)); 122c4762a1bSJed Brown } 1239566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(coordSection)); 1249566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize)); 1259566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates)); 1269566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates")); 1279566063dSJacob Faibussowitsch PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE)); 1289566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(coordinates)); 1299566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordinates, &coords)); 130c4762a1bSJed Brown for (v = vStart; v < vEnd; ++v) { 131c4762a1bSJed Brown PetscInt off; 132c4762a1bSJed Brown 1339566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(coordSection, v, &off)); 134ad540459SPierre Jolivet for (d = 0; d < spaceDim; ++d) coords[off + d] = vertexCoords[(v - vStart) * spaceDim + d]; 135c4762a1bSJed Brown } 1369566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coordinates, &coords)); 1379566063dSJacob Faibussowitsch PetscCall(DMSetCoordinateDim(dm, spaceDim)); 1389566063dSJacob Faibussowitsch PetscCall(DMSetCoordinatesLocal(dm, coordinates)); 1399566063dSJacob Faibussowitsch PetscCall(VecDestroy(&coordinates)); 1409566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 1413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 142c4762a1bSJed Brown } 143c4762a1bSJed Brown 144c4762a1bSJed Brown #define RelativeError(a, b) PetscAbs(a - b) / (1.0 + PetscMax(PetscAbs(a), PetscAbs(b))) 145c4762a1bSJed Brown 146d71ae5a4SJacob Faibussowitsch static PetscErrorCode CheckFEMGeometry(DM dm, PetscInt cell, PetscInt spaceDim, PetscReal v0Ex[], PetscReal JEx[], PetscReal invJEx[], PetscReal detJEx) 147d71ae5a4SJacob Faibussowitsch { 148c4762a1bSJed Brown PetscReal v0[3], J[9], invJ[9], detJ; 149c4762a1bSJed Brown PetscInt d, i, j; 150c4762a1bSJed Brown 151c4762a1bSJed Brown PetscFunctionBegin; 1529566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFEM(dm, cell, NULL, v0, J, invJ, &detJ)); 153c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) { 154c4762a1bSJed Brown if (v0[d] != v0Ex[d]) { 155c4762a1bSJed Brown switch (spaceDim) { 156d71ae5a4SJacob Faibussowitsch case 2: 157d71ae5a4SJacob Faibussowitsch 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]); 158d71ae5a4SJacob Faibussowitsch case 3: 159d71ae5a4SJacob Faibussowitsch 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]); 160d71ae5a4SJacob Faibussowitsch default: 161d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid space dimension %" PetscInt_FMT, spaceDim); 162c4762a1bSJed Brown } 163c4762a1bSJed Brown } 164c4762a1bSJed Brown } 165c4762a1bSJed Brown for (i = 0; i < spaceDim; ++i) { 166c4762a1bSJed Brown for (j = 0; j < spaceDim; ++j) { 16763a3b9bcSJacob 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]); 16863a3b9bcSJacob 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]); 169c4762a1bSJed Brown } 170c4762a1bSJed Brown } 171700a43edSPierre 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)); 1723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 173c4762a1bSJed Brown } 174c4762a1bSJed Brown 175d71ae5a4SJacob Faibussowitsch static PetscErrorCode CheckFVMGeometry(DM dm, PetscInt cell, PetscInt spaceDim, PetscReal centroidEx[], PetscReal normalEx[], PetscReal volEx) 176d71ae5a4SJacob Faibussowitsch { 1774f99dae5SMatthew G. Knepley PetscReal tol = PetscMax(10 * PETSC_SMALL, 1e-10); 178c4762a1bSJed Brown PetscReal centroid[3], normal[3], vol; 179c4762a1bSJed Brown PetscInt d; 180c4762a1bSJed Brown 181c4762a1bSJed Brown PetscFunctionBegin; 1829566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFVM(dm, cell, volEx ? &vol : NULL, centroidEx ? centroid : NULL, normalEx ? normal : NULL)); 183c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) { 1849bf2564aSMatt McGurn if (centroidEx) 18563a3b9bcSJacob 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])); 1869371c9d4SSatish Balay if (normalEx) 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]); 187c4762a1bSJed Brown } 1889371c9d4SSatish Balay if (volEx) 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)); 1893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1904f99dae5SMatthew G. Knepley } 1914f99dae5SMatthew G. Knepley 192d71ae5a4SJacob Faibussowitsch static PetscErrorCode CheckGaussLaw(DM dm, PetscInt cell) 193d71ae5a4SJacob Faibussowitsch { 1944f99dae5SMatthew G. Knepley DMPolytopeType ct; 1954f99dae5SMatthew G. Knepley PetscReal tol = PetscMax(10 * PETSC_SMALL, 1e-10); 1964f99dae5SMatthew G. Knepley PetscReal normal[3], integral[3] = {0., 0., 0.}, area; 1974f99dae5SMatthew G. Knepley const PetscInt *cone, *ornt; 1984f99dae5SMatthew G. Knepley PetscInt coneSize, f, dim, cdim, d; 1994f99dae5SMatthew G. Knepley 2004f99dae5SMatthew G. Knepley PetscFunctionBegin; 2019566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 2029566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &cdim)); 2033ba16761SJacob Faibussowitsch if (dim != cdim) PetscFunctionReturn(PETSC_SUCCESS); 2049566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, cell, &ct)); 2053ba16761SJacob Faibussowitsch if (ct == DM_POLYTOPE_TRI_PRISM_TENSOR) PetscFunctionReturn(PETSC_SUCCESS); 2069566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, cell, &coneSize)); 2079566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, cell, &cone)); 2089566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeOrientation(dm, cell, &ornt)); 2094f99dae5SMatthew G. Knepley for (f = 0; f < coneSize; ++f) { 2109bf2564aSMatt McGurn const PetscInt sgn = dim == 1 ? (f == 0 ? -1 : 1) : (ornt[f] < 0 ? -1 : 1); 2114f99dae5SMatthew G. Knepley 2129566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFVM(dm, cone[f], &area, NULL, normal)); 2134f99dae5SMatthew G. Knepley for (d = 0; d < cdim; ++d) integral[d] += sgn * area * normal[d]; 2144f99dae5SMatthew G. Knepley } 2159371c9d4SSatish Balay for (d = 0; d < cdim; ++d) 2169371c9d4SSatish Balay 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]); 2173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 218c4762a1bSJed Brown } 219c4762a1bSJed Brown 220d71ae5a4SJacob Faibussowitsch 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[]) 221d71ae5a4SJacob Faibussowitsch { 222c4762a1bSJed Brown const PetscInt *cone; 223c4762a1bSJed Brown PetscInt coneSize, c; 224c4762a1bSJed Brown PetscInt dim, depth, cdim; 225c4762a1bSJed Brown 226c4762a1bSJed Brown PetscFunctionBegin; 2279566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepth(dm, &depth)); 2289566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 2299566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &cdim)); 2301baa6e33SBarry Smith if (v0Ex) PetscCall(CheckFEMGeometry(dm, cell, cdim, v0Ex, JEx, invJEx, detJEx)); 231c4762a1bSJed Brown if (dim == depth && centroidEx) { 2329566063dSJacob Faibussowitsch PetscCall(CheckFVMGeometry(dm, cell, cdim, centroidEx, normalEx, volEx)); 2339566063dSJacob Faibussowitsch PetscCall(CheckGaussLaw(dm, cell)); 234c4762a1bSJed Brown if (faceCentroidEx) { 2359566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, cell, &coneSize)); 2369566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, cell, &cone)); 23748a46eb9SPierre Jolivet for (c = 0; c < coneSize; ++c) PetscCall(CheckFVMGeometry(dm, cone[c], dim, &faceCentroidEx[c * dim], &faceNormalEx[c * dim], faceVolEx[c])); 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)); 27848a46eb9SPierre Jolivet for (d = 0; d < cdim; ++d) PetscCall(PetscRandomGetValue(r, &trans[d])); 279c4762a1bSJed Brown switch (cdim) { 280c4762a1bSJed Brown case 2: 2819371c9d4SSatish Balay R[0] = PetscCosReal(phi); 2829371c9d4SSatish Balay R[1] = -PetscSinReal(phi); 2839371c9d4SSatish Balay R[2] = PetscSinReal(phi); 2849371c9d4SSatish Balay R[3] = PetscCosReal(phi); 285c4762a1bSJed Brown break; 2869371c9d4SSatish Balay case 3: { 287c4762a1bSJed Brown const PetscReal ct = PetscCosReal(theta), st = PetscSinReal(theta); 288c4762a1bSJed Brown const PetscReal cp = PetscCosReal(phi), sp = PetscSinReal(phi); 289c4762a1bSJed Brown const PetscReal cs = PetscCosReal(psi), ss = PetscSinReal(psi); 2909371c9d4SSatish Balay R[0] = ct * cs; 2919371c9d4SSatish Balay R[1] = sp * st * cs - cp * ss; 2929371c9d4SSatish Balay R[2] = sp * ss + cp * st * cs; 2939371c9d4SSatish Balay R[3] = ct * ss; 2949371c9d4SSatish Balay R[4] = cp * cs + sp * st * ss; 2959371c9d4SSatish Balay R[5] = cp * st * ss - sp * cs; 2969371c9d4SSatish Balay R[6] = -st; 2979371c9d4SSatish Balay R[7] = sp * ct; 2989371c9d4SSatish Balay R[8] = cp * ct; 299c4762a1bSJed Brown break; 300c4762a1bSJed Brown } 301d71ae5a4SJacob Faibussowitsch default: 302d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid coordinate dimension %" PetscInt_FMT, cdim); 303c4762a1bSJed Brown } 304c4762a1bSJed Brown if (v0Ex) { 305c4762a1bSJed Brown detJExT = detJEx; 306c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 307c4762a1bSJed Brown v0ExT[d] = v0Ex[d]; 308c4762a1bSJed Brown for (e = 0; e < cdim; ++e) { 309c4762a1bSJed Brown JExT[d * cdim + e] = JEx[d * cdim + e]; 310c4762a1bSJed Brown invJExT[d * cdim + e] = invJEx[d * cdim + e]; 311c4762a1bSJed Brown } 312c4762a1bSJed Brown } 313c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 314c4762a1bSJed Brown v0ExT[d] *= scale; 315c4762a1bSJed Brown v0ExT[d] += PetscRealPart(trans[d]); 316c4762a1bSJed Brown /* Only scale dimensions in the manifold */ 317c4762a1bSJed Brown for (e = 0; e < dim; ++e) { 318c4762a1bSJed Brown JExT[d * cdim + e] *= scale; 319c4762a1bSJed Brown invJExT[d * cdim + e] /= scale; 320c4762a1bSJed Brown } 321c4762a1bSJed Brown if (d < dim) detJExT *= scale; 322c4762a1bSJed Brown } 323c4762a1bSJed Brown /* Do scaling and translation before rotation, so that we can leave out the normal dimension for lower dimensional manifolds */ 324c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 325ad540459SPierre Jolivet for (e = 0, rot[d] = 0.0; e < cdim; ++e) rot[d] += R[d * cdim + e] * v0ExT[e]; 326c4762a1bSJed Brown } 327c4762a1bSJed Brown for (d = 0; d < cdim; ++d) v0ExT[d] = rot[d]; 328c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 329c4762a1bSJed Brown for (e = 0; e < cdim; ++e) { 330ad540459SPierre Jolivet for (f = 0, rotM[d * cdim + e] = 0.0; f < cdim; ++f) rotM[d * cdim + e] += R[d * cdim + f] * JExT[f * cdim + e]; 331c4762a1bSJed Brown } 332c4762a1bSJed Brown } 333c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 334ad540459SPierre Jolivet for (e = 0; e < cdim; ++e) JExT[d * cdim + e] = rotM[d * cdim + e]; 335c4762a1bSJed Brown } 336c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 337c4762a1bSJed Brown for (e = 0; e < cdim; ++e) { 338ad540459SPierre Jolivet for (f = 0, rotM[d * cdim + e] = 0.0; f < cdim; ++f) rotM[d * cdim + e] += invJExT[d * cdim + f] * R[e * cdim + f]; 339c4762a1bSJed Brown } 340c4762a1bSJed Brown } 3419371c9d4SSatish Balay for (d = 0; d < cdim; ++d) { 342ad540459SPierre Jolivet for (e = 0; e < cdim; ++e) invJExT[d * cdim + e] = rotM[d * cdim + e]; 3439371c9d4SSatish Balay } 344c4762a1bSJed Brown } 345c4762a1bSJed Brown if (centroidEx) { 346c4762a1bSJed Brown volExT = volEx; 347c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 348c4762a1bSJed Brown centroidExT[d] = centroidEx[d]; 349c4762a1bSJed Brown normalExT[d] = normalEx[d]; 350c4762a1bSJed Brown } 351c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 352c4762a1bSJed Brown centroidExT[d] *= scale; 353c4762a1bSJed Brown centroidExT[d] += PetscRealPart(trans[d]); 354c4762a1bSJed Brown normalExT[d] /= scale; 355c4762a1bSJed Brown /* Only scale dimensions in the manifold */ 356c4762a1bSJed Brown if (d < dim) volExT *= scale; 357c4762a1bSJed Brown } 358c4762a1bSJed Brown /* Do scaling and translation before rotation, so that we can leave out the normal dimension for lower dimensional manifolds */ 359c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 360ad540459SPierre Jolivet for (e = 0, rot[d] = 0.0; e < cdim; ++e) rot[d] += R[d * cdim + e] * centroidExT[e]; 361c4762a1bSJed Brown } 362c4762a1bSJed Brown for (d = 0; d < cdim; ++d) centroidExT[d] = rot[d]; 363c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 364ad540459SPierre Jolivet for (e = 0, rot[d] = 0.0; e < cdim; ++e) rot[d] += R[d * cdim + e] * normalExT[e]; 365c4762a1bSJed Brown } 366c4762a1bSJed Brown for (d = 0; d < cdim; ++d) normalExT[d] = rot[d]; 367c4762a1bSJed Brown for (d = 0, norm = 0.0; d < cdim; ++d) norm += PetscSqr(normalExT[d]); 368c4762a1bSJed Brown norm = PetscSqrtReal(norm); 3699371c9d4SSatish Balay if (norm != 0.) 3709371c9d4SSatish Balay for (d = 0; d < cdim; ++d) normalExT[d] /= norm; 371c4762a1bSJed Brown } 372c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 373c4762a1bSJed Brown for (p = 0; p < numCorners; ++p) { 374c4762a1bSJed Brown newCoords[p * cdim + d] *= scale; 375c4762a1bSJed Brown newCoords[p * cdim + d] += trans[d]; 376c4762a1bSJed Brown } 377c4762a1bSJed Brown } 378c4762a1bSJed Brown for (p = 0; p < numCorners; ++p) { 379c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 380ad540459SPierre Jolivet for (e = 0, rot[d] = 0.0; e < cdim; ++e) rot[d] += R[d * cdim + e] * PetscRealPart(newCoords[p * cdim + e]); 381c4762a1bSJed Brown } 382c4762a1bSJed Brown for (d = 0; d < cdim; ++d) newCoords[p * cdim + d] = rot[d]; 383c4762a1bSJed Brown } 384c4762a1bSJed Brown 3859566063dSJacob Faibussowitsch PetscCall(ChangeCoordinates(dm, cdim, newCoords)); 3861baa6e33SBarry Smith if (v0Ex) PetscCall(CheckFEMGeometry(dm, 0, cdim, v0ExT, JExT, invJExT, detJExT)); 387c4762a1bSJed Brown if (dim == depth && centroidEx) { 3889566063dSJacob Faibussowitsch PetscCall(CheckFVMGeometry(dm, cell, cdim, centroidExT, normalExT, volExT)); 3899566063dSJacob Faibussowitsch PetscCall(CheckGaussLaw(dm, cell)); 390c4762a1bSJed Brown if (faceCentroidEx) { 3919566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, cell, &coneSize)); 3929566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, cell, &cone)); 393c4762a1bSJed Brown for (c = 0; c < coneSize; ++c) { 394c4762a1bSJed Brown PetscInt off = c * cdim; 395c4762a1bSJed Brown 396c4762a1bSJed Brown faceVolExT = faceVolEx[c]; 397c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 398c4762a1bSJed Brown faceCentroidExT[d] = faceCentroidEx[off + d]; 399c4762a1bSJed Brown faceNormalExT[d] = faceNormalEx[off + d]; 400c4762a1bSJed Brown } 401c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 402c4762a1bSJed Brown faceCentroidExT[d] *= scale; 403c4762a1bSJed Brown faceCentroidExT[d] += PetscRealPart(trans[d]); 404c4762a1bSJed Brown faceNormalExT[d] /= scale; 405c4762a1bSJed Brown /* Only scale dimensions in the manifold */ 406ad540459SPierre Jolivet if (d < dim - 1) faceVolExT *= scale; 407c4762a1bSJed Brown } 408c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 409ad540459SPierre Jolivet for (e = 0, rot[d] = 0.0; e < cdim; ++e) rot[d] += R[d * cdim + e] * faceCentroidExT[e]; 410c4762a1bSJed Brown } 411c4762a1bSJed Brown for (d = 0; d < cdim; ++d) faceCentroidExT[d] = rot[d]; 412c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 413ad540459SPierre Jolivet for (e = 0, rot[d] = 0.0; e < cdim; ++e) rot[d] += R[d * cdim + e] * faceNormalExT[e]; 414c4762a1bSJed Brown } 415c4762a1bSJed Brown for (d = 0; d < cdim; ++d) faceNormalExT[d] = rot[d]; 416c4762a1bSJed Brown for (d = 0, norm = 0.0; d < cdim; ++d) norm += PetscSqr(faceNormalExT[d]); 417c4762a1bSJed Brown norm = PetscSqrtReal(norm); 418c4762a1bSJed Brown for (d = 0; d < cdim; ++d) faceNormalExT[d] /= norm; 419c4762a1bSJed Brown 4209566063dSJacob Faibussowitsch PetscCall(CheckFVMGeometry(dm, cone[c], cdim, faceCentroidExT, faceNormalExT, faceVolExT)); 421c4762a1bSJed Brown } 422c4762a1bSJed Brown } 423c4762a1bSJed Brown } 424c4762a1bSJed Brown } 4259566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&r)); 4269566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&ang)); 4279566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&ang2)); 4289566063dSJacob Faibussowitsch PetscCall(PetscFree2(origCoords, newCoords)); 4299566063dSJacob Faibussowitsch PetscCall(PetscFree5(v0ExT, JExT, invJExT, centroidExT, normalExT)); 4309566063dSJacob Faibussowitsch PetscCall(PetscFree2(faceCentroidExT, faceNormalExT)); 431c4762a1bSJed Brown } 4323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 433c4762a1bSJed Brown } 434c4762a1bSJed Brown 435d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestTriangle(MPI_Comm comm, PetscBool transform) 436d71ae5a4SJacob Faibussowitsch { 437c4762a1bSJed Brown DM dm; 438c4762a1bSJed Brown 439c4762a1bSJed Brown PetscFunctionBegin; 4409566063dSJacob Faibussowitsch PetscCall(DMPlexCreateReferenceCell(comm, DM_POLYTOPE_TRIANGLE, &dm)); 4419566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 442c4762a1bSJed Brown /* Check reference geometry: determinant is scaled by reference volume (2.0) */ 443c4762a1bSJed Brown { 444c4762a1bSJed Brown PetscReal v0Ex[2] = {-1.0, -1.0}; 445c4762a1bSJed Brown PetscReal JEx[4] = {1.0, 0.0, 0.0, 1.0}; 446c4762a1bSJed Brown PetscReal invJEx[4] = {1.0, 0.0, 0.0, 1.0}; 447c4762a1bSJed Brown PetscReal detJEx = 1.0; 448c4762a1bSJed Brown PetscReal centroidEx[2] = {-((PetscReal)1.) / ((PetscReal)3.), -((PetscReal)1.) / ((PetscReal)3.)}; 449c4762a1bSJed Brown PetscReal normalEx[2] = {0.0, 0.0}; 450c4762a1bSJed Brown PetscReal volEx = 2.0; 451c4762a1bSJed Brown 4529566063dSJacob Faibussowitsch PetscCall(CheckCell(dm, 0, transform, v0Ex, JEx, invJEx, detJEx, centroidEx, normalEx, volEx, NULL, NULL, NULL)); 453c4762a1bSJed Brown } 454c4762a1bSJed Brown /* Move to 3D: Check reference geometry: determinant is scaled by reference volume (2.0) */ 455c4762a1bSJed Brown { 456c4762a1bSJed Brown PetscScalar vertexCoords[9] = {-1.0, -1.0, 0.0, 1.0, -1.0, 0.0, -1.0, 1.0, 0.0}; 457c4762a1bSJed Brown PetscReal v0Ex[3] = {-1.0, -1.0, 0.0}; 458c4762a1bSJed Brown PetscReal JEx[9] = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0}; 459c4762a1bSJed Brown PetscReal invJEx[9] = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0}; 460c4762a1bSJed Brown PetscReal detJEx = 1.0; 461c4762a1bSJed Brown PetscReal centroidEx[3] = {-((PetscReal)1.) / ((PetscReal)3.), -((PetscReal)1.) / ((PetscReal)3.), 0.0}; 462c4762a1bSJed Brown PetscReal normalEx[3] = {0.0, 0.0, 1.0}; 463c4762a1bSJed Brown PetscReal volEx = 2.0; 464c4762a1bSJed Brown 4659566063dSJacob Faibussowitsch PetscCall(ChangeCoordinates(dm, 3, vertexCoords)); 4669566063dSJacob Faibussowitsch PetscCall(CheckCell(dm, 0, transform, v0Ex, JEx, invJEx, detJEx, centroidEx, normalEx, volEx, NULL, NULL, NULL)); 467c4762a1bSJed Brown } 468c4762a1bSJed Brown /* Cleanup */ 4699566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 4703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 471c4762a1bSJed Brown } 472c4762a1bSJed Brown 473d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestQuadrilateral(MPI_Comm comm, PetscBool transform) 474d71ae5a4SJacob Faibussowitsch { 475c4762a1bSJed Brown DM dm; 476c4762a1bSJed Brown 477c4762a1bSJed Brown PetscFunctionBegin; 4789566063dSJacob Faibussowitsch PetscCall(DMPlexCreateReferenceCell(comm, DM_POLYTOPE_QUADRILATERAL, &dm)); 4799566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 480c4762a1bSJed Brown /* Check reference geometry: determinant is scaled by reference volume (2.0) */ 481c4762a1bSJed Brown { 482c4762a1bSJed Brown PetscReal v0Ex[2] = {-1.0, -1.0}; 483c4762a1bSJed Brown PetscReal JEx[4] = {1.0, 0.0, 0.0, 1.0}; 484c4762a1bSJed Brown PetscReal invJEx[4] = {1.0, 0.0, 0.0, 1.0}; 485c4762a1bSJed Brown PetscReal detJEx = 1.0; 486c4762a1bSJed Brown PetscReal centroidEx[2] = {0.0, 0.0}; 487c4762a1bSJed Brown PetscReal normalEx[2] = {0.0, 0.0}; 488c4762a1bSJed Brown PetscReal volEx = 4.0; 489c4762a1bSJed Brown 4909566063dSJacob Faibussowitsch PetscCall(CheckCell(dm, 0, transform, v0Ex, JEx, invJEx, detJEx, centroidEx, normalEx, volEx, NULL, NULL, NULL)); 491c4762a1bSJed Brown } 492c4762a1bSJed Brown /* Move to 3D: Check reference geometry: determinant is scaled by reference volume (4.0) */ 493c4762a1bSJed Brown { 494c4762a1bSJed 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}; 495c4762a1bSJed Brown PetscReal v0Ex[3] = {-1.0, -1.0, 0.0}; 496c4762a1bSJed Brown PetscReal JEx[9] = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0}; 497c4762a1bSJed Brown PetscReal invJEx[9] = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0}; 498c4762a1bSJed Brown PetscReal detJEx = 1.0; 499c4762a1bSJed Brown PetscReal centroidEx[3] = {0.0, 0.0, 0.0}; 500c4762a1bSJed Brown PetscReal normalEx[3] = {0.0, 0.0, 1.0}; 501c4762a1bSJed Brown PetscReal volEx = 4.0; 502c4762a1bSJed Brown 5039566063dSJacob Faibussowitsch PetscCall(ChangeCoordinates(dm, 3, vertexCoords)); 5049566063dSJacob Faibussowitsch PetscCall(CheckCell(dm, 0, transform, v0Ex, JEx, invJEx, detJEx, centroidEx, normalEx, volEx, NULL, NULL, NULL)); 505c4762a1bSJed Brown } 506c4762a1bSJed Brown /* Cleanup */ 5079566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 5083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 509c4762a1bSJed Brown } 510c4762a1bSJed Brown 511d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestTetrahedron(MPI_Comm comm, PetscBool transform) 512d71ae5a4SJacob Faibussowitsch { 513c4762a1bSJed Brown DM dm; 514c4762a1bSJed Brown 515c4762a1bSJed Brown PetscFunctionBegin; 5169566063dSJacob Faibussowitsch PetscCall(DMPlexCreateReferenceCell(comm, DM_POLYTOPE_TETRAHEDRON, &dm)); 5179566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 518c4762a1bSJed Brown /* Check reference geometry: determinant is scaled by reference volume (4/3) */ 519c4762a1bSJed Brown { 520c4762a1bSJed Brown PetscReal v0Ex[3] = {-1.0, -1.0, -1.0}; 521c4762a1bSJed Brown PetscReal JEx[9] = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0}; 522c4762a1bSJed Brown PetscReal invJEx[9] = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0}; 523c4762a1bSJed Brown PetscReal detJEx = 1.0; 524c4762a1bSJed Brown PetscReal centroidEx[3] = {-0.5, -0.5, -0.5}; 525c4762a1bSJed Brown PetscReal normalEx[3] = {0.0, 0.0, 0.0}; 526c4762a1bSJed Brown PetscReal volEx = (PetscReal)4.0 / (PetscReal)3.0; 527c4762a1bSJed Brown 5289566063dSJacob Faibussowitsch PetscCall(CheckCell(dm, 0, transform, v0Ex, JEx, invJEx, detJEx, centroidEx, normalEx, volEx, NULL, NULL, NULL)); 529c4762a1bSJed Brown } 530c4762a1bSJed Brown /* Cleanup */ 5319566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 5323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 533c4762a1bSJed Brown } 534c4762a1bSJed Brown 535d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestHexahedron(MPI_Comm comm, PetscBool transform) 536d71ae5a4SJacob Faibussowitsch { 537c4762a1bSJed Brown DM dm; 538c4762a1bSJed Brown 539c4762a1bSJed Brown PetscFunctionBegin; 5409566063dSJacob Faibussowitsch PetscCall(DMPlexCreateReferenceCell(comm, DM_POLYTOPE_HEXAHEDRON, &dm)); 5419566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 542c4762a1bSJed Brown /* Check reference geometry: determinant is scaled by reference volume 8.0 */ 543c4762a1bSJed Brown { 544c4762a1bSJed Brown PetscReal v0Ex[3] = {-1.0, -1.0, -1.0}; 545c4762a1bSJed Brown PetscReal JEx[9] = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0}; 546c4762a1bSJed Brown PetscReal invJEx[9] = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0}; 547c4762a1bSJed Brown PetscReal detJEx = 1.0; 548c4762a1bSJed Brown PetscReal centroidEx[3] = {0.0, 0.0, 0.0}; 549c4762a1bSJed Brown PetscReal normalEx[3] = {0.0, 0.0, 0.0}; 550c4762a1bSJed Brown PetscReal volEx = 8.0; 551c4762a1bSJed Brown 5529566063dSJacob Faibussowitsch PetscCall(CheckCell(dm, 0, transform, v0Ex, JEx, invJEx, detJEx, centroidEx, normalEx, volEx, NULL, NULL, NULL)); 553c4762a1bSJed Brown } 554c4762a1bSJed Brown /* Cleanup */ 5559566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 5563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 557c4762a1bSJed Brown } 558c4762a1bSJed Brown 559d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestHexahedronCurved(MPI_Comm comm) 560d71ae5a4SJacob Faibussowitsch { 561c4762a1bSJed Brown DM dm; 5629371c9d4SSatish Balay 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, -1.0, -1.0, 1.1, 1.0, -1.0, 1.0, 1.0, 1.0, 1.1, -1.0, 1.0, 1.0}; 563c4762a1bSJed Brown 564c4762a1bSJed Brown PetscFunctionBegin; 5659566063dSJacob Faibussowitsch PetscCall(DMPlexCreateReferenceCell(comm, DM_POLYTOPE_HEXAHEDRON, &dm)); 5669566063dSJacob Faibussowitsch PetscCall(ChangeCoordinates(dm, 3, coords)); 5679566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 5684f99dae5SMatthew G. Knepley { 5694f99dae5SMatthew G. Knepley PetscReal centroidEx[3] = {0.0, 0.0, 0.016803278688524603}; 5704f99dae5SMatthew G. Knepley PetscReal normalEx[3] = {0.0, 0.0, 0.0}; 5714f99dae5SMatthew G. Knepley PetscReal volEx = 8.1333333333333346; 5724f99dae5SMatthew G. Knepley 5739566063dSJacob Faibussowitsch PetscCall(CheckCell(dm, 0, PETSC_FALSE, NULL, NULL, NULL, 0.0, centroidEx, normalEx, volEx, NULL, NULL, NULL)); 574c4762a1bSJed Brown } 5759566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 5763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5774f99dae5SMatthew G. Knepley } 5784f99dae5SMatthew G. Knepley 5794f99dae5SMatthew G. Knepley /* This wedge is a tensor product cell, rather than a normal wedge */ 580d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestWedge(MPI_Comm comm, PetscBool transform) 581d71ae5a4SJacob Faibussowitsch { 5824f99dae5SMatthew G. Knepley DM dm; 5834f99dae5SMatthew G. Knepley 5844f99dae5SMatthew G. Knepley PetscFunctionBegin; 5859566063dSJacob Faibussowitsch PetscCall(DMPlexCreateReferenceCell(comm, DM_POLYTOPE_TRI_PRISM_TENSOR, &dm)); 5869566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 587c4762a1bSJed Brown /* Check reference geometry: determinant is scaled by reference volume 4.0 */ 588c4762a1bSJed Brown { 589c4762a1bSJed Brown #if 0 590c4762a1bSJed Brown /* FEM geometry is not functional for wedges */ 591c4762a1bSJed Brown PetscReal v0Ex[3] = {-1.0, -1.0, -1.0}; 592c4762a1bSJed Brown PetscReal JEx[9] = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0}; 593c4762a1bSJed Brown PetscReal invJEx[9] = {1.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0}; 594c4762a1bSJed Brown PetscReal detJEx = 1.0; 595c4762a1bSJed Brown #endif 596c4762a1bSJed Brown 5974f99dae5SMatthew G. Knepley { 598c4762a1bSJed Brown PetscReal centroidEx[3] = {-((PetscReal)1.) / ((PetscReal)3.), -((PetscReal)1.) / ((PetscReal)3.), 0.0}; 599c4762a1bSJed Brown PetscReal normalEx[3] = {0.0, 0.0, 0.0}; 600c4762a1bSJed Brown PetscReal volEx = 4.0; 601c4762a1bSJed Brown PetscReal faceVolEx[5] = {2.0, 2.0, 4.0, PETSC_SQRT2 * 4.0, 4.0}; 602c4762a1bSJed 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}; 6039371c9d4SSatish Balay PetscReal faceCentroidEx[15] = {-((PetscReal)1.) / ((PetscReal)3.), -((PetscReal)1.) / ((PetscReal)3.), -1.0, -((PetscReal)1.) / ((PetscReal)3.), -((PetscReal)1.) / ((PetscReal)3.), 1.0, 0.0, -1.0, 0.0, 0.0, 0.0, 0.0, -1.0, 0.0, 0.0}; 604c4762a1bSJed Brown 6059566063dSJacob Faibussowitsch PetscCall(CheckCell(dm, 0, transform, NULL, NULL, NULL, 0.0, centroidEx, normalEx, volEx, faceCentroidEx, faceNormalEx, faceVolEx)); 606c4762a1bSJed Brown } 607c4762a1bSJed Brown } 608c4762a1bSJed Brown /* Cleanup */ 6099566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 6103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 611c4762a1bSJed Brown } 612c4762a1bSJed Brown 613d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 614d71ae5a4SJacob Faibussowitsch { 615c4762a1bSJed Brown AppCtx user; 616c4762a1bSJed Brown 617327415f7SBarry Smith PetscFunctionBeginUser; 6189566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 6199566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 620c4762a1bSJed Brown if (user.runType == RUN_REFERENCE) { 6219566063dSJacob Faibussowitsch PetscCall(TestTriangle(PETSC_COMM_SELF, user.transform)); 6229566063dSJacob Faibussowitsch PetscCall(TestQuadrilateral(PETSC_COMM_SELF, user.transform)); 6239566063dSJacob Faibussowitsch PetscCall(TestTetrahedron(PETSC_COMM_SELF, user.transform)); 6249566063dSJacob Faibussowitsch PetscCall(TestHexahedron(PETSC_COMM_SELF, user.transform)); 6259566063dSJacob Faibussowitsch PetscCall(TestWedge(PETSC_COMM_SELF, user.transform)); 6264f99dae5SMatthew G. Knepley } else if (user.runType == RUN_HEX_CURVED) { 6279566063dSJacob Faibussowitsch PetscCall(TestHexahedronCurved(PETSC_COMM_SELF)); 628c4762a1bSJed Brown } else if (user.runType == RUN_FILE) { 629c4762a1bSJed Brown PetscInt dim, cStart, cEnd, c; 630c4762a1bSJed Brown 6319566063dSJacob Faibussowitsch PetscCall(DMGetDimension(user.dm, &dim)); 6329566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(user.dm, 0, &cStart, &cEnd)); 633c4762a1bSJed Brown for (c = 0; c < cEnd - cStart; ++c) { 6344f99dae5SMatthew G. Knepley PetscReal *v0 = user.v0 ? &user.v0[c * dim] : NULL; 6354f99dae5SMatthew G. Knepley PetscReal *J = user.J ? &user.J[c * dim * dim] : NULL; 6364f99dae5SMatthew G. Knepley PetscReal *invJ = user.invJ ? &user.invJ[c * dim * dim] : NULL; 6374f99dae5SMatthew G. Knepley PetscReal detJ = user.detJ ? user.detJ[c] : 0.0; 6384f99dae5SMatthew G. Knepley PetscReal *centroid = user.centroid ? &user.centroid[c * dim] : NULL; 6394f99dae5SMatthew G. Knepley PetscReal *normal = user.normal ? &user.normal[c * dim] : NULL; 6404f99dae5SMatthew G. Knepley PetscReal vol = user.vol ? user.vol[c] : 0.0; 6414f99dae5SMatthew G. Knepley 6429566063dSJacob Faibussowitsch PetscCall(CheckCell(user.dm, c + cStart, PETSC_FALSE, v0, J, invJ, detJ, centroid, normal, vol, NULL, NULL, NULL)); 643c4762a1bSJed Brown } 6449566063dSJacob Faibussowitsch PetscCall(PetscFree4(user.v0, user.J, user.invJ, user.detJ)); 6459566063dSJacob Faibussowitsch PetscCall(PetscFree(user.centroid)); 6469566063dSJacob Faibussowitsch PetscCall(PetscFree(user.normal)); 6479566063dSJacob Faibussowitsch PetscCall(PetscFree(user.vol)); 6489566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.dm)); 649c4762a1bSJed Brown } else if (user.runType == RUN_DISPLAY) { 650c4762a1bSJed Brown DM gdm, dmCell; 651c4762a1bSJed Brown Vec cellgeom, facegeom; 652c4762a1bSJed Brown const PetscScalar *cgeom; 653c4762a1bSJed Brown PetscInt dim, d, cStart, cEnd, cEndInterior, c; 654c4762a1bSJed Brown 6559566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(user.dm, &dim)); 6569566063dSJacob Faibussowitsch PetscCall(DMPlexConstructGhostCells(user.dm, NULL, NULL, &gdm)); 657c4762a1bSJed Brown if (gdm) { 6589566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.dm)); 659c4762a1bSJed Brown user.dm = gdm; 660c4762a1bSJed Brown } 6619566063dSJacob Faibussowitsch PetscCall(DMPlexComputeGeometryFVM(user.dm, &cellgeom, &facegeom)); 6629566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(user.dm, 0, &cStart, &cEnd)); 663*2827ebadSStefano Zampini PetscCall(DMPlexGetCellTypeStratum(user.dm, DM_POLYTOPE_FV_GHOST, &cEndInterior, NULL)); 664c4762a1bSJed Brown if (cEndInterior >= 0) cEnd = cEndInterior; 6659566063dSJacob Faibussowitsch PetscCall(VecGetDM(cellgeom, &dmCell)); 6669566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(cellgeom, &cgeom)); 667c4762a1bSJed Brown for (c = 0; c < cEnd - cStart; ++c) { 668c4762a1bSJed Brown PetscFVCellGeom *cg; 669c4762a1bSJed Brown 6709566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRead(dmCell, c, cgeom, &cg)); 67163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Cell %4" PetscInt_FMT ": Centroid (", c)); 672c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 6739566063dSJacob Faibussowitsch if (d > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", ")); 67463a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "%12.2g", (double)cg->centroid[d])); 675c4762a1bSJed Brown } 67663a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, ") Vol %12.2g\n", (double)cg->volume)); 677c4762a1bSJed Brown } 6789566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(cellgeom, &cgeom)); 6799566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cellgeom)); 6809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&facegeom)); 6819566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.dm)); 682c4762a1bSJed Brown } 6839566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 684b122ec5aSJacob Faibussowitsch return 0; 685c4762a1bSJed Brown } 686c4762a1bSJed Brown 687c4762a1bSJed Brown /*TEST 688c4762a1bSJed Brown 689c4762a1bSJed Brown test: 6904f99dae5SMatthew G. Knepley suffix: 1 691c4762a1bSJed Brown args: -dm_view ascii::ascii_info_detail 692c4762a1bSJed Brown test: 693c4762a1bSJed Brown suffix: 2 6944f99dae5SMatthew G. Knepley args: -run_type hex_curved 695c4762a1bSJed Brown test: 696c4762a1bSJed Brown suffix: 3 6974f99dae5SMatthew G. Knepley args: -transform 698c4762a1bSJed Brown test: 699c4762a1bSJed Brown suffix: 4 700c4762a1bSJed Brown requires: exodusii 7014f99dae5SMatthew 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 702c4762a1bSJed Brown test: 703c4762a1bSJed Brown suffix: 5 7044f99dae5SMatthew 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 7059bf2564aSMatt McGurn test: 7069bf2564aSMatt McGurn suffix: 6 7079bf2564aSMatt 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 708c4762a1bSJed Brown TEST*/ 709