xref: /petsc/src/dm/impls/plex/tests/ex33.c (revision 30602db00db74b7e41a0c75e517aefe6711423f0)
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 
68b438e21SMatthew G. Knepley typedef enum {TRANSFORM_NONE, TRANSFORM_SHEAR, TRANSFORM_ANNULUS, TRANSFORM_SHELL} Transform;
78b438e21SMatthew G. Knepley const char * const TransformTypes[] = {"none", "shear", "annulus", "shell", "Mesh Transform", "TRANSFORM_", NULL};
88b438e21SMatthew G. Knepley 
98b438e21SMatthew G. Knepley typedef struct {
10*30602db0SMatthew 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   PetscErrorCode ierr;
228b438e21SMatthew G. Knepley 
238b438e21SMatthew G. Knepley   PetscFunctionBegin;
24*30602db0SMatthew G. Knepley   options->coordSpace        = PETSC_TRUE;
258b438e21SMatthew G. Knepley   options->meshTransform     = TRANSFORM_NONE;
268b438e21SMatthew G. Knepley   options->transformDataReal = NULL;
278b438e21SMatthew G. Knepley   options->transformData     = NULL;
288b438e21SMatthew G. Knepley   options->volume            = -1.0;
298b438e21SMatthew G. Knepley   options->tol               = PETSC_SMALL;
308b438e21SMatthew G. Knepley 
318b438e21SMatthew G. Knepley   ierr = PetscOptionsBegin(comm, "", "Meshing Interpolation Test Options", "DMPLEX");CHKERRQ(ierr);
32*30602db0SMatthew G. Knepley   ierr = PetscOptionsBool("-coord_space", "Flag to create a coordinate space", "ex33.c", options->coordSpace, &options->coordSpace, NULL);CHKERRQ(ierr);
338b438e21SMatthew G. Knepley   ierr = PetscOptionsEnum("-mesh_transform", "Method to transform initial box mesh <none, shear, annulus, shell>", "ex33.c", TransformTypes, (PetscEnum) options->meshTransform, (PetscEnum *) &options->meshTransform, NULL);CHKERRQ(ierr);
348b438e21SMatthew G. Knepley   switch (options->meshTransform) {
358b438e21SMatthew G. Knepley     case TRANSFORM_NONE: break;
368b438e21SMatthew G. Knepley     case TRANSFORM_SHEAR:
3739290ff6SMatthew G. Knepley       n = 2;
388b438e21SMatthew G. Knepley       ierr = PetscMalloc1(n, &options->transformDataReal);CHKERRQ(ierr);
398b438e21SMatthew G. Knepley       for (i = 0; i < n; ++i) options->transformDataReal[i] = 1.0;
408b438e21SMatthew G. Knepley       ierr = PetscOptionsRealArray("-transform_data", "Parameters for mesh transforms", "ex33.c", options->transformDataReal, &n, NULL);CHKERRQ(ierr);
418b438e21SMatthew G. Knepley       break;
428b438e21SMatthew G. Knepley     case TRANSFORM_ANNULUS:
438b438e21SMatthew G. Knepley       n = 2;
448b438e21SMatthew G. Knepley       ierr = PetscMalloc1(n, &options->transformData);CHKERRQ(ierr);
458b438e21SMatthew G. Knepley       options->transformData[0] = 1.0;
468b438e21SMatthew G. Knepley       options->transformData[1] = 2.0;
478b438e21SMatthew G. Knepley       ierr = PetscOptionsScalarArray("-transform_data", "Parameters for mesh transforms", "ex33.c", options->transformData, &n, NULL);CHKERRQ(ierr);
488b438e21SMatthew G. Knepley       break;
498b438e21SMatthew G. Knepley     case TRANSFORM_SHELL:
508b438e21SMatthew G. Knepley       n = 2;
518b438e21SMatthew G. Knepley       ierr = PetscMalloc1(n, &options->transformData);CHKERRQ(ierr);
528b438e21SMatthew G. Knepley       options->transformData[0] = 1.0;
538b438e21SMatthew G. Knepley       options->transformData[1] = 2.0;
548b438e21SMatthew G. Knepley       ierr = PetscOptionsScalarArray("-transform_data", "Parameters for mesh transforms", "ex33.c", options->transformData, &n, NULL);CHKERRQ(ierr);
558b438e21SMatthew G. Knepley       break;
568b438e21SMatthew G. Knepley     default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Unknown mesh transform %D", options->meshTransform);
578b438e21SMatthew G. Knepley   }
588b438e21SMatthew G. Knepley   ierr = PetscOptionsReal("-volume", "The analytical volume of the mesh", "ex33.c", options->volume, &options->volume, NULL);CHKERRQ(ierr);
598b438e21SMatthew G. Knepley   ierr = PetscOptionsReal("-tol", "The tolerance for the volume check", "ex33.c", options->tol, &options->tol, NULL);CHKERRQ(ierr);
608b438e21SMatthew G. Knepley   ierr = PetscOptionsEnd();
618b438e21SMatthew G. Knepley   PetscFunctionReturn(0);
628b438e21SMatthew G. Knepley }
638b438e21SMatthew G. Knepley 
648b438e21SMatthew G. Knepley static void identity(PetscInt dim, PetscInt Nf, PetscInt NfAux,
658b438e21SMatthew G. Knepley                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
668b438e21SMatthew G. Knepley                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
678b438e21SMatthew G. Knepley                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
688b438e21SMatthew G. Knepley {
698b438e21SMatthew G. Knepley   const PetscInt Nc = uOff[1] - uOff[0];
708b438e21SMatthew G. Knepley   PetscInt       c;
718b438e21SMatthew G. Knepley 
728b438e21SMatthew G. Knepley   for (c = 0; c < Nc; ++c) f0[c] = u[c];
738b438e21SMatthew G. Knepley }
748b438e21SMatthew G. Knepley 
758b438e21SMatthew G. Knepley /*
768b438e21SMatthew 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
778b438e21SMatthew G. Knepley   will correspond to the top and bottom of our square. So
788b438e21SMatthew G. Knepley 
798b438e21SMatthew G. Knepley     (0,0)--(1,0)  ==>  (1,0)--(2,0)      Just a shift of (1,0)
808b438e21SMatthew G. Knepley     (0,1)--(1,1)  ==>  (0,1)--(0,2)      Switch x and y
818b438e21SMatthew G. Knepley 
828b438e21SMatthew 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:
838b438e21SMatthew G. Knepley 
848b438e21SMatthew G. Knepley     (x, y)  ==>  (x+1, \pi/2 y)                           in (r', \theta') space
858b438e21SMatthew G. Knepley             ==>  ((x+1) cos(\pi/2 y), (x+1) sin(\pi/2 y)) in (x', y') space
868b438e21SMatthew G. Knepley */
878b438e21SMatthew G. Knepley static void f0_annulus(PetscInt dim, PetscInt Nf, PetscInt NfAux,
888b438e21SMatthew G. Knepley                        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
898b438e21SMatthew G. Knepley                        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
908b438e21SMatthew G. Knepley                        PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar xp[])
918b438e21SMatthew G. Knepley {
928b438e21SMatthew G. Knepley   const PetscReal ri = PetscRealPart(constants[0]);
938b438e21SMatthew G. Knepley   const PetscReal ro = PetscRealPart(constants[1]);
948b438e21SMatthew G. Knepley 
958b438e21SMatthew G. Knepley   xp[0] = (x[0] * (ro-ri) + ri) * PetscCosReal(0.5*PETSC_PI*x[1]);
968b438e21SMatthew G. Knepley   xp[1] = (x[0] * (ro-ri) + ri) * PetscSinReal(0.5*PETSC_PI*x[1]);
978b438e21SMatthew G. Knepley }
988b438e21SMatthew G. Knepley 
998b438e21SMatthew G. Knepley /*
1008b438e21SMatthew 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
1018b438e21SMatthew G. Knepley   lower hemisphere and the upper surface onto the top, letting z be the radius.
1028b438e21SMatthew G. Knepley 
1038b438e21SMatthew G. Knepley     (x, y)  ==>  ((z+3)/2, \pi/2 (|x| or |y|), arctan y/x)                                                  in (r', \theta', \phi') space
1048b438e21SMatthew 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
1058b438e21SMatthew G. Knepley */
1068b438e21SMatthew G. Knepley static void f0_shell(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1078b438e21SMatthew G. Knepley                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1088b438e21SMatthew G. Knepley                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1098b438e21SMatthew G. Knepley                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar xp[])
1108b438e21SMatthew G. Knepley {
1118b438e21SMatthew G. Knepley   const PetscReal pi4    = PETSC_PI/4.0;
1128b438e21SMatthew G. Knepley   const PetscReal ri     = PetscRealPart(constants[0]);
1138b438e21SMatthew G. Knepley   const PetscReal ro     = PetscRealPart(constants[1]);
1148b438e21SMatthew G. Knepley   const PetscReal rp     = (x[2]+1) * 0.5*(ro-ri) + ri;
1158b438e21SMatthew G. Knepley   const PetscReal phip   = PetscAtan2Real(x[1], x[0]);
1168b438e21SMatthew 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])));
1178b438e21SMatthew G. Knepley 
1188b438e21SMatthew G. Knepley   xp[0] = rp * PetscCosReal(thetap) * PetscCosReal(phip);
1198b438e21SMatthew G. Knepley   xp[1] = rp * PetscCosReal(thetap) * PetscSinReal(phip);
1208b438e21SMatthew G. Knepley   xp[2] = rp * PetscSinReal(thetap);
1218b438e21SMatthew G. Knepley }
1228b438e21SMatthew G. Knepley 
12339290ff6SMatthew G. Knepley static PetscErrorCode DMCreateCoordinateDisc(DM dm)
1248b438e21SMatthew G. Knepley {
1258b438e21SMatthew G. Knepley   DM             cdm;
1268b438e21SMatthew G. Knepley   PetscFE        fe;
12739290ff6SMatthew G. Knepley   DMPolytopeType ct;
12839290ff6SMatthew G. Knepley   PetscInt       dim, dE, cStart;
12939290ff6SMatthew G. Knepley   PetscBool      simplex;
1308b438e21SMatthew G. Knepley   PetscErrorCode ierr;
1318b438e21SMatthew G. Knepley 
1328b438e21SMatthew G. Knepley   PetscFunctionBegin;
1338b438e21SMatthew G. Knepley   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
1348b438e21SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1358b438e21SMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dE);CHKERRQ(ierr);
13639290ff6SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(cdm, 0, &cStart, NULL);CHKERRQ(ierr);
13739290ff6SMatthew G. Knepley   ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr);
13839290ff6SMatthew G. Knepley   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
139*30602db0SMatthew G. Knepley   ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, dE, simplex, "dm_coord_", -1, &fe);CHKERRQ(ierr);
14039290ff6SMatthew G. Knepley   ierr = DMProjectCoordinates(dm, fe);CHKERRQ(ierr);
1418b438e21SMatthew G. Knepley   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
1428b438e21SMatthew G. Knepley   PetscFunctionReturn(0);
1438b438e21SMatthew G. Knepley }
1448b438e21SMatthew G. Knepley 
1458b438e21SMatthew G. Knepley PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm)
1468b438e21SMatthew G. Knepley {
1478b438e21SMatthew G. Knepley   DM             cdm;
1488b438e21SMatthew G. Knepley   PetscDS        cds;
1498b438e21SMatthew G. Knepley   PetscErrorCode ierr;
1508b438e21SMatthew G. Knepley 
1518b438e21SMatthew G. Knepley   PetscFunctionBegin;
152*30602db0SMatthew G. Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
153*30602db0SMatthew G. Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
1548b438e21SMatthew G. Knepley   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
155*30602db0SMatthew G. Knepley 
156*30602db0SMatthew G. Knepley   if (ctx->coordSpace) {ierr = DMCreateCoordinateDisc(*dm);CHKERRQ(ierr);}
1578b438e21SMatthew G. Knepley   switch (ctx->meshTransform) {
1588b438e21SMatthew G. Knepley     case TRANSFORM_NONE:
1598b438e21SMatthew G. Knepley       ierr = DMPlexRemapGeometry(*dm, 0.0, identity);CHKERRQ(ierr);
1608b438e21SMatthew G. Knepley       break;
1618b438e21SMatthew G. Knepley     case TRANSFORM_SHEAR:
1628b438e21SMatthew G. Knepley       ierr = DMPlexShearGeometry(*dm, DM_X, ctx->transformDataReal);CHKERRQ(ierr);
1638b438e21SMatthew G. Knepley       break;
1648b438e21SMatthew G. Knepley     case TRANSFORM_ANNULUS:
1658b438e21SMatthew G. Knepley       ierr = DMGetCoordinateDM(*dm, &cdm);CHKERRQ(ierr);
1668b438e21SMatthew G. Knepley       ierr = DMGetDS(cdm, &cds);CHKERRQ(ierr);
1678b438e21SMatthew G. Knepley       ierr = PetscDSSetConstants(cds, 2, ctx->transformData);CHKERRQ(ierr);
1688b438e21SMatthew G. Knepley       ierr = DMPlexRemapGeometry(*dm, 0.0, f0_annulus);CHKERRQ(ierr);
1698b438e21SMatthew G. Knepley       break;
1708b438e21SMatthew G. Knepley     case TRANSFORM_SHELL:
1718b438e21SMatthew G. Knepley       ierr = DMGetCoordinateDM(*dm, &cdm);CHKERRQ(ierr);
1728b438e21SMatthew G. Knepley       ierr = DMGetDS(cdm, &cds);CHKERRQ(ierr);
1738b438e21SMatthew G. Knepley       ierr = PetscDSSetConstants(cds, 2, ctx->transformData);CHKERRQ(ierr);
1748b438e21SMatthew G. Knepley       ierr = DMPlexRemapGeometry(*dm, 0.0, f0_shell);CHKERRQ(ierr);
1758b438e21SMatthew G. Knepley       break;
1768b438e21SMatthew G. Knepley     default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Unknown mesh transform %D", ctx->meshTransform);
1778b438e21SMatthew G. Knepley   }
1788b438e21SMatthew G. Knepley   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
1798b438e21SMatthew G. Knepley   PetscFunctionReturn(0);
1808b438e21SMatthew G. Knepley }
1818b438e21SMatthew G. Knepley 
1828b438e21SMatthew G. Knepley static void volume(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1838b438e21SMatthew G. Knepley                    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1848b438e21SMatthew G. Knepley                    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1858b438e21SMatthew G. Knepley                    PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar vol[])
1868b438e21SMatthew G. Knepley {
1878b438e21SMatthew G. Knepley   vol[0] = 1.;
1888b438e21SMatthew G. Knepley }
1898b438e21SMatthew G. Knepley 
1908b438e21SMatthew G. Knepley static PetscErrorCode CreateDiscretization(DM dm, AppCtx *ctx)
1918b438e21SMatthew G. Knepley {
1928b438e21SMatthew G. Knepley   PetscDS        ds;
1938b438e21SMatthew G. Knepley   PetscFE        fe;
19439290ff6SMatthew G. Knepley   DMPolytopeType ct;
19539290ff6SMatthew G. Knepley   PetscInt       dim, cStart;
19639290ff6SMatthew G. Knepley   PetscBool      simplex;
1978b438e21SMatthew G. Knepley   PetscErrorCode ierr;
1988b438e21SMatthew G. Knepley 
1998b438e21SMatthew G. Knepley   PetscFunctionBeginUser;
20039290ff6SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
20139290ff6SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr);
20239290ff6SMatthew G. Knepley   ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr);
20339290ff6SMatthew G. Knepley   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
20439290ff6SMatthew G. Knepley   ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, NULL, PETSC_DETERMINE, &fe);CHKERRQ(ierr);
2058b438e21SMatthew G. Knepley   ierr = PetscFESetName(fe, "scalar");CHKERRQ(ierr);
2068b438e21SMatthew G. Knepley   ierr = DMAddField(dm, NULL, (PetscObject) fe);
2078b438e21SMatthew G. Knepley   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
2088b438e21SMatthew G. Knepley   ierr = DMCreateDS(dm);CHKERRQ(ierr);
2098b438e21SMatthew G. Knepley   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
2108b438e21SMatthew G. Knepley   ierr = PetscDSSetObjective(ds, 0, volume);CHKERRQ(ierr);
2118b438e21SMatthew G. Knepley   PetscFunctionReturn(0);
2128b438e21SMatthew G. Knepley }
2138b438e21SMatthew G. Knepley 
2148b438e21SMatthew G. Knepley static PetscErrorCode CheckVolume(DM dm, AppCtx *ctx)
2158b438e21SMatthew G. Knepley {
2168b438e21SMatthew G. Knepley   Vec            u;
2178b438e21SMatthew G. Knepley   PetscScalar    result;
2188b438e21SMatthew G. Knepley   PetscReal      vol, tol = ctx->tol;
2198b438e21SMatthew G. Knepley   PetscErrorCode ierr;
2208b438e21SMatthew G. Knepley 
2218b438e21SMatthew G. Knepley   PetscFunctionBeginUser;
2228b438e21SMatthew G. Knepley   ierr = DMGetGlobalVector(dm, &u);CHKERRQ(ierr);
2238b438e21SMatthew G. Knepley   ierr = DMPlexComputeIntegralFEM(dm, u, &result, ctx);CHKERRQ(ierr);
2248b438e21SMatthew G. Knepley   vol  = PetscRealPart(result);
2258b438e21SMatthew G. Knepley   ierr = DMRestoreGlobalVector(dm, &u);CHKERRQ(ierr);
2268b438e21SMatthew G. Knepley   ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Volume: %g\n", (double) vol);CHKERRQ(ierr);
2278b438e21SMatthew G. Knepley   if (ctx->volume > 0.0 && PetscAbsReal(ctx->volume - vol) > tol) {
2288b438e21SMatthew G. Knepley     SETERRQ4(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);
2298b438e21SMatthew G. Knepley   }
2308b438e21SMatthew G. Knepley   PetscFunctionReturn(0);
2318b438e21SMatthew G. Knepley }
2328b438e21SMatthew G. Knepley 
2338b438e21SMatthew G. Knepley int main(int argc, char **argv)
2348b438e21SMatthew G. Knepley {
23539290ff6SMatthew G. Knepley   DM             dm;
2368b438e21SMatthew G. Knepley   AppCtx         user;
2378b438e21SMatthew G. Knepley   PetscErrorCode ierr;
2388b438e21SMatthew G. Knepley 
2398b438e21SMatthew G. Knepley   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
2408b438e21SMatthew G. Knepley   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
24139290ff6SMatthew G. Knepley   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
24239290ff6SMatthew G. Knepley   ierr = CreateDiscretization(dm, &user);CHKERRQ(ierr);
24339290ff6SMatthew G. Knepley   ierr = CheckVolume(dm, &user);CHKERRQ(ierr);
24439290ff6SMatthew G. Knepley   ierr = DMDestroy(&dm);CHKERRQ(ierr);
2458b438e21SMatthew G. Knepley   ierr = PetscFree(user.transformDataReal);CHKERRQ(ierr);
2468b438e21SMatthew G. Knepley   ierr = PetscFree(user.transformData);CHKERRQ(ierr);
2478b438e21SMatthew G. Knepley   ierr = PetscFinalize();
2488b438e21SMatthew G. Knepley   return ierr;
2498b438e21SMatthew G. Knepley }
2508b438e21SMatthew G. Knepley 
2518b438e21SMatthew G. Knepley /*TEST
2528b438e21SMatthew G. Knepley 
253*30602db0SMatthew G. Knepley   testset:
254*30602db0SMatthew 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
255*30602db0SMatthew G. Knepley 
2568b438e21SMatthew G. Knepley     test:
2578b438e21SMatthew G. Knepley       suffix: square_0
258*30602db0SMatthew G. Knepley       args: -dm_coord_petscspace_degree 1
2598b438e21SMatthew G. Knepley 
2608b438e21SMatthew G. Knepley     test:
2618b438e21SMatthew G. Knepley       suffix: square_1
262*30602db0SMatthew G. Knepley       args: -dm_coord_petscspace_degree 2
2638b438e21SMatthew G. Knepley 
2648b438e21SMatthew G. Knepley     test:
2658b438e21SMatthew G. Knepley       suffix: square_2
266*30602db0SMatthew G. Knepley       args: -dm_refine 1 -dm_coord_petscspace_degree 1
2678b438e21SMatthew G. Knepley 
2688b438e21SMatthew G. Knepley     test:
2698b438e21SMatthew G. Knepley       suffix: square_3
270*30602db0SMatthew G. Knepley       args: -dm_refine 1 -dm_coord_petscspace_degree 2
271*30602db0SMatthew G. Knepley 
272*30602db0SMatthew G. Knepley   testset:
273*30602db0SMatthew 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
277*30602db0SMatthew G. Knepley       args: -dm_coord_petscspace_degree 1
2788b438e21SMatthew G. Knepley 
2798b438e21SMatthew G. Knepley     test:
2808b438e21SMatthew G. Knepley       suffix: cube_1
281*30602db0SMatthew G. Knepley       args: -dm_coord_petscspace_degree 2
2828b438e21SMatthew G. Knepley 
2838b438e21SMatthew G. Knepley     test:
2848b438e21SMatthew G. Knepley       suffix: cube_2
285*30602db0SMatthew G. Knepley       args: -dm_refine 1 -dm_coord_petscspace_degree 1
2868b438e21SMatthew G. Knepley 
2878b438e21SMatthew G. Knepley     test:
2888b438e21SMatthew G. Knepley       suffix: cube_3
289*30602db0SMatthew G. Knepley       args: -dm_refine 1 -dm_coord_petscspace_degree 2
290*30602db0SMatthew G. Knepley 
291*30602db0SMatthew G. Knepley   testset:
292*30602db0SMatthew 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
296*30602db0SMatthew 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
300*30602db0SMatthew 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
304*30602db0SMatthew 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
308*30602db0SMatthew G. Knepley       args: -dm_refine 1 -dm_coord_petscspace_degree 2 -mesh_transform shear -transform_data 3.0
309*30602db0SMatthew G. Knepley 
310*30602db0SMatthew G. Knepley   testset:
311*30602db0SMatthew 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
315*30602db0SMatthew 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
319*30602db0SMatthew 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
323*30602db0SMatthew 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
327*30602db0SMatthew G. Knepley       args: -dm_refine 1 -dm_coord_petscspace_degree 2 -mesh_transform shear -transform_data 3.0,4.0
328*30602db0SMatthew G. Knepley 
329*30602db0SMatthew G. Knepley   testset:
330*30602db0SMatthew G. Knepley     # Area: 3/4 \pi = 2.3562
331*30602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -mesh_transform annulus -volume 2.35619449019235 -dm_coord_space 0
3328b438e21SMatthew G. Knepley 
3338b438e21SMatthew G. Knepley     test:
3348b438e21SMatthew G. Knepley       # Area: (a+b)/2 h = 3/\sqrt{2} (sqrt{2} - 1/\sqrt{2}) = 3/2
3358b438e21SMatthew G. Knepley       suffix: annulus_0
3368b438e21SMatthew G. Knepley       requires: double
337*30602db0SMatthew G. Knepley       args: -dm_coord_petscspace_degree 1 -volume 1.5
3388b438e21SMatthew G. Knepley 
3398b438e21SMatthew G. Knepley     test:
3408b438e21SMatthew G. Knepley       suffix: annulus_1
3418b438e21SMatthew G. Knepley       requires: double
342*30602db0SMatthew G. Knepley       args: -dm_refine 3 -dm_coord_petscspace_degree 1 -tol .016
3438b438e21SMatthew G. Knepley 
3448b438e21SMatthew G. Knepley     test:
3458b438e21SMatthew G. Knepley       suffix: annulus_2
3468b438e21SMatthew G. Knepley       requires: double
347*30602db0SMatthew G. Knepley       args: -dm_refine 3 -dm_coord_petscspace_degree 2 -tol .0038
3488b438e21SMatthew G. Knepley 
3498b438e21SMatthew G. Knepley     test:
3508b438e21SMatthew G. Knepley       suffix: annulus_3
3518b438e21SMatthew G. Knepley       requires: double
352*30602db0SMatthew G. Knepley       args: -dm_refine 3 -dm_coord_petscspace_degree 3 -tol 2.2e-6
3538b438e21SMatthew G. Knepley 
3548b438e21SMatthew G. Knepley     test:
3558b438e21SMatthew G. Knepley       suffix: annulus_4
3568b438e21SMatthew G. Knepley       requires: double
357*30602db0SMatthew G. Knepley       args: -dm_refine 2 -dm_coord_petscspace_degree 2 -petscfe_default_quadrature_order 2 -tol .00012
3588b438e21SMatthew G. Knepley 
3598b438e21SMatthew G. Knepley     test:
3608b438e21SMatthew G. Knepley       suffix: annulus_5
3618b438e21SMatthew G. Knepley       requires: double
362*30602db0SMatthew G. Knepley       args: -dm_refine 2 -dm_coord_petscspace_degree 3 -petscfe_default_quadrature_order 3 -tol 1.2e-7
363*30602db0SMatthew G. Knepley 
364*30602db0SMatthew G. Knepley   testset:
365*30602db0SMatthew G. Knepley     # Volume: 4/3 \pi (8 - 1)/2 = 14/3 \pi = 14.66076571675238
366*30602db0SMatthew 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
3678b438e21SMatthew G. Knepley 
3688b438e21SMatthew G. Knepley     test:
3698b438e21SMatthew G. Knepley       suffix: shell_0
3708b438e21SMatthew G. Knepley       requires: double
371*30602db0SMatthew G. Knepley       args: -dm_refine 1 -dm_coord_petscspace_degree 1 -petscfe_default_quadrature_order 1 -volume 5.633164922 -tol 1.0e-7
3728b438e21SMatthew G. Knepley 
3738b438e21SMatthew G. Knepley     test:
3748b438e21SMatthew G. Knepley       suffix: shell_1
3758b438e21SMatthew G. Knepley       requires: double
376*30602db0SMatthew G. Knepley       args: -dm_refine 2 -dm_coord_petscspace_degree 1 -petscfe_default_quadrature_order 1 -tol 3.1
3778b438e21SMatthew G. Knepley 
3788b438e21SMatthew G. Knepley     test:
3798b438e21SMatthew G. Knepley       suffix: shell_2
3808b438e21SMatthew G. Knepley       requires: double
381*30602db0SMatthew G. Knepley       args: -dm_refine 2 -dm_coord_petscspace_degree 2 -petscfe_default_quadrature_order 2 -tol .1
3828b438e21SMatthew G. Knepley 
3838b438e21SMatthew G. Knepley     test:
3848b438e21SMatthew G. Knepley       suffix: shell_3
3858b438e21SMatthew G. Knepley       requires: double
386*30602db0SMatthew G. Knepley       args: -dm_refine 2 -dm_coord_petscspace_degree 3 -petscfe_default_quadrature_order 3 -tol .02
3878b438e21SMatthew G. Knepley 
3888b438e21SMatthew G. Knepley     test:
3898b438e21SMatthew G. Knepley       suffix: shell_4
3908b438e21SMatthew G. Knepley       requires: double
391*30602db0SMatthew G. Knepley       args: -dm_refine 2 -dm_coord_petscspace_degree 4 -petscfe_default_quadrature_order 4 -tol .006
39239290ff6SMatthew G. Knepley 
39339290ff6SMatthew G. Knepley   test:
39439290ff6SMatthew G. Knepley     # Volume: 1.0
39539290ff6SMatthew G. Knepley     suffix: gmsh_q2
39639290ff6SMatthew G. Knepley     requires: double
397*30602db0SMatthew 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
39839290ff6SMatthew G. Knepley 
39939290ff6SMatthew G. Knepley   test:
40039290ff6SMatthew G. Knepley     # Volume: 1.0
40139290ff6SMatthew G. Knepley     suffix: gmsh_q3
40239290ff6SMatthew G. Knepley     requires: double
403*30602db0SMatthew 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
4048b438e21SMatthew G. Knepley 
4058b438e21SMatthew G. Knepley TEST*/
406