18b438e21SMatthew G. Knepley static char help[] = "Tests for high order geometry\n\n"; 28b438e21SMatthew G. Knepley 38b438e21SMatthew G. Knepley #include <petscdmplex.h> 48b438e21SMatthew G. Knepley #include <petscds.h> 58b438e21SMatthew G. Knepley 6*9371c9d4SSatish Balay typedef enum { 7*9371c9d4SSatish Balay TRANSFORM_NONE, 8*9371c9d4SSatish Balay TRANSFORM_SHEAR, 9*9371c9d4SSatish Balay TRANSFORM_FLARE, 10*9371c9d4SSatish Balay TRANSFORM_ANNULUS, 11*9371c9d4SSatish Balay TRANSFORM_SHELL 12*9371c9d4SSatish Balay } Transform; 1360c1a66aSMatthew G. Knepley const char *const TransformTypes[] = {"none", "shear", "flare", "annulus", "shell", "Mesh Transform", "TRANSFORM_", NULL}; 148b438e21SMatthew G. Knepley 158b438e21SMatthew G. Knepley typedef struct { 1630602db0SMatthew G. Knepley PetscBool coordSpace; /* Flag to create coordinate space */ 178b438e21SMatthew G. Knepley Transform meshTransform; /* Transform for initial box mesh */ 188b438e21SMatthew G. Knepley PetscReal *transformDataReal; /* Parameters for mesh transform */ 198b438e21SMatthew G. Knepley PetscScalar *transformData; /* Parameters for mesh transform */ 208b438e21SMatthew G. Knepley PetscReal volume; /* Analytical volume of the mesh */ 218b438e21SMatthew G. Knepley PetscReal tol; /* Tolerance for volume check */ 228b438e21SMatthew G. Knepley } AppCtx; 238b438e21SMatthew G. Knepley 24*9371c9d4SSatish Balay PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) { 258b438e21SMatthew G. Knepley PetscInt n = 0, i; 268b438e21SMatthew G. Knepley 278b438e21SMatthew G. Knepley PetscFunctionBegin; 2830602db0SMatthew G. Knepley options->coordSpace = PETSC_TRUE; 298b438e21SMatthew G. Knepley options->meshTransform = TRANSFORM_NONE; 308b438e21SMatthew G. Knepley options->transformDataReal = NULL; 318b438e21SMatthew G. Knepley options->transformData = NULL; 328b438e21SMatthew G. Knepley options->volume = -1.0; 338b438e21SMatthew G. Knepley options->tol = PETSC_SMALL; 348b438e21SMatthew G. Knepley 35d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Meshing Interpolation Test Options", "DMPLEX"); 369566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-coord_space", "Flag to create a coordinate space", "ex33.c", options->coordSpace, &options->coordSpace, NULL)); 379566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-mesh_transform", "Method to transform initial box mesh <none, shear, annulus, shell>", "ex33.c", TransformTypes, (PetscEnum)options->meshTransform, (PetscEnum *)&options->meshTransform, NULL)); 388b438e21SMatthew G. Knepley switch (options->meshTransform) { 398b438e21SMatthew G. Knepley case TRANSFORM_NONE: break; 408b438e21SMatthew G. Knepley case TRANSFORM_SHEAR: 4139290ff6SMatthew G. Knepley n = 2; 429566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &options->transformDataReal)); 438b438e21SMatthew G. Knepley for (i = 0; i < n; ++i) options->transformDataReal[i] = 1.0; 449566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-transform_data", "Parameters for mesh transforms", "ex33.c", options->transformDataReal, &n, NULL)); 458b438e21SMatthew G. Knepley break; 4660c1a66aSMatthew G. Knepley case TRANSFORM_FLARE: 4760c1a66aSMatthew G. Knepley n = 4; 489566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &options->transformData)); 4960c1a66aSMatthew G. Knepley for (i = 0; i < n; ++i) options->transformData[i] = 1.0; 5060c1a66aSMatthew G. Knepley options->transformData[0] = (PetscScalar)0; 519566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalarArray("-transform_data", "Parameters for mesh transforms", "ex33.c", options->transformData, &n, NULL)); 5260c1a66aSMatthew G. Knepley break; 538b438e21SMatthew G. Knepley case TRANSFORM_ANNULUS: 548b438e21SMatthew G. Knepley n = 2; 559566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &options->transformData)); 568b438e21SMatthew G. Knepley options->transformData[0] = 1.0; 578b438e21SMatthew G. Knepley options->transformData[1] = 2.0; 589566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalarArray("-transform_data", "Parameters for mesh transforms", "ex33.c", options->transformData, &n, NULL)); 598b438e21SMatthew G. Knepley break; 608b438e21SMatthew G. Knepley case TRANSFORM_SHELL: 618b438e21SMatthew G. Knepley n = 2; 629566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &options->transformData)); 638b438e21SMatthew G. Knepley options->transformData[0] = 1.0; 648b438e21SMatthew G. Knepley options->transformData[1] = 2.0; 659566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalarArray("-transform_data", "Parameters for mesh transforms", "ex33.c", options->transformData, &n, NULL)); 668b438e21SMatthew G. Knepley break; 6763a3b9bcSJacob Faibussowitsch default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Unknown mesh transform %d", options->meshTransform); 688b438e21SMatthew G. Knepley } 699566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-volume", "The analytical volume of the mesh", "ex33.c", options->volume, &options->volume, NULL)); 709566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tol", "The tolerance for the volume check", "ex33.c", options->tol, &options->tol, NULL)); 71d0609cedSBarry Smith PetscOptionsEnd(); 728b438e21SMatthew G. Knepley PetscFunctionReturn(0); 738b438e21SMatthew G. Knepley } 748b438e21SMatthew G. Knepley 75*9371c9d4SSatish Balay static void identity(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) { 768b438e21SMatthew G. Knepley const PetscInt Nc = uOff[1] - uOff[0]; 778b438e21SMatthew G. Knepley PetscInt c; 788b438e21SMatthew G. Knepley 798b438e21SMatthew G. Knepley for (c = 0; c < Nc; ++c) f0[c] = u[c]; 808b438e21SMatthew G. Knepley } 818b438e21SMatthew G. Knepley 8260c1a66aSMatthew G. Knepley /* Flare applies the transformation, assuming we fix x_f, 8360c1a66aSMatthew G. Knepley 8460c1a66aSMatthew G. Knepley x_i = x_i * alpha_i x_f 8560c1a66aSMatthew G. Knepley */ 86*9371c9d4SSatish Balay static void f0_flare(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar coords[]) { 8760c1a66aSMatthew G. Knepley const PetscInt Nc = uOff[1] - uOff[0]; 8860c1a66aSMatthew G. Knepley const PetscInt cf = (PetscInt)PetscRealPart(constants[0]); 8960c1a66aSMatthew G. Knepley PetscInt c; 9060c1a66aSMatthew G. Knepley 91*9371c9d4SSatish Balay for (c = 0; c < Nc; ++c) { coords[c] = u[c] * (c == cf ? 1.0 : constants[c + 1] * u[cf]); } 9260c1a66aSMatthew G. Knepley } 9360c1a66aSMatthew G. Knepley 948b438e21SMatthew G. Knepley /* 958b438e21SMatthew G. Knepley We would like to map the unit square to a quarter of the annulus between circles of radius 1 and 2. We start by mapping the straight sections, which 968b438e21SMatthew G. Knepley will correspond to the top and bottom of our square. So 978b438e21SMatthew G. Knepley 988b438e21SMatthew G. Knepley (0,0)--(1,0) ==> (1,0)--(2,0) Just a shift of (1,0) 998b438e21SMatthew G. Knepley (0,1)--(1,1) ==> (0,1)--(0,2) Switch x and y 1008b438e21SMatthew G. Knepley 1018b438e21SMatthew G. Knepley So it looks like we want to map each layer in y to a ray, so x is the radius and y is the angle: 1028b438e21SMatthew G. Knepley 1038b438e21SMatthew G. Knepley (x, y) ==> (x+1, \pi/2 y) in (r', \theta') space 1048b438e21SMatthew G. Knepley ==> ((x+1) cos(\pi/2 y), (x+1) sin(\pi/2 y)) in (x', y') space 1058b438e21SMatthew G. Knepley */ 106*9371c9d4SSatish Balay static void f0_annulus(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar xp[]) { 1078b438e21SMatthew G. Knepley const PetscReal ri = PetscRealPart(constants[0]); 1088b438e21SMatthew G. Knepley const PetscReal ro = PetscRealPart(constants[1]); 1098b438e21SMatthew G. Knepley 1108b438e21SMatthew G. Knepley xp[0] = (x[0] * (ro - ri) + ri) * PetscCosReal(0.5 * PETSC_PI * x[1]); 1118b438e21SMatthew G. Knepley xp[1] = (x[0] * (ro - ri) + ri) * PetscSinReal(0.5 * PETSC_PI * x[1]); 1128b438e21SMatthew G. Knepley } 1138b438e21SMatthew G. Knepley 1148b438e21SMatthew G. Knepley /* 1158b438e21SMatthew G. Knepley We would like to map the unit cube to a hemisphere of the spherical shell between balls of radius 1 and 2. We want to map the bottom surface onto the 1168b438e21SMatthew G. Knepley lower hemisphere and the upper surface onto the top, letting z be the radius. 1178b438e21SMatthew G. Knepley 1188b438e21SMatthew G. Knepley (x, y) ==> ((z+3)/2, \pi/2 (|x| or |y|), arctan y/x) in (r', \theta', \phi') space 1198b438e21SMatthew G. Knepley ==> ((z+3)/2 \cos(\theta') cos(\phi'), (z+3)/2 \cos(\theta') sin(\phi'), (z+3)/2 sin(\theta')) in (x', y', z') space 1208b438e21SMatthew G. Knepley */ 121*9371c9d4SSatish Balay static void f0_shell(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar xp[]) { 1228b438e21SMatthew G. Knepley const PetscReal pi4 = PETSC_PI / 4.0; 1238b438e21SMatthew G. Knepley const PetscReal ri = PetscRealPart(constants[0]); 1248b438e21SMatthew G. Knepley const PetscReal ro = PetscRealPart(constants[1]); 1258b438e21SMatthew G. Knepley const PetscReal rp = (x[2] + 1) * 0.5 * (ro - ri) + ri; 1268b438e21SMatthew G. Knepley const PetscReal phip = PetscAtan2Real(x[1], x[0]); 1278b438e21SMatthew G. Knepley const PetscReal thetap = 0.5 * PETSC_PI * (1.0 - ((((phip <= pi4) && (phip >= -pi4)) || ((phip >= 3.0 * pi4) || (phip <= -3.0 * pi4))) ? PetscAbsReal(x[0]) : PetscAbsReal(x[1]))); 1288b438e21SMatthew G. Knepley 1298b438e21SMatthew G. Knepley xp[0] = rp * PetscCosReal(thetap) * PetscCosReal(phip); 1308b438e21SMatthew G. Knepley xp[1] = rp * PetscCosReal(thetap) * PetscSinReal(phip); 1318b438e21SMatthew G. Knepley xp[2] = rp * PetscSinReal(thetap); 1328b438e21SMatthew G. Knepley } 1338b438e21SMatthew G. Knepley 134*9371c9d4SSatish Balay static PetscErrorCode DMCreateCoordinateDisc(DM dm) { 1358b438e21SMatthew G. Knepley DM cdm; 1368b438e21SMatthew G. Knepley PetscFE fe; 13739290ff6SMatthew G. Knepley DMPolytopeType ct; 13839290ff6SMatthew G. Knepley PetscInt dim, dE, cStart; 13939290ff6SMatthew G. Knepley PetscBool simplex; 1408b438e21SMatthew G. Knepley 1418b438e21SMatthew G. Knepley PetscFunctionBegin; 1429566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &cdm)); 1439566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1449566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &dE)); 1459566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(cdm, 0, &cStart, NULL)); 1469566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 14739290ff6SMatthew G. Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE; 1489566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, dE, simplex, "dm_coord_", -1, &fe)); 1499566063dSJacob Faibussowitsch PetscCall(DMProjectCoordinates(dm, fe)); 1509566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 1518b438e21SMatthew G. Knepley PetscFunctionReturn(0); 1528b438e21SMatthew G. Knepley } 1538b438e21SMatthew G. Knepley 154*9371c9d4SSatish Balay PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm) { 1558b438e21SMatthew G. Knepley DM cdm; 1568b438e21SMatthew G. Knepley PetscDS cds; 1578b438e21SMatthew G. Knepley 1588b438e21SMatthew G. Knepley PetscFunctionBegin; 1599566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 1609566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 1619566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 16230602db0SMatthew G. Knepley 1639566063dSJacob Faibussowitsch if (ctx->coordSpace) PetscCall(DMCreateCoordinateDisc(*dm)); 1648b438e21SMatthew G. Knepley switch (ctx->meshTransform) { 165*9371c9d4SSatish Balay case TRANSFORM_NONE: PetscCall(DMPlexRemapGeometry(*dm, 0.0, identity)); break; 166*9371c9d4SSatish Balay case TRANSFORM_SHEAR: PetscCall(DMPlexShearGeometry(*dm, DM_X, ctx->transformDataReal)); break; 16760c1a66aSMatthew G. Knepley case TRANSFORM_FLARE: 1689566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(*dm, &cdm)); 1699566063dSJacob Faibussowitsch PetscCall(DMGetDS(cdm, &cds)); 1709566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(cds, 4, ctx->transformData)); 1719566063dSJacob Faibussowitsch PetscCall(DMPlexRemapGeometry(*dm, 0.0, f0_flare)); 17260c1a66aSMatthew G. Knepley break; 1738b438e21SMatthew G. Knepley case TRANSFORM_ANNULUS: 1749566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(*dm, &cdm)); 1759566063dSJacob Faibussowitsch PetscCall(DMGetDS(cdm, &cds)); 1769566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(cds, 2, ctx->transformData)); 1779566063dSJacob Faibussowitsch PetscCall(DMPlexRemapGeometry(*dm, 0.0, f0_annulus)); 1788b438e21SMatthew G. Knepley break; 1798b438e21SMatthew G. Knepley case TRANSFORM_SHELL: 1809566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(*dm, &cdm)); 1819566063dSJacob Faibussowitsch PetscCall(DMGetDS(cdm, &cds)); 1829566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(cds, 2, ctx->transformData)); 1839566063dSJacob Faibussowitsch PetscCall(DMPlexRemapGeometry(*dm, 0.0, f0_shell)); 1848b438e21SMatthew G. Knepley break; 18563a3b9bcSJacob Faibussowitsch default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Unknown mesh transform %d", ctx->meshTransform); 1868b438e21SMatthew G. Knepley } 1879566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 1888b438e21SMatthew G. Knepley PetscFunctionReturn(0); 1898b438e21SMatthew G. Knepley } 1908b438e21SMatthew G. Knepley 191*9371c9d4SSatish Balay static void volume(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar vol[]) { 1928b438e21SMatthew G. Knepley vol[0] = 1.; 1938b438e21SMatthew G. Knepley } 1948b438e21SMatthew G. Knepley 195*9371c9d4SSatish Balay static PetscErrorCode CreateDiscretization(DM dm, AppCtx *ctx) { 1968b438e21SMatthew G. Knepley PetscDS ds; 1978b438e21SMatthew G. Knepley PetscFE fe; 19839290ff6SMatthew G. Knepley DMPolytopeType ct; 19939290ff6SMatthew G. Knepley PetscInt dim, cStart; 20039290ff6SMatthew G. Knepley PetscBool simplex; 2018b438e21SMatthew G. Knepley 2028b438e21SMatthew G. Knepley PetscFunctionBeginUser; 2039566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 2049566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 2059566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 20639290ff6SMatthew G. Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE; 2079566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, NULL, PETSC_DETERMINE, &fe)); 2089566063dSJacob Faibussowitsch PetscCall(PetscFESetName(fe, "scalar")); 2099566063dSJacob Faibussowitsch PetscCall(DMAddField(dm, NULL, (PetscObject)fe)); 2109566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 2119566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 2129566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 2139566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(ds, 0, volume)); 2148b438e21SMatthew G. Knepley PetscFunctionReturn(0); 2158b438e21SMatthew G. Knepley } 2168b438e21SMatthew G. Knepley 217*9371c9d4SSatish Balay static PetscErrorCode CheckVolume(DM dm, AppCtx *ctx) { 2188b438e21SMatthew G. Knepley Vec u; 2198b438e21SMatthew G. Knepley PetscScalar result; 2208b438e21SMatthew G. Knepley PetscReal vol, tol = ctx->tol; 2218b438e21SMatthew G. Knepley 2228b438e21SMatthew G. Knepley PetscFunctionBeginUser; 2239566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &u)); 2249566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(dm, u, &result, ctx)); 2258b438e21SMatthew G. Knepley vol = PetscRealPart(result); 2269566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &u)); 2279566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Volume: %g\n", (double)vol)); 2288b438e21SMatthew G. Knepley if (ctx->volume > 0.0 && PetscAbsReal(ctx->volume - vol) > tol) { 22998921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Calculated volume %g != %g actual volume (error %g > %g tol)", (double)vol, (double)ctx->volume, (double)PetscAbsReal(ctx->volume - vol), (double)tol); 2308b438e21SMatthew G. Knepley } 2318b438e21SMatthew G. Knepley PetscFunctionReturn(0); 2328b438e21SMatthew G. Knepley } 2338b438e21SMatthew G. Knepley 234*9371c9d4SSatish Balay int main(int argc, char **argv) { 23539290ff6SMatthew G. Knepley DM dm; 2368b438e21SMatthew G. Knepley AppCtx user; 2378b438e21SMatthew G. Knepley 238327415f7SBarry Smith PetscFunctionBeginUser; 2399566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 2409566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 2419566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 2429566063dSJacob Faibussowitsch PetscCall(CreateDiscretization(dm, &user)); 2439566063dSJacob Faibussowitsch PetscCall(CheckVolume(dm, &user)); 2449566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 2459566063dSJacob Faibussowitsch PetscCall(PetscFree(user.transformDataReal)); 2469566063dSJacob Faibussowitsch PetscCall(PetscFree(user.transformData)); 2479566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 248b122ec5aSJacob Faibussowitsch return 0; 2498b438e21SMatthew G. Knepley } 2508b438e21SMatthew G. Knepley 2518b438e21SMatthew G. Knepley /*TEST 2528b438e21SMatthew G. Knepley 25330602db0SMatthew G. Knepley testset: 25430602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower -1.,-1. -dm_plex_box_upper 1.,1. -volume 4. -dm_coord_space 0 25530602db0SMatthew G. Knepley 2568b438e21SMatthew G. Knepley test: 2578b438e21SMatthew G. Knepley suffix: square_0 25830602db0SMatthew G. Knepley args: -dm_coord_petscspace_degree 1 2598b438e21SMatthew G. Knepley 2608b438e21SMatthew G. Knepley test: 2618b438e21SMatthew G. Knepley suffix: square_1 26230602db0SMatthew G. Knepley args: -dm_coord_petscspace_degree 2 2638b438e21SMatthew G. Knepley 2648b438e21SMatthew G. Knepley test: 2658b438e21SMatthew G. Knepley suffix: square_2 26630602db0SMatthew G. Knepley args: -dm_refine 1 -dm_coord_petscspace_degree 1 2678b438e21SMatthew G. Knepley 2688b438e21SMatthew G. Knepley test: 2698b438e21SMatthew G. Knepley suffix: square_3 27030602db0SMatthew G. Knepley args: -dm_refine 1 -dm_coord_petscspace_degree 2 27130602db0SMatthew G. Knepley 27230602db0SMatthew G. Knepley testset: 27330602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_faces 1,1,1 -dm_plex_box_lower -1.,-1.,-1. -dm_plex_box_upper 1.,1.,1. -volume 8. -dm_coord_space 0 2748b438e21SMatthew G. Knepley 2758b438e21SMatthew G. Knepley test: 2768b438e21SMatthew G. Knepley suffix: cube_0 27730602db0SMatthew G. Knepley args: -dm_coord_petscspace_degree 1 2788b438e21SMatthew G. Knepley 2798b438e21SMatthew G. Knepley test: 2808b438e21SMatthew G. Knepley suffix: cube_1 28130602db0SMatthew G. Knepley args: -dm_coord_petscspace_degree 2 2828b438e21SMatthew G. Knepley 2838b438e21SMatthew G. Knepley test: 2848b438e21SMatthew G. Knepley suffix: cube_2 28530602db0SMatthew G. Knepley args: -dm_refine 1 -dm_coord_petscspace_degree 1 2868b438e21SMatthew G. Knepley 2878b438e21SMatthew G. Knepley test: 2888b438e21SMatthew G. Knepley suffix: cube_3 28930602db0SMatthew G. Knepley args: -dm_refine 1 -dm_coord_petscspace_degree 2 29030602db0SMatthew G. Knepley 29130602db0SMatthew G. Knepley testset: 29230602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower -1.,-1. -dm_plex_box_upper 1.,1. -volume 4. -dm_coord_space 0 2938b438e21SMatthew G. Knepley 2948b438e21SMatthew G. Knepley test: 2958b438e21SMatthew G. Knepley suffix: shear_0 29630602db0SMatthew G. Knepley args: -dm_coord_petscspace_degree 1 -mesh_transform shear -transform_data 3.0 2978b438e21SMatthew G. Knepley 2988b438e21SMatthew G. Knepley test: 2998b438e21SMatthew G. Knepley suffix: shear_1 30030602db0SMatthew G. Knepley args: -dm_coord_petscspace_degree 2 -mesh_transform shear -transform_data 3.0 3018b438e21SMatthew G. Knepley 3028b438e21SMatthew G. Knepley test: 3038b438e21SMatthew G. Knepley suffix: shear_2 30430602db0SMatthew G. Knepley args: -dm_refine 1 -dm_coord_petscspace_degree 1 -mesh_transform shear -transform_data 3.0 3058b438e21SMatthew G. Knepley 3068b438e21SMatthew G. Knepley test: 3078b438e21SMatthew G. Knepley suffix: shear_3 30830602db0SMatthew G. Knepley args: -dm_refine 1 -dm_coord_petscspace_degree 2 -mesh_transform shear -transform_data 3.0 30930602db0SMatthew G. Knepley 31030602db0SMatthew G. Knepley testset: 31130602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_faces 1,1,1 -dm_plex_box_lower -1.,-1.,-1. -dm_plex_box_upper 1.,1.,1. -volume 8. -dm_coord_space 0 3128b438e21SMatthew G. Knepley 3138b438e21SMatthew G. Knepley test: 3148b438e21SMatthew G. Knepley suffix: shear_4 31530602db0SMatthew G. Knepley args: -dm_coord_petscspace_degree 1 -mesh_transform shear -transform_data 3.0 3168b438e21SMatthew G. Knepley 3178b438e21SMatthew G. Knepley test: 3188b438e21SMatthew G. Knepley suffix: shear_5 31930602db0SMatthew G. Knepley args: -dm_coord_petscspace_degree 2 -mesh_transform shear -transform_data 3.0 3208b438e21SMatthew G. Knepley 3218b438e21SMatthew G. Knepley test: 3228b438e21SMatthew G. Knepley suffix: shear_6 32330602db0SMatthew G. Knepley args: -dm_refine 1 -dm_coord_petscspace_degree 1 -mesh_transform shear -transform_data 3.0,4.0 3248b438e21SMatthew G. Knepley 3258b438e21SMatthew G. Knepley test: 3268b438e21SMatthew G. Knepley suffix: shear_7 32730602db0SMatthew G. Knepley args: -dm_refine 1 -dm_coord_petscspace_degree 2 -mesh_transform shear -transform_data 3.0,4.0 32830602db0SMatthew G. Knepley 32930602db0SMatthew G. Knepley testset: 33060c1a66aSMatthew G. Knepley args: -dm_coord_space 0 -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower 1.,-1. -dm_plex_box_upper 3.,1. \ 33160c1a66aSMatthew G. Knepley -dm_coord_petscspace_degree 1 -mesh_transform flare -volume 8. 33260c1a66aSMatthew G. Knepley 33360c1a66aSMatthew G. Knepley test: 33460c1a66aSMatthew G. Knepley suffix: flare_0 33560c1a66aSMatthew G. Knepley args: 33660c1a66aSMatthew G. Knepley 33760c1a66aSMatthew G. Knepley test: 33860c1a66aSMatthew G. Knepley suffix: flare_1 33960c1a66aSMatthew G. Knepley args: -dm_refine 2 34060c1a66aSMatthew G. Knepley 34160c1a66aSMatthew G. Knepley testset: 34230602db0SMatthew G. Knepley # Area: 3/4 \pi = 2.3562 34330602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -mesh_transform annulus -volume 2.35619449019235 -dm_coord_space 0 3448b438e21SMatthew G. Knepley 3458b438e21SMatthew G. Knepley test: 3468b438e21SMatthew G. Knepley # Area: (a+b)/2 h = 3/\sqrt{2} (sqrt{2} - 1/\sqrt{2}) = 3/2 3478b438e21SMatthew G. Knepley suffix: annulus_0 3488b438e21SMatthew G. Knepley requires: double 34930602db0SMatthew G. Knepley args: -dm_coord_petscspace_degree 1 -volume 1.5 3508b438e21SMatthew G. Knepley 3518b438e21SMatthew G. Knepley test: 3528b438e21SMatthew G. Knepley suffix: annulus_1 3538b438e21SMatthew G. Knepley requires: double 35430602db0SMatthew G. Knepley args: -dm_refine 3 -dm_coord_petscspace_degree 1 -tol .016 3558b438e21SMatthew G. Knepley 3568b438e21SMatthew G. Knepley test: 3578b438e21SMatthew G. Knepley suffix: annulus_2 3588b438e21SMatthew G. Knepley requires: double 35930602db0SMatthew G. Knepley args: -dm_refine 3 -dm_coord_petscspace_degree 2 -tol .0038 3608b438e21SMatthew G. Knepley 3618b438e21SMatthew G. Knepley test: 3628b438e21SMatthew G. Knepley suffix: annulus_3 3638b438e21SMatthew G. Knepley requires: double 36430602db0SMatthew G. Knepley args: -dm_refine 3 -dm_coord_petscspace_degree 3 -tol 2.2e-6 3658b438e21SMatthew G. Knepley 3668b438e21SMatthew G. Knepley test: 3678b438e21SMatthew G. Knepley suffix: annulus_4 3688b438e21SMatthew G. Knepley requires: double 36930602db0SMatthew G. Knepley args: -dm_refine 2 -dm_coord_petscspace_degree 2 -petscfe_default_quadrature_order 2 -tol .00012 3708b438e21SMatthew G. Knepley 3718b438e21SMatthew G. Knepley test: 3728b438e21SMatthew G. Knepley suffix: annulus_5 3738b438e21SMatthew G. Knepley requires: double 37430602db0SMatthew G. Knepley args: -dm_refine 2 -dm_coord_petscspace_degree 3 -petscfe_default_quadrature_order 3 -tol 1.2e-7 37530602db0SMatthew G. Knepley 37630602db0SMatthew G. Knepley testset: 37730602db0SMatthew G. Knepley # Volume: 4/3 \pi (8 - 1)/2 = 14/3 \pi = 14.66076571675238 37830602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_faces 1,1,1 -dm_plex_box_lower -1.,-1.,-1. -dm_plex_box_upper 1.,1.,1. -mesh_transform shell -volume 14.66076571675238 -dm_coord_space 0 3798b438e21SMatthew G. Knepley 3808b438e21SMatthew G. Knepley test: 3818b438e21SMatthew G. Knepley suffix: shell_0 3828b438e21SMatthew G. Knepley requires: double 38330602db0SMatthew G. Knepley args: -dm_refine 1 -dm_coord_petscspace_degree 1 -petscfe_default_quadrature_order 1 -volume 5.633164922 -tol 1.0e-7 3848b438e21SMatthew G. Knepley 3858b438e21SMatthew G. Knepley test: 3868b438e21SMatthew G. Knepley suffix: shell_1 3878b438e21SMatthew G. Knepley requires: double 38830602db0SMatthew G. Knepley args: -dm_refine 2 -dm_coord_petscspace_degree 1 -petscfe_default_quadrature_order 1 -tol 3.1 3898b438e21SMatthew G. Knepley 3908b438e21SMatthew G. Knepley test: 3918b438e21SMatthew G. Knepley suffix: shell_2 3928b438e21SMatthew G. Knepley requires: double 39330602db0SMatthew G. Knepley args: -dm_refine 2 -dm_coord_petscspace_degree 2 -petscfe_default_quadrature_order 2 -tol .1 3948b438e21SMatthew G. Knepley 3958b438e21SMatthew G. Knepley test: 3968b438e21SMatthew G. Knepley suffix: shell_3 3978b438e21SMatthew G. Knepley requires: double 39830602db0SMatthew G. Knepley args: -dm_refine 2 -dm_coord_petscspace_degree 3 -petscfe_default_quadrature_order 3 -tol .02 3998b438e21SMatthew G. Knepley 4008b438e21SMatthew G. Knepley test: 4018b438e21SMatthew G. Knepley suffix: shell_4 4028b438e21SMatthew G. Knepley requires: double 40330602db0SMatthew G. Knepley args: -dm_refine 2 -dm_coord_petscspace_degree 4 -petscfe_default_quadrature_order 4 -tol .006 40439290ff6SMatthew G. Knepley 40539290ff6SMatthew G. Knepley test: 40639290ff6SMatthew G. Knepley # Volume: 1.0 40739290ff6SMatthew G. Knepley suffix: gmsh_q2 40839290ff6SMatthew G. Knepley requires: double 40930602db0SMatthew G. Knepley args: -coord_space 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/quads-q2.msh -dm_plex_gmsh_project -volume 1.0 -tol 1e-6 41039290ff6SMatthew G. Knepley 41139290ff6SMatthew G. Knepley test: 41239290ff6SMatthew G. Knepley # Volume: 1.0 41339290ff6SMatthew G. Knepley suffix: gmsh_q3 41439290ff6SMatthew G. Knepley requires: double 41524b9a4b1SJed Brown nsize: {{1 2}} 41630602db0SMatthew G. Knepley args: -coord_space 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/quads-q3.msh -dm_plex_gmsh_project -volume 1.0 -tol 1e-6 4178b438e21SMatthew G. Knepley 4188b438e21SMatthew G. Knepley TEST*/ 419