xref: /petsc/src/dm/impls/plex/tests/ex33.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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 
660c1a66aSMatthew G. Knepley typedef enum {TRANSFORM_NONE, TRANSFORM_SHEAR, TRANSFORM_FLARE, TRANSFORM_ANNULUS, TRANSFORM_SHELL} Transform;
760c1a66aSMatthew G. Knepley const char * const TransformTypes[] = {"none", "shear", "flare", "annulus", "shell", "Mesh Transform", "TRANSFORM_", NULL};
88b438e21SMatthew G. Knepley 
98b438e21SMatthew G. Knepley typedef struct {
1030602db0SMatthew G. Knepley   PetscBool   coordSpace;         /* Flag to create coordinate space */
118b438e21SMatthew G. Knepley   Transform   meshTransform;      /* Transform for initial box mesh */
128b438e21SMatthew G. Knepley   PetscReal   *transformDataReal; /* Parameters for mesh transform */
138b438e21SMatthew G. Knepley   PetscScalar *transformData;     /* Parameters for mesh transform */
148b438e21SMatthew G. Knepley   PetscReal   volume;             /* Analytical volume of the mesh */
158b438e21SMatthew G. Knepley   PetscReal   tol;                /* Tolerance for volume check */
168b438e21SMatthew G. Knepley } AppCtx;
178b438e21SMatthew G. Knepley 
188b438e21SMatthew G. Knepley PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
198b438e21SMatthew G. Knepley {
208b438e21SMatthew G. Knepley   PetscInt       n = 0, i;
218b438e21SMatthew G. Knepley 
228b438e21SMatthew G. Knepley   PetscFunctionBegin;
2330602db0SMatthew G. Knepley   options->coordSpace        = PETSC_TRUE;
248b438e21SMatthew G. Knepley   options->meshTransform     = TRANSFORM_NONE;
258b438e21SMatthew G. Knepley   options->transformDataReal = NULL;
268b438e21SMatthew G. Knepley   options->transformData     = NULL;
278b438e21SMatthew G. Knepley   options->volume            = -1.0;
288b438e21SMatthew G. Knepley   options->tol               = PETSC_SMALL;
298b438e21SMatthew G. Knepley 
30d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Meshing Interpolation Test Options", "DMPLEX");
319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-coord_space", "Flag to create a coordinate space", "ex33.c", options->coordSpace, &options->coordSpace, NULL));
329566063dSJacob 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));
338b438e21SMatthew G. Knepley   switch (options->meshTransform) {
348b438e21SMatthew G. Knepley     case TRANSFORM_NONE: break;
358b438e21SMatthew G. Knepley     case TRANSFORM_SHEAR:
3639290ff6SMatthew G. Knepley       n = 2;
379566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(n, &options->transformDataReal));
388b438e21SMatthew G. Knepley       for (i = 0; i < n; ++i) options->transformDataReal[i] = 1.0;
399566063dSJacob Faibussowitsch       PetscCall(PetscOptionsRealArray("-transform_data", "Parameters for mesh transforms", "ex33.c", options->transformDataReal, &n, NULL));
408b438e21SMatthew G. Knepley       break;
4160c1a66aSMatthew G. Knepley     case TRANSFORM_FLARE:
4260c1a66aSMatthew G. Knepley       n = 4;
439566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(n, &options->transformData));
4460c1a66aSMatthew G. Knepley       for (i = 0; i < n; ++i) options->transformData[i] = 1.0;
4560c1a66aSMatthew G. Knepley       options->transformData[0] = (PetscScalar) 0;
469566063dSJacob Faibussowitsch       PetscCall(PetscOptionsScalarArray("-transform_data", "Parameters for mesh transforms", "ex33.c", options->transformData, &n, NULL));
4760c1a66aSMatthew G. Knepley       break;
488b438e21SMatthew G. Knepley     case TRANSFORM_ANNULUS:
498b438e21SMatthew G. Knepley       n = 2;
509566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(n, &options->transformData));
518b438e21SMatthew G. Knepley       options->transformData[0] = 1.0;
528b438e21SMatthew G. Knepley       options->transformData[1] = 2.0;
539566063dSJacob Faibussowitsch       PetscCall(PetscOptionsScalarArray("-transform_data", "Parameters for mesh transforms", "ex33.c", options->transformData, &n, NULL));
548b438e21SMatthew G. Knepley       break;
558b438e21SMatthew G. Knepley     case TRANSFORM_SHELL:
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;
6263a3b9bcSJacob Faibussowitsch     default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Unknown mesh transform %d", options->meshTransform);
638b438e21SMatthew G. Knepley   }
649566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-volume", "The analytical volume of the mesh", "ex33.c", options->volume, &options->volume, NULL));
659566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tol", "The tolerance for the volume check", "ex33.c", options->tol, &options->tol, NULL));
66d0609cedSBarry Smith   PetscOptionsEnd();
678b438e21SMatthew G. Knepley   PetscFunctionReturn(0);
688b438e21SMatthew G. Knepley }
698b438e21SMatthew G. Knepley 
708b438e21SMatthew G. Knepley static void identity(PetscInt dim, PetscInt Nf, PetscInt NfAux,
718b438e21SMatthew G. Knepley                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
728b438e21SMatthew G. Knepley                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
738b438e21SMatthew G. Knepley                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
748b438e21SMatthew G. Knepley {
758b438e21SMatthew G. Knepley   const PetscInt Nc = uOff[1] - uOff[0];
768b438e21SMatthew G. Knepley   PetscInt       c;
778b438e21SMatthew G. Knepley 
788b438e21SMatthew G. Knepley   for (c = 0; c < Nc; ++c) f0[c] = u[c];
798b438e21SMatthew G. Knepley }
808b438e21SMatthew G. Knepley 
8160c1a66aSMatthew G. Knepley /* Flare applies the transformation, assuming we fix x_f,
8260c1a66aSMatthew G. Knepley 
8360c1a66aSMatthew G. Knepley    x_i = x_i * alpha_i x_f
8460c1a66aSMatthew G. Knepley */
8560c1a66aSMatthew G. Knepley static void f0_flare(PetscInt dim, PetscInt Nf, PetscInt NfAux,
8660c1a66aSMatthew G. Knepley                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
8760c1a66aSMatthew G. Knepley                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
8860c1a66aSMatthew G. Knepley                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar coords[])
8960c1a66aSMatthew G. Knepley {
9060c1a66aSMatthew G. Knepley   const PetscInt Nc = uOff[1]-uOff[0];
9160c1a66aSMatthew G. Knepley   const PetscInt cf = (PetscInt) PetscRealPart(constants[0]);
9260c1a66aSMatthew G. Knepley   PetscInt       c;
9360c1a66aSMatthew G. Knepley 
9460c1a66aSMatthew G. Knepley   for (c = 0; c < Nc; ++c) {
9560c1a66aSMatthew G. Knepley     coords[c] = u[c] * (c == cf ? 1.0 : constants[c+1]*u[cf]);
9660c1a66aSMatthew G. Knepley   }
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 */
1118b438e21SMatthew G. Knepley static void f0_annulus(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1128b438e21SMatthew G. Knepley                        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1138b438e21SMatthew G. Knepley                        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1148b438e21SMatthew G. Knepley                        PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar xp[])
1158b438e21SMatthew G. Knepley {
1168b438e21SMatthew G. Knepley   const PetscReal ri = PetscRealPart(constants[0]);
1178b438e21SMatthew G. Knepley   const PetscReal ro = PetscRealPart(constants[1]);
1188b438e21SMatthew G. Knepley 
1198b438e21SMatthew G. Knepley   xp[0] = (x[0] * (ro-ri) + ri) * PetscCosReal(0.5*PETSC_PI*x[1]);
1208b438e21SMatthew G. Knepley   xp[1] = (x[0] * (ro-ri) + ri) * PetscSinReal(0.5*PETSC_PI*x[1]);
1218b438e21SMatthew G. Knepley }
1228b438e21SMatthew G. Knepley 
1238b438e21SMatthew G. Knepley /*
1248b438e21SMatthew 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
1258b438e21SMatthew G. Knepley   lower hemisphere and the upper surface onto the top, letting z be the radius.
1268b438e21SMatthew G. Knepley 
1278b438e21SMatthew G. Knepley     (x, y)  ==>  ((z+3)/2, \pi/2 (|x| or |y|), arctan y/x)                                                  in (r', \theta', \phi') space
1288b438e21SMatthew 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
1298b438e21SMatthew G. Knepley */
1308b438e21SMatthew G. Knepley static void f0_shell(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1318b438e21SMatthew G. Knepley                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1328b438e21SMatthew G. Knepley                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1338b438e21SMatthew G. Knepley                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar xp[])
1348b438e21SMatthew G. Knepley {
1358b438e21SMatthew G. Knepley   const PetscReal pi4    = PETSC_PI/4.0;
1368b438e21SMatthew G. Knepley   const PetscReal ri     = PetscRealPart(constants[0]);
1378b438e21SMatthew G. Knepley   const PetscReal ro     = PetscRealPart(constants[1]);
1388b438e21SMatthew G. Knepley   const PetscReal rp     = (x[2]+1) * 0.5*(ro-ri) + ri;
1398b438e21SMatthew G. Knepley   const PetscReal phip   = PetscAtan2Real(x[1], x[0]);
1408b438e21SMatthew 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])));
1418b438e21SMatthew G. Knepley 
1428b438e21SMatthew G. Knepley   xp[0] = rp * PetscCosReal(thetap) * PetscCosReal(phip);
1438b438e21SMatthew G. Knepley   xp[1] = rp * PetscCosReal(thetap) * PetscSinReal(phip);
1448b438e21SMatthew G. Knepley   xp[2] = rp * PetscSinReal(thetap);
1458b438e21SMatthew G. Knepley }
1468b438e21SMatthew G. Knepley 
14739290ff6SMatthew G. Knepley static PetscErrorCode DMCreateCoordinateDisc(DM dm)
1488b438e21SMatthew G. Knepley {
1498b438e21SMatthew G. Knepley   DM             cdm;
1508b438e21SMatthew G. Knepley   PetscFE        fe;
15139290ff6SMatthew G. Knepley   DMPolytopeType ct;
15239290ff6SMatthew G. Knepley   PetscInt       dim, dE, cStart;
15339290ff6SMatthew G. Knepley   PetscBool      simplex;
1548b438e21SMatthew G. Knepley 
1558b438e21SMatthew G. Knepley   PetscFunctionBegin;
1569566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(dm, &cdm));
1579566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1589566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &dE));
1599566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(cdm, 0, &cStart, NULL));
1609566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
16139290ff6SMatthew G. Knepley   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
1629566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, dE, simplex, "dm_coord_", -1, &fe));
1639566063dSJacob Faibussowitsch   PetscCall(DMProjectCoordinates(dm, fe));
1649566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
1658b438e21SMatthew G. Knepley   PetscFunctionReturn(0);
1668b438e21SMatthew G. Knepley }
1678b438e21SMatthew G. Knepley 
1688b438e21SMatthew G. Knepley PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm)
1698b438e21SMatthew G. Knepley {
1708b438e21SMatthew G. Knepley   DM             cdm;
1718b438e21SMatthew G. Knepley   PetscDS        cds;
1728b438e21SMatthew G. Knepley 
1738b438e21SMatthew G. Knepley   PetscFunctionBegin;
1749566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
1759566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
1769566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
17730602db0SMatthew G. Knepley 
1789566063dSJacob Faibussowitsch   if (ctx->coordSpace) PetscCall(DMCreateCoordinateDisc(*dm));
1798b438e21SMatthew G. Knepley   switch (ctx->meshTransform) {
1808b438e21SMatthew G. Knepley     case TRANSFORM_NONE:
1819566063dSJacob Faibussowitsch       PetscCall(DMPlexRemapGeometry(*dm, 0.0, identity));
1828b438e21SMatthew G. Knepley       break;
1838b438e21SMatthew G. Knepley     case TRANSFORM_SHEAR:
1849566063dSJacob Faibussowitsch       PetscCall(DMPlexShearGeometry(*dm, DM_X, ctx->transformDataReal));
1858b438e21SMatthew G. Knepley       break;
18660c1a66aSMatthew G. Knepley     case TRANSFORM_FLARE:
1879566063dSJacob Faibussowitsch       PetscCall(DMGetCoordinateDM(*dm, &cdm));
1889566063dSJacob Faibussowitsch       PetscCall(DMGetDS(cdm, &cds));
1899566063dSJacob Faibussowitsch       PetscCall(PetscDSSetConstants(cds, 4, ctx->transformData));
1909566063dSJacob Faibussowitsch       PetscCall(DMPlexRemapGeometry(*dm, 0.0, f0_flare));
19160c1a66aSMatthew G. Knepley       break;
1928b438e21SMatthew G. Knepley     case TRANSFORM_ANNULUS:
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_annulus));
1978b438e21SMatthew G. Knepley       break;
1988b438e21SMatthew G. Knepley     case TRANSFORM_SHELL:
1999566063dSJacob Faibussowitsch       PetscCall(DMGetCoordinateDM(*dm, &cdm));
2009566063dSJacob Faibussowitsch       PetscCall(DMGetDS(cdm, &cds));
2019566063dSJacob Faibussowitsch       PetscCall(PetscDSSetConstants(cds, 2, ctx->transformData));
2029566063dSJacob Faibussowitsch       PetscCall(DMPlexRemapGeometry(*dm, 0.0, f0_shell));
2038b438e21SMatthew G. Knepley       break;
20463a3b9bcSJacob Faibussowitsch     default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Unknown mesh transform %d", ctx->meshTransform);
2058b438e21SMatthew G. Knepley   }
2069566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
2078b438e21SMatthew G. Knepley   PetscFunctionReturn(0);
2088b438e21SMatthew G. Knepley }
2098b438e21SMatthew G. Knepley 
2108b438e21SMatthew G. Knepley static void volume(PetscInt dim, PetscInt Nf, PetscInt NfAux,
2118b438e21SMatthew G. Knepley                    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
2128b438e21SMatthew G. Knepley                    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
2138b438e21SMatthew G. Knepley                    PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar vol[])
2148b438e21SMatthew G. Knepley {
2158b438e21SMatthew G. Knepley   vol[0] = 1.;
2168b438e21SMatthew G. Knepley }
2178b438e21SMatthew G. Knepley 
2188b438e21SMatthew G. Knepley static PetscErrorCode CreateDiscretization(DM dm, AppCtx *ctx)
2198b438e21SMatthew G. Knepley {
2208b438e21SMatthew G. Knepley   PetscDS        ds;
2218b438e21SMatthew G. Knepley   PetscFE        fe;
22239290ff6SMatthew G. Knepley   DMPolytopeType ct;
22339290ff6SMatthew G. Knepley   PetscInt       dim, cStart;
22439290ff6SMatthew G. Knepley   PetscBool      simplex;
2258b438e21SMatthew G. Knepley 
2268b438e21SMatthew G. Knepley   PetscFunctionBeginUser;
2279566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
2289566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
2299566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
23039290ff6SMatthew G. Knepley   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
2319566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, NULL, PETSC_DETERMINE, &fe));
2329566063dSJacob Faibussowitsch   PetscCall(PetscFESetName(fe, "scalar"));
2339566063dSJacob Faibussowitsch   PetscCall(DMAddField(dm, NULL, (PetscObject) fe));
2349566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
2359566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
2369566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
2379566063dSJacob Faibussowitsch   PetscCall(PetscDSSetObjective(ds, 0, volume));
2388b438e21SMatthew G. Knepley   PetscFunctionReturn(0);
2398b438e21SMatthew G. Knepley }
2408b438e21SMatthew G. Knepley 
2418b438e21SMatthew G. Knepley static PetscErrorCode CheckVolume(DM dm, AppCtx *ctx)
2428b438e21SMatthew G. Knepley {
2438b438e21SMatthew G. Knepley   Vec            u;
2448b438e21SMatthew G. Knepley   PetscScalar    result;
2458b438e21SMatthew G. Knepley   PetscReal      vol, tol = ctx->tol;
2468b438e21SMatthew G. Knepley 
2478b438e21SMatthew G. Knepley   PetscFunctionBeginUser;
2489566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dm, &u));
2499566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeIntegralFEM(dm, u, &result, ctx));
2508b438e21SMatthew G. Knepley   vol  = PetscRealPart(result);
2519566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dm, &u));
2529566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PetscObjectComm((PetscObject) dm), "Volume: %g\n", (double) vol));
2538b438e21SMatthew G. Knepley   if (ctx->volume > 0.0 && PetscAbsReal(ctx->volume - vol) > tol) {
25498921bdaSJacob 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);
2558b438e21SMatthew G. Knepley   }
2568b438e21SMatthew G. Knepley   PetscFunctionReturn(0);
2578b438e21SMatthew G. Knepley }
2588b438e21SMatthew G. Knepley 
2598b438e21SMatthew G. Knepley int main(int argc, char **argv)
2608b438e21SMatthew G. Knepley {
26139290ff6SMatthew G. Knepley   DM             dm;
2628b438e21SMatthew G. Knepley   AppCtx         user;
2638b438e21SMatthew G. Knepley 
264*327415f7SBarry Smith   PetscFunctionBeginUser;
2659566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL,help));
2669566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
2679566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
2689566063dSJacob Faibussowitsch   PetscCall(CreateDiscretization(dm, &user));
2699566063dSJacob Faibussowitsch   PetscCall(CheckVolume(dm, &user));
2709566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
2719566063dSJacob Faibussowitsch   PetscCall(PetscFree(user.transformDataReal));
2729566063dSJacob Faibussowitsch   PetscCall(PetscFree(user.transformData));
2739566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
274b122ec5aSJacob Faibussowitsch   return 0;
2758b438e21SMatthew G. Knepley }
2768b438e21SMatthew G. Knepley 
2778b438e21SMatthew G. Knepley /*TEST
2788b438e21SMatthew G. Knepley 
27930602db0SMatthew G. Knepley   testset:
28030602db0SMatthew 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
28130602db0SMatthew G. Knepley 
2828b438e21SMatthew G. Knepley     test:
2838b438e21SMatthew G. Knepley       suffix: square_0
28430602db0SMatthew G. Knepley       args: -dm_coord_petscspace_degree 1
2858b438e21SMatthew G. Knepley 
2868b438e21SMatthew G. Knepley     test:
2878b438e21SMatthew G. Knepley       suffix: square_1
28830602db0SMatthew G. Knepley       args: -dm_coord_petscspace_degree 2
2898b438e21SMatthew G. Knepley 
2908b438e21SMatthew G. Knepley     test:
2918b438e21SMatthew G. Knepley       suffix: square_2
29230602db0SMatthew G. Knepley       args: -dm_refine 1 -dm_coord_petscspace_degree 1
2938b438e21SMatthew G. Knepley 
2948b438e21SMatthew G. Knepley     test:
2958b438e21SMatthew G. Knepley       suffix: square_3
29630602db0SMatthew G. Knepley       args: -dm_refine 1 -dm_coord_petscspace_degree 2
29730602db0SMatthew G. Knepley 
29830602db0SMatthew G. Knepley   testset:
29930602db0SMatthew 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
3008b438e21SMatthew G. Knepley 
3018b438e21SMatthew G. Knepley     test:
3028b438e21SMatthew G. Knepley       suffix: cube_0
30330602db0SMatthew G. Knepley       args: -dm_coord_petscspace_degree 1
3048b438e21SMatthew G. Knepley 
3058b438e21SMatthew G. Knepley     test:
3068b438e21SMatthew G. Knepley       suffix: cube_1
30730602db0SMatthew G. Knepley       args: -dm_coord_petscspace_degree 2
3088b438e21SMatthew G. Knepley 
3098b438e21SMatthew G. Knepley     test:
3108b438e21SMatthew G. Knepley       suffix: cube_2
31130602db0SMatthew G. Knepley       args: -dm_refine 1 -dm_coord_petscspace_degree 1
3128b438e21SMatthew G. Knepley 
3138b438e21SMatthew G. Knepley     test:
3148b438e21SMatthew G. Knepley       suffix: cube_3
31530602db0SMatthew G. Knepley       args: -dm_refine 1 -dm_coord_petscspace_degree 2
31630602db0SMatthew G. Knepley 
31730602db0SMatthew G. Knepley   testset:
31830602db0SMatthew 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
3198b438e21SMatthew G. Knepley 
3208b438e21SMatthew G. Knepley     test:
3218b438e21SMatthew G. Knepley       suffix: shear_0
32230602db0SMatthew G. Knepley       args: -dm_coord_petscspace_degree 1 -mesh_transform shear -transform_data 3.0
3238b438e21SMatthew G. Knepley 
3248b438e21SMatthew G. Knepley     test:
3258b438e21SMatthew G. Knepley       suffix: shear_1
32630602db0SMatthew G. Knepley       args: -dm_coord_petscspace_degree 2 -mesh_transform shear -transform_data 3.0
3278b438e21SMatthew G. Knepley 
3288b438e21SMatthew G. Knepley     test:
3298b438e21SMatthew G. Knepley       suffix: shear_2
33030602db0SMatthew G. Knepley       args: -dm_refine 1 -dm_coord_petscspace_degree 1 -mesh_transform shear -transform_data 3.0
3318b438e21SMatthew G. Knepley 
3328b438e21SMatthew G. Knepley     test:
3338b438e21SMatthew G. Knepley       suffix: shear_3
33430602db0SMatthew G. Knepley       args: -dm_refine 1 -dm_coord_petscspace_degree 2 -mesh_transform shear -transform_data 3.0
33530602db0SMatthew G. Knepley 
33630602db0SMatthew G. Knepley   testset:
33730602db0SMatthew 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
3388b438e21SMatthew G. Knepley 
3398b438e21SMatthew G. Knepley     test:
3408b438e21SMatthew G. Knepley       suffix: shear_4
34130602db0SMatthew G. Knepley       args: -dm_coord_petscspace_degree 1 -mesh_transform shear -transform_data 3.0
3428b438e21SMatthew G. Knepley 
3438b438e21SMatthew G. Knepley     test:
3448b438e21SMatthew G. Knepley       suffix: shear_5
34530602db0SMatthew G. Knepley       args: -dm_coord_petscspace_degree 2 -mesh_transform shear -transform_data 3.0
3468b438e21SMatthew G. Knepley 
3478b438e21SMatthew G. Knepley     test:
3488b438e21SMatthew G. Knepley       suffix: shear_6
34930602db0SMatthew G. Knepley       args: -dm_refine 1 -dm_coord_petscspace_degree 1 -mesh_transform shear -transform_data 3.0,4.0
3508b438e21SMatthew G. Knepley 
3518b438e21SMatthew G. Knepley     test:
3528b438e21SMatthew G. Knepley       suffix: shear_7
35330602db0SMatthew G. Knepley       args: -dm_refine 1 -dm_coord_petscspace_degree 2 -mesh_transform shear -transform_data 3.0,4.0
35430602db0SMatthew G. Knepley 
35530602db0SMatthew G. Knepley   testset:
35660c1a66aSMatthew 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. \
35760c1a66aSMatthew G. Knepley           -dm_coord_petscspace_degree 1 -mesh_transform flare -volume 8.
35860c1a66aSMatthew G. Knepley 
35960c1a66aSMatthew G. Knepley     test:
36060c1a66aSMatthew G. Knepley       suffix: flare_0
36160c1a66aSMatthew G. Knepley       args:
36260c1a66aSMatthew G. Knepley 
36360c1a66aSMatthew G. Knepley     test:
36460c1a66aSMatthew G. Knepley       suffix: flare_1
36560c1a66aSMatthew G. Knepley       args: -dm_refine 2
36660c1a66aSMatthew G. Knepley 
36760c1a66aSMatthew G. Knepley   testset:
36830602db0SMatthew G. Knepley     # Area: 3/4 \pi = 2.3562
36930602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -mesh_transform annulus -volume 2.35619449019235 -dm_coord_space 0
3708b438e21SMatthew G. Knepley 
3718b438e21SMatthew G. Knepley     test:
3728b438e21SMatthew G. Knepley       # Area: (a+b)/2 h = 3/\sqrt{2} (sqrt{2} - 1/\sqrt{2}) = 3/2
3738b438e21SMatthew G. Knepley       suffix: annulus_0
3748b438e21SMatthew G. Knepley       requires: double
37530602db0SMatthew G. Knepley       args: -dm_coord_petscspace_degree 1 -volume 1.5
3768b438e21SMatthew G. Knepley 
3778b438e21SMatthew G. Knepley     test:
3788b438e21SMatthew G. Knepley       suffix: annulus_1
3798b438e21SMatthew G. Knepley       requires: double
38030602db0SMatthew G. Knepley       args: -dm_refine 3 -dm_coord_petscspace_degree 1 -tol .016
3818b438e21SMatthew G. Knepley 
3828b438e21SMatthew G. Knepley     test:
3838b438e21SMatthew G. Knepley       suffix: annulus_2
3848b438e21SMatthew G. Knepley       requires: double
38530602db0SMatthew G. Knepley       args: -dm_refine 3 -dm_coord_petscspace_degree 2 -tol .0038
3868b438e21SMatthew G. Knepley 
3878b438e21SMatthew G. Knepley     test:
3888b438e21SMatthew G. Knepley       suffix: annulus_3
3898b438e21SMatthew G. Knepley       requires: double
39030602db0SMatthew G. Knepley       args: -dm_refine 3 -dm_coord_petscspace_degree 3 -tol 2.2e-6
3918b438e21SMatthew G. Knepley 
3928b438e21SMatthew G. Knepley     test:
3938b438e21SMatthew G. Knepley       suffix: annulus_4
3948b438e21SMatthew G. Knepley       requires: double
39530602db0SMatthew G. Knepley       args: -dm_refine 2 -dm_coord_petscspace_degree 2 -petscfe_default_quadrature_order 2 -tol .00012
3968b438e21SMatthew G. Knepley 
3978b438e21SMatthew G. Knepley     test:
3988b438e21SMatthew G. Knepley       suffix: annulus_5
3998b438e21SMatthew G. Knepley       requires: double
40030602db0SMatthew G. Knepley       args: -dm_refine 2 -dm_coord_petscspace_degree 3 -petscfe_default_quadrature_order 3 -tol 1.2e-7
40130602db0SMatthew G. Knepley 
40230602db0SMatthew G. Knepley   testset:
40330602db0SMatthew G. Knepley     # Volume: 4/3 \pi (8 - 1)/2 = 14/3 \pi = 14.66076571675238
40430602db0SMatthew 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
4058b438e21SMatthew G. Knepley 
4068b438e21SMatthew G. Knepley     test:
4078b438e21SMatthew G. Knepley       suffix: shell_0
4088b438e21SMatthew G. Knepley       requires: double
40930602db0SMatthew G. Knepley       args: -dm_refine 1 -dm_coord_petscspace_degree 1 -petscfe_default_quadrature_order 1 -volume 5.633164922 -tol 1.0e-7
4108b438e21SMatthew G. Knepley 
4118b438e21SMatthew G. Knepley     test:
4128b438e21SMatthew G. Knepley       suffix: shell_1
4138b438e21SMatthew G. Knepley       requires: double
41430602db0SMatthew G. Knepley       args: -dm_refine 2 -dm_coord_petscspace_degree 1 -petscfe_default_quadrature_order 1 -tol 3.1
4158b438e21SMatthew G. Knepley 
4168b438e21SMatthew G. Knepley     test:
4178b438e21SMatthew G. Knepley       suffix: shell_2
4188b438e21SMatthew G. Knepley       requires: double
41930602db0SMatthew G. Knepley       args: -dm_refine 2 -dm_coord_petscspace_degree 2 -petscfe_default_quadrature_order 2 -tol .1
4208b438e21SMatthew G. Knepley 
4218b438e21SMatthew G. Knepley     test:
4228b438e21SMatthew G. Knepley       suffix: shell_3
4238b438e21SMatthew G. Knepley       requires: double
42430602db0SMatthew G. Knepley       args: -dm_refine 2 -dm_coord_petscspace_degree 3 -petscfe_default_quadrature_order 3 -tol .02
4258b438e21SMatthew G. Knepley 
4268b438e21SMatthew G. Knepley     test:
4278b438e21SMatthew G. Knepley       suffix: shell_4
4288b438e21SMatthew G. Knepley       requires: double
42930602db0SMatthew G. Knepley       args: -dm_refine 2 -dm_coord_petscspace_degree 4 -petscfe_default_quadrature_order 4 -tol .006
43039290ff6SMatthew G. Knepley 
43139290ff6SMatthew G. Knepley   test:
43239290ff6SMatthew G. Knepley     # Volume: 1.0
43339290ff6SMatthew G. Knepley     suffix: gmsh_q2
43439290ff6SMatthew G. Knepley     requires: double
43530602db0SMatthew 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
43639290ff6SMatthew G. Knepley 
43739290ff6SMatthew G. Knepley   test:
43839290ff6SMatthew G. Knepley     # Volume: 1.0
43939290ff6SMatthew G. Knepley     suffix: gmsh_q3
44039290ff6SMatthew G. Knepley     requires: double
44124b9a4b1SJed Brown     nsize: {{1 2}}
44230602db0SMatthew 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
4438b438e21SMatthew G. Knepley 
4448b438e21SMatthew G. Knepley TEST*/
445