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 69371c9d4SSatish Balay typedef enum { 79371c9d4SSatish Balay TRANSFORM_NONE, 89371c9d4SSatish Balay TRANSFORM_SHEAR, 99371c9d4SSatish Balay TRANSFORM_FLARE, 109371c9d4SSatish Balay TRANSFORM_ANNULUS, 119371c9d4SSatish Balay TRANSFORM_SHELL 129371c9d4SSatish 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 24d71ae5a4SJacob Faibussowitsch PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 25d71ae5a4SJacob Faibussowitsch { 268b438e21SMatthew G. Knepley PetscInt n = 0, i; 278b438e21SMatthew G. Knepley 288b438e21SMatthew G. Knepley PetscFunctionBegin; 2930602db0SMatthew G. Knepley options->coordSpace = PETSC_TRUE; 308b438e21SMatthew G. Knepley options->meshTransform = TRANSFORM_NONE; 318b438e21SMatthew G. Knepley options->transformDataReal = NULL; 328b438e21SMatthew G. Knepley options->transformData = NULL; 338b438e21SMatthew G. Knepley options->volume = -1.0; 348b438e21SMatthew G. Knepley options->tol = PETSC_SMALL; 358b438e21SMatthew G. Knepley 36d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Meshing Interpolation Test Options", "DMPLEX"); 379566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-coord_space", "Flag to create a coordinate space", "ex33.c", options->coordSpace, &options->coordSpace, NULL)); 389566063dSJacob 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)); 398b438e21SMatthew G. Knepley switch (options->meshTransform) { 40d71ae5a4SJacob Faibussowitsch case TRANSFORM_NONE: 41d71ae5a4SJacob Faibussowitsch break; 428b438e21SMatthew G. Knepley case TRANSFORM_SHEAR: 4339290ff6SMatthew G. Knepley n = 2; 449566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &options->transformDataReal)); 458b438e21SMatthew G. Knepley for (i = 0; i < n; ++i) options->transformDataReal[i] = 1.0; 469566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-transform_data", "Parameters for mesh transforms", "ex33.c", options->transformDataReal, &n, NULL)); 478b438e21SMatthew G. Knepley break; 4860c1a66aSMatthew G. Knepley case TRANSFORM_FLARE: 4960c1a66aSMatthew G. Knepley n = 4; 509566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &options->transformData)); 5160c1a66aSMatthew G. Knepley for (i = 0; i < n; ++i) options->transformData[i] = 1.0; 5260c1a66aSMatthew G. Knepley options->transformData[0] = (PetscScalar)0; 539566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalarArray("-transform_data", "Parameters for mesh transforms", "ex33.c", options->transformData, &n, NULL)); 5460c1a66aSMatthew G. Knepley break; 558b438e21SMatthew G. Knepley case TRANSFORM_ANNULUS: 568b438e21SMatthew G. Knepley n = 2; 579566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &options->transformData)); 588b438e21SMatthew G. Knepley options->transformData[0] = 1.0; 598b438e21SMatthew G. Knepley options->transformData[1] = 2.0; 609566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalarArray("-transform_data", "Parameters for mesh transforms", "ex33.c", options->transformData, &n, NULL)); 618b438e21SMatthew G. Knepley break; 628b438e21SMatthew G. Knepley case TRANSFORM_SHELL: 638b438e21SMatthew G. Knepley n = 2; 649566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &options->transformData)); 658b438e21SMatthew G. Knepley options->transformData[0] = 1.0; 668b438e21SMatthew G. Knepley options->transformData[1] = 2.0; 679566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalarArray("-transform_data", "Parameters for mesh transforms", "ex33.c", options->transformData, &n, NULL)); 688b438e21SMatthew G. Knepley break; 69d71ae5a4SJacob Faibussowitsch default: 70d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Unknown mesh transform %d", options->meshTransform); 718b438e21SMatthew G. Knepley } 729566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-volume", "The analytical volume of the mesh", "ex33.c", options->volume, &options->volume, NULL)); 739566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tol", "The tolerance for the volume check", "ex33.c", options->tol, &options->tol, NULL)); 74d0609cedSBarry Smith PetscOptionsEnd(); 753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 768b438e21SMatthew G. Knepley } 778b438e21SMatthew G. Knepley 78d71ae5a4SJacob Faibussowitsch 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[]) 79d71ae5a4SJacob Faibussowitsch { 808b438e21SMatthew G. Knepley const PetscInt Nc = uOff[1] - uOff[0]; 818b438e21SMatthew G. Knepley PetscInt c; 828b438e21SMatthew G. Knepley 838b438e21SMatthew G. Knepley for (c = 0; c < Nc; ++c) f0[c] = u[c]; 848b438e21SMatthew G. Knepley } 858b438e21SMatthew G. Knepley 8660c1a66aSMatthew G. Knepley /* Flare applies the transformation, assuming we fix x_f, 8760c1a66aSMatthew G. Knepley 8860c1a66aSMatthew G. Knepley x_i = x_i * alpha_i x_f 8960c1a66aSMatthew G. Knepley */ 90d71ae5a4SJacob Faibussowitsch 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[]) 91d71ae5a4SJacob Faibussowitsch { 9260c1a66aSMatthew G. Knepley const PetscInt Nc = uOff[1] - uOff[0]; 9360c1a66aSMatthew G. Knepley const PetscInt cf = (PetscInt)PetscRealPart(constants[0]); 9460c1a66aSMatthew G. Knepley PetscInt c; 9560c1a66aSMatthew G. Knepley 96ad540459SPierre Jolivet for (c = 0; c < Nc; ++c) coords[c] = u[c] * (c == cf ? 1.0 : constants[c + 1] * u[cf]); 9760c1a66aSMatthew G. Knepley } 9860c1a66aSMatthew G. Knepley 998b438e21SMatthew G. Knepley /* 1008b438e21SMatthew 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 1018b438e21SMatthew G. Knepley will correspond to the top and bottom of our square. So 1028b438e21SMatthew G. Knepley 1038b438e21SMatthew G. Knepley (0,0)--(1,0) ==> (1,0)--(2,0) Just a shift of (1,0) 1048b438e21SMatthew G. Knepley (0,1)--(1,1) ==> (0,1)--(0,2) Switch x and y 1058b438e21SMatthew G. Knepley 1068b438e21SMatthew 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: 1078b438e21SMatthew G. Knepley 1088b438e21SMatthew G. Knepley (x, y) ==> (x+1, \pi/2 y) in (r', \theta') space 1098b438e21SMatthew G. Knepley ==> ((x+1) cos(\pi/2 y), (x+1) sin(\pi/2 y)) in (x', y') space 1108b438e21SMatthew G. Knepley */ 111d71ae5a4SJacob Faibussowitsch 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[]) 112d71ae5a4SJacob Faibussowitsch { 1138b438e21SMatthew G. Knepley const PetscReal ri = PetscRealPart(constants[0]); 1148b438e21SMatthew G. Knepley const PetscReal ro = PetscRealPart(constants[1]); 1158b438e21SMatthew G. Knepley 1168b438e21SMatthew G. Knepley xp[0] = (x[0] * (ro - ri) + ri) * PetscCosReal(0.5 * PETSC_PI * x[1]); 1178b438e21SMatthew G. Knepley xp[1] = (x[0] * (ro - ri) + ri) * PetscSinReal(0.5 * PETSC_PI * x[1]); 1188b438e21SMatthew G. Knepley } 1198b438e21SMatthew G. Knepley 1208b438e21SMatthew G. Knepley /* 1218b438e21SMatthew 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 1228b438e21SMatthew G. Knepley lower hemisphere and the upper surface onto the top, letting z be the radius. 1238b438e21SMatthew G. Knepley 1248b438e21SMatthew G. Knepley (x, y) ==> ((z+3)/2, \pi/2 (|x| or |y|), arctan y/x) in (r', \theta', \phi') space 1258b438e21SMatthew 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 1268b438e21SMatthew G. Knepley */ 127d71ae5a4SJacob Faibussowitsch 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[]) 128d71ae5a4SJacob Faibussowitsch { 1298b438e21SMatthew G. Knepley const PetscReal pi4 = PETSC_PI / 4.0; 1308b438e21SMatthew G. Knepley const PetscReal ri = PetscRealPart(constants[0]); 1318b438e21SMatthew G. Knepley const PetscReal ro = PetscRealPart(constants[1]); 1328b438e21SMatthew G. Knepley const PetscReal rp = (x[2] + 1) * 0.5 * (ro - ri) + ri; 1338b438e21SMatthew G. Knepley const PetscReal phip = PetscAtan2Real(x[1], x[0]); 1348b438e21SMatthew 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]))); 1358b438e21SMatthew G. Knepley 1368b438e21SMatthew G. Knepley xp[0] = rp * PetscCosReal(thetap) * PetscCosReal(phip); 1378b438e21SMatthew G. Knepley xp[1] = rp * PetscCosReal(thetap) * PetscSinReal(phip); 1388b438e21SMatthew G. Knepley xp[2] = rp * PetscSinReal(thetap); 1398b438e21SMatthew G. Knepley } 1408b438e21SMatthew G. Knepley 141d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCreateCoordinateDisc(DM dm) 142d71ae5a4SJacob Faibussowitsch { 1438b438e21SMatthew G. Knepley DM cdm; 1448b438e21SMatthew G. Knepley PetscFE fe; 14539290ff6SMatthew G. Knepley DMPolytopeType ct; 14639290ff6SMatthew G. Knepley PetscInt dim, dE, cStart; 14739290ff6SMatthew G. Knepley PetscBool simplex; 1488b438e21SMatthew G. Knepley 1498b438e21SMatthew G. Knepley PetscFunctionBegin; 1509566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &cdm)); 1519566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1529566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &dE)); 1539566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(cdm, 0, &cStart, NULL)); 1549566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 15539290ff6SMatthew G. Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE; 1569566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, dE, simplex, "dm_coord_", -1, &fe)); 1579566063dSJacob Faibussowitsch PetscCall(DMProjectCoordinates(dm, fe)); 1589566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 1593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1608b438e21SMatthew G. Knepley } 1618b438e21SMatthew G. Knepley 162d71ae5a4SJacob Faibussowitsch PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm) 163d71ae5a4SJacob Faibussowitsch { 1648b438e21SMatthew G. Knepley DM cdm; 1658b438e21SMatthew G. Knepley PetscDS cds; 1668b438e21SMatthew G. Knepley 1678b438e21SMatthew G. Knepley PetscFunctionBegin; 1689566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 1699566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 1709566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 17130602db0SMatthew G. Knepley 1729566063dSJacob Faibussowitsch if (ctx->coordSpace) PetscCall(DMCreateCoordinateDisc(*dm)); 1738b438e21SMatthew G. Knepley switch (ctx->meshTransform) { 174d71ae5a4SJacob Faibussowitsch case TRANSFORM_NONE: 175d71ae5a4SJacob Faibussowitsch PetscCall(DMPlexRemapGeometry(*dm, 0.0, identity)); 176d71ae5a4SJacob Faibussowitsch break; 177d71ae5a4SJacob Faibussowitsch case TRANSFORM_SHEAR: 178d71ae5a4SJacob Faibussowitsch PetscCall(DMPlexShearGeometry(*dm, DM_X, ctx->transformDataReal)); 179d71ae5a4SJacob Faibussowitsch break; 18060c1a66aSMatthew G. Knepley case TRANSFORM_FLARE: 1819566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(*dm, &cdm)); 1829566063dSJacob Faibussowitsch PetscCall(DMGetDS(cdm, &cds)); 1839566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(cds, 4, ctx->transformData)); 1849566063dSJacob Faibussowitsch PetscCall(DMPlexRemapGeometry(*dm, 0.0, f0_flare)); 18560c1a66aSMatthew G. Knepley break; 1868b438e21SMatthew G. Knepley case TRANSFORM_ANNULUS: 1879566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(*dm, &cdm)); 1889566063dSJacob Faibussowitsch PetscCall(DMGetDS(cdm, &cds)); 1899566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(cds, 2, ctx->transformData)); 1909566063dSJacob Faibussowitsch PetscCall(DMPlexRemapGeometry(*dm, 0.0, f0_annulus)); 1918b438e21SMatthew G. Knepley break; 1928b438e21SMatthew G. Knepley case TRANSFORM_SHELL: 1939566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(*dm, &cdm)); 1949566063dSJacob Faibussowitsch PetscCall(DMGetDS(cdm, &cds)); 1959566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(cds, 2, ctx->transformData)); 1969566063dSJacob Faibussowitsch PetscCall(DMPlexRemapGeometry(*dm, 0.0, f0_shell)); 1978b438e21SMatthew G. Knepley break; 198d71ae5a4SJacob Faibussowitsch default: 199d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Unknown mesh transform %d", ctx->meshTransform); 2008b438e21SMatthew G. Knepley } 2019566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 2023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2038b438e21SMatthew G. Knepley } 2048b438e21SMatthew G. Knepley 205d71ae5a4SJacob Faibussowitsch 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[]) 206d71ae5a4SJacob Faibussowitsch { 2078b438e21SMatthew G. Knepley vol[0] = 1.; 2088b438e21SMatthew G. Knepley } 2098b438e21SMatthew G. Knepley 210d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateDiscretization(DM dm, AppCtx *ctx) 211d71ae5a4SJacob Faibussowitsch { 2128b438e21SMatthew G. Knepley PetscDS ds; 2138b438e21SMatthew G. Knepley PetscFE fe; 21439290ff6SMatthew G. Knepley DMPolytopeType ct; 21539290ff6SMatthew G. Knepley PetscInt dim, cStart; 21639290ff6SMatthew G. Knepley PetscBool simplex; 2178b438e21SMatthew G. Knepley 2188b438e21SMatthew G. Knepley PetscFunctionBeginUser; 2199566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 2209566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 2219566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 22239290ff6SMatthew G. Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE; 2239566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, NULL, PETSC_DETERMINE, &fe)); 2249566063dSJacob Faibussowitsch PetscCall(PetscFESetName(fe, "scalar")); 2259566063dSJacob Faibussowitsch PetscCall(DMAddField(dm, NULL, (PetscObject)fe)); 2269566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 2279566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 2289566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 2299566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(ds, 0, volume)); 2303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2318b438e21SMatthew G. Knepley } 2328b438e21SMatthew G. Knepley 233d71ae5a4SJacob Faibussowitsch static PetscErrorCode CheckVolume(DM dm, AppCtx *ctx) 234d71ae5a4SJacob Faibussowitsch { 2358b438e21SMatthew G. Knepley Vec u; 2368b438e21SMatthew G. Knepley PetscScalar result; 2378b438e21SMatthew G. Knepley PetscReal vol, tol = ctx->tol; 2388b438e21SMatthew G. Knepley 2398b438e21SMatthew G. Knepley PetscFunctionBeginUser; 2409566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &u)); 2419566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(dm, u, &result, ctx)); 2428b438e21SMatthew G. Knepley vol = PetscRealPart(result); 2439566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &u)); 2449566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Volume: %g\n", (double)vol)); 2458b438e21SMatthew G. Knepley if (ctx->volume > 0.0 && PetscAbsReal(ctx->volume - vol) > tol) { 24698921bdaSJacob 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); 2478b438e21SMatthew G. Knepley } 2483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2498b438e21SMatthew G. Knepley } 2508b438e21SMatthew G. Knepley 251d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 252d71ae5a4SJacob Faibussowitsch { 25339290ff6SMatthew G. Knepley DM dm; 2548b438e21SMatthew G. Knepley AppCtx user; 2558b438e21SMatthew G. Knepley 256327415f7SBarry Smith PetscFunctionBeginUser; 2579566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 2589566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 2599566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 2609566063dSJacob Faibussowitsch PetscCall(CreateDiscretization(dm, &user)); 2619566063dSJacob Faibussowitsch PetscCall(CheckVolume(dm, &user)); 2629566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 2639566063dSJacob Faibussowitsch PetscCall(PetscFree(user.transformDataReal)); 2649566063dSJacob Faibussowitsch PetscCall(PetscFree(user.transformData)); 2659566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 266b122ec5aSJacob Faibussowitsch return 0; 2678b438e21SMatthew G. Knepley } 2688b438e21SMatthew G. Knepley 2698b438e21SMatthew G. Knepley /*TEST 2708b438e21SMatthew G. Knepley 27130602db0SMatthew G. Knepley testset: 27230602db0SMatthew 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 27330602db0SMatthew G. Knepley 2748b438e21SMatthew G. Knepley test: 2758b438e21SMatthew G. Knepley suffix: square_0 27630602db0SMatthew G. Knepley args: -dm_coord_petscspace_degree 1 2778b438e21SMatthew G. Knepley 2788b438e21SMatthew G. Knepley test: 2798b438e21SMatthew G. Knepley suffix: square_1 28030602db0SMatthew G. Knepley args: -dm_coord_petscspace_degree 2 2818b438e21SMatthew G. Knepley 2828b438e21SMatthew G. Knepley test: 2838b438e21SMatthew G. Knepley suffix: square_2 28430602db0SMatthew G. Knepley args: -dm_refine 1 -dm_coord_petscspace_degree 1 2858b438e21SMatthew G. Knepley 2868b438e21SMatthew G. Knepley test: 2878b438e21SMatthew G. Knepley suffix: square_3 28830602db0SMatthew G. Knepley args: -dm_refine 1 -dm_coord_petscspace_degree 2 28930602db0SMatthew G. Knepley 29030602db0SMatthew G. Knepley testset: 29130602db0SMatthew 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 2928b438e21SMatthew G. Knepley 2938b438e21SMatthew G. Knepley test: 2948b438e21SMatthew G. Knepley suffix: cube_0 29530602db0SMatthew G. Knepley args: -dm_coord_petscspace_degree 1 2968b438e21SMatthew G. Knepley 2978b438e21SMatthew G. Knepley test: 2988b438e21SMatthew G. Knepley suffix: cube_1 29930602db0SMatthew G. Knepley args: -dm_coord_petscspace_degree 2 3008b438e21SMatthew G. Knepley 3018b438e21SMatthew G. Knepley test: 3028b438e21SMatthew G. Knepley suffix: cube_2 30330602db0SMatthew G. Knepley args: -dm_refine 1 -dm_coord_petscspace_degree 1 3048b438e21SMatthew G. Knepley 3058b438e21SMatthew G. Knepley test: 3068b438e21SMatthew G. Knepley suffix: cube_3 30730602db0SMatthew G. Knepley args: -dm_refine 1 -dm_coord_petscspace_degree 2 30830602db0SMatthew G. Knepley 30930602db0SMatthew G. Knepley testset: 31030602db0SMatthew 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 3118b438e21SMatthew G. Knepley 3128b438e21SMatthew G. Knepley test: 3138b438e21SMatthew G. Knepley suffix: shear_0 31430602db0SMatthew G. Knepley args: -dm_coord_petscspace_degree 1 -mesh_transform shear -transform_data 3.0 3158b438e21SMatthew G. Knepley 3168b438e21SMatthew G. Knepley test: 3178b438e21SMatthew G. Knepley suffix: shear_1 31830602db0SMatthew G. Knepley args: -dm_coord_petscspace_degree 2 -mesh_transform shear -transform_data 3.0 3198b438e21SMatthew G. Knepley 3208b438e21SMatthew G. Knepley test: 3218b438e21SMatthew G. Knepley suffix: shear_2 32230602db0SMatthew G. Knepley args: -dm_refine 1 -dm_coord_petscspace_degree 1 -mesh_transform shear -transform_data 3.0 3238b438e21SMatthew G. Knepley 3248b438e21SMatthew G. Knepley test: 3258b438e21SMatthew G. Knepley suffix: shear_3 32630602db0SMatthew G. Knepley args: -dm_refine 1 -dm_coord_petscspace_degree 2 -mesh_transform shear -transform_data 3.0 32730602db0SMatthew G. Knepley 32830602db0SMatthew G. Knepley testset: 32930602db0SMatthew 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 3308b438e21SMatthew G. Knepley 3318b438e21SMatthew G. Knepley test: 3328b438e21SMatthew G. Knepley suffix: shear_4 33330602db0SMatthew G. Knepley args: -dm_coord_petscspace_degree 1 -mesh_transform shear -transform_data 3.0 3348b438e21SMatthew G. Knepley 3358b438e21SMatthew G. Knepley test: 3368b438e21SMatthew G. Knepley suffix: shear_5 33730602db0SMatthew G. Knepley args: -dm_coord_petscspace_degree 2 -mesh_transform shear -transform_data 3.0 3388b438e21SMatthew G. Knepley 3398b438e21SMatthew G. Knepley test: 3408b438e21SMatthew G. Knepley suffix: shear_6 34130602db0SMatthew G. Knepley args: -dm_refine 1 -dm_coord_petscspace_degree 1 -mesh_transform shear -transform_data 3.0,4.0 3428b438e21SMatthew G. Knepley 3438b438e21SMatthew G. Knepley test: 3448b438e21SMatthew G. Knepley suffix: shear_7 34530602db0SMatthew G. Knepley args: -dm_refine 1 -dm_coord_petscspace_degree 2 -mesh_transform shear -transform_data 3.0,4.0 34630602db0SMatthew G. Knepley 34730602db0SMatthew G. Knepley testset: 34860c1a66aSMatthew 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. \ 34960c1a66aSMatthew G. Knepley -dm_coord_petscspace_degree 1 -mesh_transform flare -volume 8. 35060c1a66aSMatthew G. Knepley 35160c1a66aSMatthew G. Knepley test: 35260c1a66aSMatthew G. Knepley suffix: flare_0 35360c1a66aSMatthew G. Knepley args: 35460c1a66aSMatthew G. Knepley 35560c1a66aSMatthew G. Knepley test: 35660c1a66aSMatthew G. Knepley suffix: flare_1 35760c1a66aSMatthew G. Knepley args: -dm_refine 2 35860c1a66aSMatthew G. Knepley 35960c1a66aSMatthew G. Knepley testset: 36030602db0SMatthew G. Knepley # Area: 3/4 \pi = 2.3562 36130602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -mesh_transform annulus -volume 2.35619449019235 -dm_coord_space 0 3628b438e21SMatthew G. Knepley 3638b438e21SMatthew G. Knepley test: 3648b438e21SMatthew G. Knepley # Area: (a+b)/2 h = 3/\sqrt{2} (sqrt{2} - 1/\sqrt{2}) = 3/2 3658b438e21SMatthew G. Knepley suffix: annulus_0 3668b438e21SMatthew G. Knepley requires: double 36730602db0SMatthew G. Knepley args: -dm_coord_petscspace_degree 1 -volume 1.5 3688b438e21SMatthew G. Knepley 3698b438e21SMatthew G. Knepley test: 3708b438e21SMatthew G. Knepley suffix: annulus_1 3718b438e21SMatthew G. Knepley requires: double 37230602db0SMatthew G. Knepley args: -dm_refine 3 -dm_coord_petscspace_degree 1 -tol .016 3738b438e21SMatthew G. Knepley 3748b438e21SMatthew G. Knepley test: 3758b438e21SMatthew G. Knepley suffix: annulus_2 3768b438e21SMatthew G. Knepley requires: double 37730602db0SMatthew G. Knepley args: -dm_refine 3 -dm_coord_petscspace_degree 2 -tol .0038 3788b438e21SMatthew G. Knepley 3798b438e21SMatthew G. Knepley test: 3808b438e21SMatthew G. Knepley suffix: annulus_3 3818b438e21SMatthew G. Knepley requires: double 38230602db0SMatthew G. Knepley args: -dm_refine 3 -dm_coord_petscspace_degree 3 -tol 2.2e-6 3838b438e21SMatthew G. Knepley 3848b438e21SMatthew G. Knepley test: 3858b438e21SMatthew G. Knepley suffix: annulus_4 3868b438e21SMatthew G. Knepley requires: double 38730602db0SMatthew G. Knepley args: -dm_refine 2 -dm_coord_petscspace_degree 2 -petscfe_default_quadrature_order 2 -tol .00012 3888b438e21SMatthew G. Knepley 3898b438e21SMatthew G. Knepley test: 3908b438e21SMatthew G. Knepley suffix: annulus_5 3918b438e21SMatthew G. Knepley requires: double 39230602db0SMatthew G. Knepley args: -dm_refine 2 -dm_coord_petscspace_degree 3 -petscfe_default_quadrature_order 3 -tol 1.2e-7 39330602db0SMatthew G. Knepley 39430602db0SMatthew G. Knepley testset: 39530602db0SMatthew G. Knepley # Volume: 4/3 \pi (8 - 1)/2 = 14/3 \pi = 14.66076571675238 39630602db0SMatthew 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 3978b438e21SMatthew G. Knepley 3988b438e21SMatthew G. Knepley test: 3998b438e21SMatthew G. Knepley suffix: shell_0 4008b438e21SMatthew G. Knepley requires: double 40130602db0SMatthew G. Knepley args: -dm_refine 1 -dm_coord_petscspace_degree 1 -petscfe_default_quadrature_order 1 -volume 5.633164922 -tol 1.0e-7 4028b438e21SMatthew G. Knepley 4038b438e21SMatthew G. Knepley test: 4048b438e21SMatthew G. Knepley suffix: shell_1 4058b438e21SMatthew G. Knepley requires: double 40630602db0SMatthew G. Knepley args: -dm_refine 2 -dm_coord_petscspace_degree 1 -petscfe_default_quadrature_order 1 -tol 3.1 4078b438e21SMatthew G. Knepley 4088b438e21SMatthew G. Knepley test: 4098b438e21SMatthew G. Knepley suffix: shell_2 4108b438e21SMatthew G. Knepley requires: double 41130602db0SMatthew G. Knepley args: -dm_refine 2 -dm_coord_petscspace_degree 2 -petscfe_default_quadrature_order 2 -tol .1 4128b438e21SMatthew G. Knepley 4138b438e21SMatthew G. Knepley test: 4148b438e21SMatthew G. Knepley suffix: shell_3 4158b438e21SMatthew G. Knepley requires: double 41630602db0SMatthew G. Knepley args: -dm_refine 2 -dm_coord_petscspace_degree 3 -petscfe_default_quadrature_order 3 -tol .02 4178b438e21SMatthew G. Knepley 4188b438e21SMatthew G. Knepley test: 4198b438e21SMatthew G. Knepley suffix: shell_4 4208b438e21SMatthew G. Knepley requires: double 42130602db0SMatthew G. Knepley args: -dm_refine 2 -dm_coord_petscspace_degree 4 -petscfe_default_quadrature_order 4 -tol .006 42239290ff6SMatthew G. Knepley 42339290ff6SMatthew G. Knepley test: 42439290ff6SMatthew G. Knepley # Volume: 1.0 42539290ff6SMatthew G. Knepley suffix: gmsh_q2 42639290ff6SMatthew G. Knepley requires: double 42730602db0SMatthew 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 42839290ff6SMatthew G. Knepley 42939290ff6SMatthew G. Knepley test: 43039290ff6SMatthew G. Knepley # Volume: 1.0 43139290ff6SMatthew G. Knepley suffix: gmsh_q3 43239290ff6SMatthew G. Knepley requires: double 43324b9a4b1SJed Brown nsize: {{1 2}} 43430602db0SMatthew 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 4358b438e21SMatthew G. Knepley 436*396e07e5SMatthew G. Knepley test: 437*396e07e5SMatthew G. Knepley # Volume: 1.0 438*396e07e5SMatthew G. Knepley suffix: gmsh_3d_q2 439*396e07e5SMatthew G. Knepley requires: double 440*396e07e5SMatthew G. Knepley args: -coord_space 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/cube_q2.msh -dm_plex_gmsh_project -volume 1.0 -tol 1e-6 441*396e07e5SMatthew G. Knepley 442*396e07e5SMatthew G. Knepley test: 443*396e07e5SMatthew G. Knepley # Volume: 1.0 444*396e07e5SMatthew G. Knepley suffix: gmsh_3d_q3 445*396e07e5SMatthew G. Knepley requires: double 446*396e07e5SMatthew G. Knepley nsize: {{1 2}} 447*396e07e5SMatthew G. Knepley args: -coord_space 0 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/cube_q3.msh -dm_plex_gmsh_project -volume 1.0 -tol 1e-6 448*396e07e5SMatthew G. Knepley 4498b438e21SMatthew G. Knepley TEST*/ 450