xref: /petsc/src/dm/impls/plex/tests/ex33.c (revision 8b438e214ee71501125a9982e700ba9ee276a16d)
1*8b438e21SMatthew G. Knepley static char help[] = "Tests for high order geometry\n\n";
2*8b438e21SMatthew G. Knepley 
3*8b438e21SMatthew G. Knepley #include <petscdmplex.h>
4*8b438e21SMatthew G. Knepley #include <petscds.h>
5*8b438e21SMatthew G. Knepley 
6*8b438e21SMatthew G. Knepley typedef enum {TRANSFORM_NONE, TRANSFORM_SHEAR, TRANSFORM_ANNULUS, TRANSFORM_SHELL} Transform;
7*8b438e21SMatthew G. Knepley const char * const TransformTypes[] = {"none", "shear", "annulus", "shell", "Mesh Transform", "TRANSFORM_", NULL};
8*8b438e21SMatthew G. Knepley 
9*8b438e21SMatthew G. Knepley typedef struct {
10*8b438e21SMatthew G. Knepley   DM          dm;
11*8b438e21SMatthew G. Knepley   PetscInt    dim;              /* The topological mesh dimension */
12*8b438e21SMatthew G. Knepley   PetscBool   simplex;          /* Use simplices or hexes */
13*8b438e21SMatthew G. Knepley   char        filename[PETSC_MAX_PATH_LEN]; /* Import mesh from file */
14*8b438e21SMatthew G. Knepley   Transform   meshTransform;    /* Transform for initial box mesh */
15*8b438e21SMatthew G. Knepley   PetscReal   *transformDataReal; /* Parameters for mesh transform */
16*8b438e21SMatthew G. Knepley   PetscScalar *transformData;   /* Parameters for mesh transform */
17*8b438e21SMatthew G. Knepley   PetscReal   volume;           /* Analytical volume of the mesh */
18*8b438e21SMatthew G. Knepley   PetscReal   tol;              /* Tolerance for volume check */
19*8b438e21SMatthew G. Knepley } AppCtx;
20*8b438e21SMatthew G. Knepley 
21*8b438e21SMatthew G. Knepley PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
22*8b438e21SMatthew G. Knepley {
23*8b438e21SMatthew G. Knepley   PetscInt       n = 0, i;
24*8b438e21SMatthew G. Knepley   PetscErrorCode ierr;
25*8b438e21SMatthew G. Knepley 
26*8b438e21SMatthew G. Knepley   PetscFunctionBegin;
27*8b438e21SMatthew G. Knepley   options->dim                = 2;
28*8b438e21SMatthew G. Knepley   options->simplex            = PETSC_TRUE;
29*8b438e21SMatthew G. Knepley   options->filename[0]        = '\0';
30*8b438e21SMatthew G. Knepley   options->meshTransform      = TRANSFORM_NONE;
31*8b438e21SMatthew G. Knepley   options->transformDataReal = NULL;
32*8b438e21SMatthew G. Knepley   options->transformData      = NULL;
33*8b438e21SMatthew G. Knepley   options->volume             = -1.0;
34*8b438e21SMatthew G. Knepley   options->tol                = PETSC_SMALL;
35*8b438e21SMatthew G. Knepley 
36*8b438e21SMatthew G. Knepley   ierr = PetscOptionsBegin(comm, "", "Meshing Interpolation Test Options", "DMPLEX");CHKERRQ(ierr);
37*8b438e21SMatthew G. Knepley   ierr = PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex33.c", options->dim, &options->dim, NULL,1,3);CHKERRQ(ierr);
38*8b438e21SMatthew G. Knepley   ierr = PetscOptionsBool("-simplex", "Use simplices if true, otherwise hexes", "ex33.c", options->simplex, &options->simplex, NULL);CHKERRQ(ierr);
39*8b438e21SMatthew G. Knepley   ierr = PetscOptionsString("-filename", "The mesh file", "ex33.c", options->filename, options->filename, sizeof(options->filename), NULL);CHKERRQ(ierr);
40*8b438e21SMatthew 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);
41*8b438e21SMatthew G. Knepley   switch (options->meshTransform) {
42*8b438e21SMatthew G. Knepley     case TRANSFORM_NONE: break;
43*8b438e21SMatthew G. Knepley     case TRANSFORM_SHEAR:
44*8b438e21SMatthew G. Knepley       n = options->dim-1;
45*8b438e21SMatthew G. Knepley       ierr = PetscMalloc1(n, &options->transformDataReal);CHKERRQ(ierr);
46*8b438e21SMatthew G. Knepley       for (i = 0; i < n; ++i) options->transformDataReal[i] = 1.0;
47*8b438e21SMatthew G. Knepley       ierr = PetscOptionsRealArray("-transform_data", "Parameters for mesh transforms", "ex33.c", options->transformDataReal, &n, NULL);CHKERRQ(ierr);
48*8b438e21SMatthew G. Knepley       break;
49*8b438e21SMatthew G. Knepley     case TRANSFORM_ANNULUS:
50*8b438e21SMatthew G. Knepley       n = 2;
51*8b438e21SMatthew G. Knepley       ierr = PetscMalloc1(n, &options->transformData);CHKERRQ(ierr);
52*8b438e21SMatthew G. Knepley       options->transformData[0] = 1.0;
53*8b438e21SMatthew G. Knepley       options->transformData[1] = 2.0;
54*8b438e21SMatthew G. Knepley       ierr = PetscOptionsScalarArray("-transform_data", "Parameters for mesh transforms", "ex33.c", options->transformData, &n, NULL);CHKERRQ(ierr);
55*8b438e21SMatthew G. Knepley       break;
56*8b438e21SMatthew G. Knepley     case TRANSFORM_SHELL:
57*8b438e21SMatthew G. Knepley       n = 2;
58*8b438e21SMatthew G. Knepley       ierr = PetscMalloc1(n, &options->transformData);CHKERRQ(ierr);
59*8b438e21SMatthew G. Knepley       options->transformData[0] = 1.0;
60*8b438e21SMatthew G. Knepley       options->transformData[1] = 2.0;
61*8b438e21SMatthew G. Knepley       ierr = PetscOptionsScalarArray("-transform_data", "Parameters for mesh transforms", "ex33.c", options->transformData, &n, NULL);CHKERRQ(ierr);
62*8b438e21SMatthew G. Knepley       break;
63*8b438e21SMatthew G. Knepley     default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Unknown mesh transform %D", options->meshTransform);
64*8b438e21SMatthew G. Knepley   }
65*8b438e21SMatthew G. Knepley   ierr = PetscOptionsReal("-volume", "The analytical volume of the mesh", "ex33.c", options->volume, &options->volume, NULL);CHKERRQ(ierr);
66*8b438e21SMatthew G. Knepley   ierr = PetscOptionsReal("-tol", "The tolerance for the volume check", "ex33.c", options->tol, &options->tol, NULL);CHKERRQ(ierr);
67*8b438e21SMatthew G. Knepley   ierr = PetscOptionsEnd();
68*8b438e21SMatthew G. Knepley   PetscFunctionReturn(0);
69*8b438e21SMatthew G. Knepley }
70*8b438e21SMatthew G. Knepley 
71*8b438e21SMatthew G. Knepley static void identity(PetscInt dim, PetscInt Nf, PetscInt NfAux,
72*8b438e21SMatthew G. Knepley                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
73*8b438e21SMatthew G. Knepley                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
74*8b438e21SMatthew G. Knepley                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
75*8b438e21SMatthew G. Knepley {
76*8b438e21SMatthew G. Knepley   const PetscInt Nc = uOff[1] - uOff[0];
77*8b438e21SMatthew G. Knepley   PetscInt       c;
78*8b438e21SMatthew G. Knepley 
79*8b438e21SMatthew G. Knepley   for (c = 0; c < Nc; ++c) f0[c] = u[c];
80*8b438e21SMatthew G. Knepley }
81*8b438e21SMatthew G. Knepley 
82*8b438e21SMatthew G. Knepley /*
83*8b438e21SMatthew 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
84*8b438e21SMatthew G. Knepley   will correspond to the top and bottom of our square. So
85*8b438e21SMatthew G. Knepley 
86*8b438e21SMatthew G. Knepley     (0,0)--(1,0)  ==>  (1,0)--(2,0)      Just a shift of (1,0)
87*8b438e21SMatthew G. Knepley     (0,1)--(1,1)  ==>  (0,1)--(0,2)      Switch x and y
88*8b438e21SMatthew G. Knepley 
89*8b438e21SMatthew 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:
90*8b438e21SMatthew G. Knepley 
91*8b438e21SMatthew G. Knepley     (x, y)  ==>  (x+1, \pi/2 y)                           in (r', \theta') space
92*8b438e21SMatthew G. Knepley             ==>  ((x+1) cos(\pi/2 y), (x+1) sin(\pi/2 y)) in (x', y') space
93*8b438e21SMatthew G. Knepley */
94*8b438e21SMatthew G. Knepley static void f0_annulus(PetscInt dim, PetscInt Nf, PetscInt NfAux,
95*8b438e21SMatthew G. Knepley                        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
96*8b438e21SMatthew G. Knepley                        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
97*8b438e21SMatthew G. Knepley                        PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar xp[])
98*8b438e21SMatthew G. Knepley {
99*8b438e21SMatthew G. Knepley   const PetscReal ri = PetscRealPart(constants[0]);
100*8b438e21SMatthew G. Knepley   const PetscReal ro = PetscRealPart(constants[1]);
101*8b438e21SMatthew G. Knepley 
102*8b438e21SMatthew G. Knepley   xp[0] = (x[0] * (ro-ri) + ri) * PetscCosReal(0.5*PETSC_PI*x[1]);
103*8b438e21SMatthew G. Knepley   xp[1] = (x[0] * (ro-ri) + ri) * PetscSinReal(0.5*PETSC_PI*x[1]);
104*8b438e21SMatthew G. Knepley }
105*8b438e21SMatthew G. Knepley 
106*8b438e21SMatthew G. Knepley /*
107*8b438e21SMatthew 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
108*8b438e21SMatthew G. Knepley   lower hemisphere and the upper surface onto the top, letting z be the radius.
109*8b438e21SMatthew G. Knepley 
110*8b438e21SMatthew G. Knepley     (x, y)  ==>  ((z+3)/2, \pi/2 (|x| or |y|), arctan y/x)                                                  in (r', \theta', \phi') space
111*8b438e21SMatthew 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
112*8b438e21SMatthew G. Knepley */
113*8b438e21SMatthew G. Knepley static void f0_shell(PetscInt dim, PetscInt Nf, PetscInt NfAux,
114*8b438e21SMatthew G. Knepley                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
115*8b438e21SMatthew G. Knepley                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
116*8b438e21SMatthew G. Knepley                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar xp[])
117*8b438e21SMatthew G. Knepley {
118*8b438e21SMatthew G. Knepley   const PetscReal pi4    = PETSC_PI/4.0;
119*8b438e21SMatthew G. Knepley   const PetscReal ri     = PetscRealPart(constants[0]);
120*8b438e21SMatthew G. Knepley   const PetscReal ro     = PetscRealPart(constants[1]);
121*8b438e21SMatthew G. Knepley   const PetscReal rp     = (x[2]+1) * 0.5*(ro-ri) + ri;
122*8b438e21SMatthew G. Knepley   const PetscReal phip   = PetscAtan2Real(x[1], x[0]);
123*8b438e21SMatthew 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])));
124*8b438e21SMatthew G. Knepley 
125*8b438e21SMatthew G. Knepley   xp[0] = rp * PetscCosReal(thetap) * PetscCosReal(phip);
126*8b438e21SMatthew G. Knepley   xp[1] = rp * PetscCosReal(thetap) * PetscSinReal(phip);
127*8b438e21SMatthew G. Knepley   xp[2] = rp * PetscSinReal(thetap);
128*8b438e21SMatthew G. Knepley }
129*8b438e21SMatthew G. Knepley 
130*8b438e21SMatthew G. Knepley static PetscErrorCode DMCreateCoordinateDisc(DM dm, PetscBool isSimplex)
131*8b438e21SMatthew G. Knepley {
132*8b438e21SMatthew G. Knepley   DM              cdm;
133*8b438e21SMatthew G. Knepley   PetscFE         fe;
134*8b438e21SMatthew G. Knepley   PetscSpace      basis;
135*8b438e21SMatthew G. Knepley   PetscDualSpace  dual;
136*8b438e21SMatthew G. Knepley   PetscInt        dim, dE;
137*8b438e21SMatthew G. Knepley   PetscErrorCode  ierr;
138*8b438e21SMatthew G. Knepley 
139*8b438e21SMatthew G. Knepley   PetscFunctionBegin;
140*8b438e21SMatthew G. Knepley   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
141*8b438e21SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
142*8b438e21SMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dE);CHKERRQ(ierr);
143*8b438e21SMatthew G. Knepley   ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, dE, isSimplex, "geom_", -1, &fe);CHKERRQ(ierr);
144*8b438e21SMatthew G. Knepley   ierr = DMSetField(cdm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr);
145*8b438e21SMatthew G. Knepley   ierr = PetscFEGetBasisSpace(fe, &basis);CHKERRQ(ierr);
146*8b438e21SMatthew G. Knepley   ierr = PetscFEGetDualSpace(fe, &dual);CHKERRQ(ierr);
147*8b438e21SMatthew G. Knepley   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
148*8b438e21SMatthew G. Knepley   ierr = DMCreateDS(cdm);CHKERRQ(ierr);
149*8b438e21SMatthew G. Knepley 
150*8b438e21SMatthew G. Knepley   {
151*8b438e21SMatthew G. Knepley     DM      cdmLinear;
152*8b438e21SMatthew G. Knepley     PetscFE feLinear;
153*8b438e21SMatthew G. Knepley     Mat     In;
154*8b438e21SMatthew G. Knepley     Vec     coordinates, coordinatesNew;
155*8b438e21SMatthew G. Knepley 
156*8b438e21SMatthew G. Knepley     /* Have to clear out previous coordinate structures */
157*8b438e21SMatthew G. Knepley     ierr = DMGetCoordinates(dm, &coordinates); CHKERRQ(ierr);
158*8b438e21SMatthew G. Knepley     ierr = DMSetCoordinateField(dm, NULL);CHKERRQ(ierr);
159*8b438e21SMatthew G. Knepley     ierr = DMSetLocalSection(cdm, NULL);CHKERRQ(ierr);
160*8b438e21SMatthew G. Knepley     ierr = DMSetGlobalSection(cdm, NULL);CHKERRQ(ierr);
161*8b438e21SMatthew G. Knepley     ierr = DMSetSectionSF(cdm, NULL);CHKERRQ(ierr);
162*8b438e21SMatthew G. Knepley     /* Project from vertex coordinates to chosen space */
163*8b438e21SMatthew G. Knepley     ierr = DMClone(cdm, &cdmLinear);CHKERRQ(ierr);
164*8b438e21SMatthew G. Knepley     ierr = PetscFECreateLagrange(PETSC_COMM_SELF, dim, dE, isSimplex, 1, -1, &feLinear);CHKERRQ(ierr);
165*8b438e21SMatthew G. Knepley     ierr = DMSetField(cdmLinear, 0, NULL, (PetscObject) feLinear);CHKERRQ(ierr);
166*8b438e21SMatthew G. Knepley     ierr = PetscFEDestroy(&feLinear);CHKERRQ(ierr);
167*8b438e21SMatthew G. Knepley     ierr = DMCreateDS(cdmLinear);CHKERRQ(ierr);
168*8b438e21SMatthew G. Knepley     ierr = DMCreateInterpolation(cdmLinear, cdm, &In, NULL); CHKERRQ(ierr);
169*8b438e21SMatthew G. Knepley     ierr = DMCreateGlobalVector(cdm, &coordinatesNew); CHKERRQ(ierr);
170*8b438e21SMatthew G. Knepley     ierr = MatInterpolate(In, coordinates, coordinatesNew); CHKERRQ(ierr);
171*8b438e21SMatthew G. Knepley     ierr = MatDestroy(&In);CHKERRQ(ierr);
172*8b438e21SMatthew G. Knepley     ierr = DMSetCoordinates(dm, coordinatesNew);CHKERRQ(ierr);
173*8b438e21SMatthew G. Knepley     ierr = VecViewFromOptions(coordinatesNew, NULL, "-coords_view");CHKERRQ(ierr);
174*8b438e21SMatthew G. Knepley     ierr = VecDestroy(&coordinatesNew);CHKERRQ(ierr);
175*8b438e21SMatthew G. Knepley     ierr = DMDestroy(&cdmLinear);CHKERRQ(ierr);
176*8b438e21SMatthew G. Knepley   }
177*8b438e21SMatthew G. Knepley   PetscFunctionReturn(0);
178*8b438e21SMatthew G. Knepley }
179*8b438e21SMatthew G. Knepley 
180*8b438e21SMatthew G. Knepley PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm)
181*8b438e21SMatthew G. Knepley {
182*8b438e21SMatthew G. Knepley   PetscInt       dim      = ctx->dim;
183*8b438e21SMatthew G. Knepley   PetscBool      simplex  = ctx->simplex;
184*8b438e21SMatthew G. Knepley   const char    *filename = ctx->filename;
185*8b438e21SMatthew G. Knepley   DM             cdm;
186*8b438e21SMatthew G. Knepley   PetscDS        cds;
187*8b438e21SMatthew G. Knepley   size_t         len;
188*8b438e21SMatthew G. Knepley   PetscMPIInt    rank;
189*8b438e21SMatthew G. Knepley   PetscErrorCode ierr;
190*8b438e21SMatthew G. Knepley 
191*8b438e21SMatthew G. Knepley   PetscFunctionBegin;
192*8b438e21SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
193*8b438e21SMatthew G. Knepley   ierr = PetscStrlen(filename, &len);CHKERRQ(ierr);
194*8b438e21SMatthew G. Knepley   if (len) {
195*8b438e21SMatthew G. Knepley     ierr = DMPlexCreateFromFile(comm, filename, PETSC_FALSE, dm);CHKERRQ(ierr);
196*8b438e21SMatthew G. Knepley     ierr = DMGetDimension(*dm, &dim);CHKERRQ(ierr);
197*8b438e21SMatthew G. Knepley   } else {
198*8b438e21SMatthew G. Knepley     ierr = DMPlexCreateBoxMesh(comm, dim, simplex, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);
199*8b438e21SMatthew G. Knepley   }
200*8b438e21SMatthew G. Knepley   {
201*8b438e21SMatthew G. Knepley     DM               pdm = NULL;
202*8b438e21SMatthew G. Knepley     PetscPartitioner part;
203*8b438e21SMatthew G. Knepley 
204*8b438e21SMatthew G. Knepley     ierr = DMPlexGetPartitioner(*dm, &part);CHKERRQ(ierr);
205*8b438e21SMatthew G. Knepley     ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
206*8b438e21SMatthew G. Knepley     /* Distribute mesh over processes */
207*8b438e21SMatthew G. Knepley     ierr = DMPlexDistribute(*dm, 0, NULL, &pdm);CHKERRQ(ierr);
208*8b438e21SMatthew G. Knepley     if (pdm) {
209*8b438e21SMatthew G. Knepley       ierr = DMDestroy(dm);CHKERRQ(ierr);
210*8b438e21SMatthew G. Knepley       *dm  = pdm;
211*8b438e21SMatthew G. Knepley     }
212*8b438e21SMatthew G. Knepley   }
213*8b438e21SMatthew G. Knepley   ierr = DMLocalizeCoordinates(*dm);CHKERRQ(ierr);
214*8b438e21SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr);
215*8b438e21SMatthew G. Knepley   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
216*8b438e21SMatthew G. Knepley   ierr = DMCreateCoordinateDisc(*dm, simplex);CHKERRQ(ierr);
217*8b438e21SMatthew G. Knepley   switch (ctx->meshTransform) {
218*8b438e21SMatthew G. Knepley     case TRANSFORM_NONE:
219*8b438e21SMatthew G. Knepley       ierr = DMPlexRemapGeometry(*dm, 0.0, identity);CHKERRQ(ierr);
220*8b438e21SMatthew G. Knepley       break;
221*8b438e21SMatthew G. Knepley     case TRANSFORM_SHEAR:
222*8b438e21SMatthew G. Knepley       ierr = DMPlexShearGeometry(*dm, DM_X, ctx->transformDataReal);CHKERRQ(ierr);
223*8b438e21SMatthew G. Knepley       break;
224*8b438e21SMatthew G. Knepley     case TRANSFORM_ANNULUS:
225*8b438e21SMatthew G. Knepley       ierr = DMGetCoordinateDM(*dm, &cdm);CHKERRQ(ierr);
226*8b438e21SMatthew G. Knepley       ierr = DMGetDS(cdm, &cds);CHKERRQ(ierr);
227*8b438e21SMatthew G. Knepley       ierr = PetscDSSetConstants(cds, 2, ctx->transformData);CHKERRQ(ierr);
228*8b438e21SMatthew G. Knepley       ierr = DMPlexRemapGeometry(*dm, 0.0, f0_annulus);CHKERRQ(ierr);
229*8b438e21SMatthew G. Knepley       break;
230*8b438e21SMatthew G. Knepley     case TRANSFORM_SHELL:
231*8b438e21SMatthew G. Knepley       ierr = DMGetCoordinateDM(*dm, &cdm);CHKERRQ(ierr);
232*8b438e21SMatthew G. Knepley       ierr = DMGetDS(cdm, &cds);CHKERRQ(ierr);
233*8b438e21SMatthew G. Knepley       ierr = PetscDSSetConstants(cds, 2, ctx->transformData);CHKERRQ(ierr);
234*8b438e21SMatthew G. Knepley       ierr = DMPlexRemapGeometry(*dm, 0.0, f0_shell);CHKERRQ(ierr);
235*8b438e21SMatthew G. Knepley       break;
236*8b438e21SMatthew G. Knepley     default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Unknown mesh transform %D", ctx->meshTransform);
237*8b438e21SMatthew G. Knepley   }
238*8b438e21SMatthew G. Knepley   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
239*8b438e21SMatthew G. Knepley   ctx->dm = *dm;
240*8b438e21SMatthew G. Knepley   PetscFunctionReturn(0);
241*8b438e21SMatthew G. Knepley }
242*8b438e21SMatthew G. Knepley 
243*8b438e21SMatthew G. Knepley static void volume(PetscInt dim, PetscInt Nf, PetscInt NfAux,
244*8b438e21SMatthew G. Knepley                    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
245*8b438e21SMatthew G. Knepley                    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
246*8b438e21SMatthew G. Knepley                    PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar vol[])
247*8b438e21SMatthew G. Knepley {
248*8b438e21SMatthew G. Knepley   vol[0] = 1.;
249*8b438e21SMatthew G. Knepley }
250*8b438e21SMatthew G. Knepley 
251*8b438e21SMatthew G. Knepley static PetscErrorCode CreateDiscretization(DM dm, AppCtx *ctx)
252*8b438e21SMatthew G. Knepley {
253*8b438e21SMatthew G. Knepley   PetscDS        ds;
254*8b438e21SMatthew G. Knepley   PetscFE        fe;
255*8b438e21SMatthew G. Knepley   PetscErrorCode ierr;
256*8b438e21SMatthew G. Knepley 
257*8b438e21SMatthew G. Knepley   PetscFunctionBeginUser;
258*8b438e21SMatthew G. Knepley   ierr = PetscFECreateDefault(PETSC_COMM_SELF, ctx->dim, 1, ctx->simplex, NULL, PETSC_DETERMINE, &fe);CHKERRQ(ierr);
259*8b438e21SMatthew G. Knepley   ierr = PetscFESetName(fe, "scalar");CHKERRQ(ierr);
260*8b438e21SMatthew G. Knepley   ierr = DMAddField(dm, NULL, (PetscObject) fe);
261*8b438e21SMatthew G. Knepley   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
262*8b438e21SMatthew G. Knepley   ierr = DMCreateDS(dm);CHKERRQ(ierr);
263*8b438e21SMatthew G. Knepley   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
264*8b438e21SMatthew G. Knepley   ierr = PetscDSSetObjective(ds, 0, volume);CHKERRQ(ierr);
265*8b438e21SMatthew G. Knepley   PetscFunctionReturn(0);
266*8b438e21SMatthew G. Knepley }
267*8b438e21SMatthew G. Knepley 
268*8b438e21SMatthew G. Knepley static PetscErrorCode CheckVolume(DM dm, AppCtx *ctx)
269*8b438e21SMatthew G. Knepley {
270*8b438e21SMatthew G. Knepley   Vec            u;
271*8b438e21SMatthew G. Knepley   PetscScalar    result;
272*8b438e21SMatthew G. Knepley   PetscReal      vol, tol = ctx->tol;
273*8b438e21SMatthew G. Knepley   PetscErrorCode ierr;
274*8b438e21SMatthew G. Knepley 
275*8b438e21SMatthew G. Knepley   PetscFunctionBeginUser;
276*8b438e21SMatthew G. Knepley   ierr = DMGetGlobalVector(dm, &u);CHKERRQ(ierr);
277*8b438e21SMatthew G. Knepley   ierr = DMPlexComputeIntegralFEM(dm, u, &result, ctx);CHKERRQ(ierr);
278*8b438e21SMatthew G. Knepley   vol  = PetscRealPart(result);
279*8b438e21SMatthew G. Knepley   ierr = DMRestoreGlobalVector(dm, &u);CHKERRQ(ierr);
280*8b438e21SMatthew G. Knepley   ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Volume: %g\n", (double) vol);CHKERRQ(ierr);
281*8b438e21SMatthew G. Knepley   if (ctx->volume > 0.0 && PetscAbsReal(ctx->volume - vol) > tol) {
282*8b438e21SMatthew 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);
283*8b438e21SMatthew G. Knepley   }
284*8b438e21SMatthew G. Knepley   PetscFunctionReturn(0);
285*8b438e21SMatthew G. Knepley }
286*8b438e21SMatthew G. Knepley 
287*8b438e21SMatthew G. Knepley int main(int argc, char **argv)
288*8b438e21SMatthew G. Knepley {
289*8b438e21SMatthew G. Knepley   AppCtx         user;
290*8b438e21SMatthew G. Knepley   PetscErrorCode ierr;
291*8b438e21SMatthew G. Knepley 
292*8b438e21SMatthew G. Knepley   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
293*8b438e21SMatthew G. Knepley   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
294*8b438e21SMatthew G. Knepley   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &user.dm);CHKERRQ(ierr);
295*8b438e21SMatthew G. Knepley   ierr = CreateDiscretization(user.dm, &user);CHKERRQ(ierr);
296*8b438e21SMatthew G. Knepley   ierr = CheckVolume(user.dm, &user);CHKERRQ(ierr);
297*8b438e21SMatthew G. Knepley   ierr = DMDestroy(&user.dm);CHKERRQ(ierr);
298*8b438e21SMatthew G. Knepley   ierr = PetscFree(user.transformDataReal);CHKERRQ(ierr);
299*8b438e21SMatthew G. Knepley   ierr = PetscFree(user.transformData);CHKERRQ(ierr);
300*8b438e21SMatthew G. Knepley   ierr = PetscFinalize();
301*8b438e21SMatthew G. Knepley   return ierr;
302*8b438e21SMatthew G. Knepley }
303*8b438e21SMatthew G. Knepley 
304*8b438e21SMatthew G. Knepley /*TEST
305*8b438e21SMatthew G. Knepley 
306*8b438e21SMatthew G. Knepley   test:
307*8b438e21SMatthew G. Knepley     suffix: square_0
308*8b438e21SMatthew G. Knepley     args: -dim 2 -simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower -1.,-1. -dm_plex_box_upper 1.,1. -geom_petscspace_degree 1 -volume 4.
309*8b438e21SMatthew G. Knepley 
310*8b438e21SMatthew G. Knepley   test:
311*8b438e21SMatthew G. Knepley     suffix: square_1
312*8b438e21SMatthew G. Knepley     args: -dim 2 -simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower -1.,-1. -dm_plex_box_upper 1.,1. -geom_petscspace_degree 2 -volume 4.
313*8b438e21SMatthew G. Knepley 
314*8b438e21SMatthew G. Knepley   test:
315*8b438e21SMatthew G. Knepley     suffix: square_2
316*8b438e21SMatthew G. Knepley     args: -dim 2 -simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower -1.,-1. -dm_plex_box_upper 1.,1. -dm_refine 1 -geom_petscspace_degree 1 -volume 4.
317*8b438e21SMatthew G. Knepley 
318*8b438e21SMatthew G. Knepley   test:
319*8b438e21SMatthew G. Knepley     suffix: square_3
320*8b438e21SMatthew G. Knepley     args: -dim 2 -simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower -1.,-1. -dm_plex_box_upper 1.,1. -dm_refine 1 -geom_petscspace_degree 2 -volume 4.
321*8b438e21SMatthew G. Knepley 
322*8b438e21SMatthew G. Knepley   test:
323*8b438e21SMatthew G. Knepley     suffix: cube_0
324*8b438e21SMatthew G. Knepley     args: -dim 3 -simplex 0 -dm_plex_box_faces 1,1,1 -dm_plex_box_lower -1.,-1.,-1. -dm_plex_box_upper 1.,1.,1. -geom_petscspace_degree 1 -volume 8.
325*8b438e21SMatthew G. Knepley 
326*8b438e21SMatthew G. Knepley   test:
327*8b438e21SMatthew G. Knepley     suffix: cube_1
328*8b438e21SMatthew G. Knepley     args: -dim 3 -simplex 0 -dm_plex_box_faces 1,1,1 -dm_plex_box_lower -1.,-1.,-1. -dm_plex_box_upper 1.,1.,1. -geom_petscspace_degree 2 -volume 8.
329*8b438e21SMatthew G. Knepley 
330*8b438e21SMatthew G. Knepley   test:
331*8b438e21SMatthew G. Knepley     suffix: cube_2
332*8b438e21SMatthew G. Knepley     args: -dim 3 -simplex 0 -dm_plex_box_faces 1,1,1 -dm_plex_box_lower -1.,-1.,-1. -dm_plex_box_upper 1.,1.,1. -dm_refine 1 -geom_petscspace_degree 1 -volume 8.
333*8b438e21SMatthew G. Knepley 
334*8b438e21SMatthew G. Knepley   test:
335*8b438e21SMatthew G. Knepley     suffix: cube_3
336*8b438e21SMatthew G. Knepley     args: -dim 3 -simplex 0 -dm_plex_box_faces 1,1,1 -dm_plex_box_lower -1.,-1.,-1. -dm_plex_box_upper 1.,1.,1. -dm_refine 1 -geom_petscspace_degree 2 -volume 8.
337*8b438e21SMatthew G. Knepley 
338*8b438e21SMatthew G. Knepley   test:
339*8b438e21SMatthew G. Knepley     suffix: shear_0
340*8b438e21SMatthew G. Knepley     args: -dim 2 -simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower -1.,-1. -dm_plex_box_upper 1.,1. -geom_petscspace_degree 1 -mesh_transform shear -transform_data 3.0 -volume 4.
341*8b438e21SMatthew G. Knepley 
342*8b438e21SMatthew G. Knepley   test:
343*8b438e21SMatthew G. Knepley     suffix: shear_1
344*8b438e21SMatthew G. Knepley     args: -dim 2 -simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower -1.,-1. -dm_plex_box_upper 1.,1. -geom_petscspace_degree 2 -mesh_transform shear -transform_data 3.0 -volume 4.
345*8b438e21SMatthew G. Knepley 
346*8b438e21SMatthew G. Knepley   test:
347*8b438e21SMatthew G. Knepley     suffix: shear_2
348*8b438e21SMatthew G. Knepley     args: -dim 2 -simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower -1.,-1. -dm_plex_box_upper 1.,1. -dm_refine 1 -geom_petscspace_degree 1 -mesh_transform shear -transform_data 3.0 -volume 4.
349*8b438e21SMatthew G. Knepley 
350*8b438e21SMatthew G. Knepley   test:
351*8b438e21SMatthew G. Knepley     suffix: shear_3
352*8b438e21SMatthew G. Knepley     args: -dim 2 -simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower -1.,-1. -dm_plex_box_upper 1.,1. -dm_refine 1 -geom_petscspace_degree 2 -mesh_transform shear -transform_data 3.0 -volume 4.
353*8b438e21SMatthew G. Knepley 
354*8b438e21SMatthew G. Knepley   test:
355*8b438e21SMatthew G. Knepley     suffix: shear_4
356*8b438e21SMatthew G. Knepley     args: -dim 3 -simplex 0 -dm_plex_box_faces 1,1,1 -dm_plex_box_lower -1.,-1.,-1. -dm_plex_box_upper 1.,1.,1. -geom_petscspace_degree 1 -mesh_transform shear -transform_data 3.0 -volume 8.
357*8b438e21SMatthew G. Knepley 
358*8b438e21SMatthew G. Knepley   test:
359*8b438e21SMatthew G. Knepley     suffix: shear_5
360*8b438e21SMatthew G. Knepley     args: -dim 3 -simplex 0 -dm_plex_box_faces 1,1,1 -dm_plex_box_lower -1.,-1.,-1. -dm_plex_box_upper 1.,1.,1. -geom_petscspace_degree 2 -mesh_transform shear -transform_data 3.0 -volume 8.
361*8b438e21SMatthew G. Knepley 
362*8b438e21SMatthew G. Knepley   test:
363*8b438e21SMatthew G. Knepley     suffix: shear_6
364*8b438e21SMatthew G. Knepley     args: -dim 3 -simplex 0 -dm_plex_box_faces 1,1,1 -dm_plex_box_lower -1.,-1,-1. -dm_plex_box_upper 1.,1.,1. -dm_refine 1 -geom_petscspace_degree 1 -mesh_transform shear -transform_data 3.0,4.0 -volume 8.
365*8b438e21SMatthew G. Knepley 
366*8b438e21SMatthew G. Knepley   test:
367*8b438e21SMatthew G. Knepley     suffix: shear_7
368*8b438e21SMatthew G. Knepley     args: -dim 3 -simplex 0 -dm_plex_box_faces 1,1,1 -dm_plex_box_lower -1.,-1.,-1. -dm_plex_box_upper 1.,1.,1. -dm_refine 1 -geom_petscspace_degree 2 -mesh_transform shear -transform_data 3.0,4.0 -volume 8.
369*8b438e21SMatthew G. Knepley 
370*8b438e21SMatthew G. Knepley   test:
371*8b438e21SMatthew G. Knepley     # Area: (a+b)/2 h = 3/\sqrt{2} (sqrt{2} - 1/\sqrt{2}) = 3/2
372*8b438e21SMatthew G. Knepley     suffix: annulus_0
373*8b438e21SMatthew G. Knepley     requires: double
374*8b438e21SMatthew G. Knepley     args: -dim 2 -simplex 0 -dm_plex_box_faces 1,1 -geom_petscspace_degree 1 -mesh_transform annulus -volume 1.5
375*8b438e21SMatthew G. Knepley 
376*8b438e21SMatthew G. Knepley   test:
377*8b438e21SMatthew G. Knepley     # Area: 3/4 \pi = 2.3562
378*8b438e21SMatthew G. Knepley     suffix: annulus_1
379*8b438e21SMatthew G. Knepley     requires: double
380*8b438e21SMatthew G. Knepley     args: -dim 2 -simplex 0 -dm_plex_box_faces 1,1 -dm_refine 3 -geom_petscspace_degree 1 -mesh_transform annulus -volume 2.35619449019235 -tol .016
381*8b438e21SMatthew G. Knepley 
382*8b438e21SMatthew G. Knepley   test:
383*8b438e21SMatthew G. Knepley     # Area: 3/4 \pi = 2.3562
384*8b438e21SMatthew G. Knepley     suffix: annulus_2
385*8b438e21SMatthew G. Knepley     requires: double
386*8b438e21SMatthew G. Knepley     args: -dim 2 -simplex 0 -dm_plex_box_faces 1,1 -dm_refine 3 -geom_petscspace_degree 2 -mesh_transform annulus -volume 2.35619449019235 -tol .0038
387*8b438e21SMatthew G. Knepley 
388*8b438e21SMatthew G. Knepley   test:
389*8b438e21SMatthew G. Knepley     # Area: 3/4 \pi = 2.3562
390*8b438e21SMatthew G. Knepley     suffix: annulus_3
391*8b438e21SMatthew G. Knepley     requires: double
392*8b438e21SMatthew G. Knepley     args: -dim 2 -simplex 0 -dm_plex_box_faces 1,1 -dm_refine 3 -geom_petscspace_degree 3 -mesh_transform annulus -volume 2.35619449019235 -tol 2.2e-6
393*8b438e21SMatthew G. Knepley 
394*8b438e21SMatthew G. Knepley   test:
395*8b438e21SMatthew G. Knepley     # Area: 3/4 \pi = 2.3562
396*8b438e21SMatthew G. Knepley     suffix: annulus_4
397*8b438e21SMatthew G. Knepley     requires: double
398*8b438e21SMatthew G. Knepley     args: -dim 2 -simplex 0 -dm_plex_box_faces 1,1 -dm_refine 2 -geom_petscspace_degree 2 -petscfe_default_quadrature_order 2 -mesh_transform annulus -volume 2.35619449019235 -tol .00012
399*8b438e21SMatthew G. Knepley 
400*8b438e21SMatthew G. Knepley   test:
401*8b438e21SMatthew G. Knepley     # Area: 3/4 \pi = 2.3562
402*8b438e21SMatthew G. Knepley     suffix: annulus_5
403*8b438e21SMatthew G. Knepley     requires: double
404*8b438e21SMatthew G. Knepley     args: -dim 2 -simplex 0 -dm_plex_box_faces 1,1 -dm_refine 2 -geom_petscspace_degree 3 -petscfe_default_quadrature_order 3 -mesh_transform annulus -volume 2.35619449019235 -tol 1.2e-7
405*8b438e21SMatthew G. Knepley 
406*8b438e21SMatthew G. Knepley   test:
407*8b438e21SMatthew G. Knepley     suffix: shell_0
408*8b438e21SMatthew G. Knepley     requires: double
409*8b438e21SMatthew G. Knepley     args: -dim 3 -simplex 0 -dm_plex_box_faces 1,1,1 -dm_plex_box_lower -1.,-1.,-1. -dm_plex_box_upper 1.,1.,1. -dm_refine 1 -geom_petscspace_degree 1 -petscfe_default_quadrature_order 1 -mesh_transform shell -volume 5.633164922 -tol 1.0e-7
410*8b438e21SMatthew G. Knepley 
411*8b438e21SMatthew G. Knepley   test:
412*8b438e21SMatthew G. Knepley     # Volume: 4/3 \pi (8 - 1)/2 = 14/3 \pi = 14.66076571675238
413*8b438e21SMatthew G. Knepley     suffix: shell_1
414*8b438e21SMatthew G. Knepley     requires: double
415*8b438e21SMatthew G. Knepley     args: -dim 3 -simplex 0 -dm_plex_box_faces 1,1,1 -dm_plex_box_lower -1.,-1.,-1. -dm_plex_box_upper 1.,1.,1. -dm_refine 2 -geom_petscspace_degree 1 -petscfe_default_quadrature_order 1 -mesh_transform shell -volume 14.66076571675238 -tol 3.1
416*8b438e21SMatthew G. Knepley 
417*8b438e21SMatthew G. Knepley   test:
418*8b438e21SMatthew G. Knepley     # Volume: 4/3 \pi (8 - 1)/2 = 14/3 \pi = 14.66076571675238
419*8b438e21SMatthew G. Knepley     suffix: shell_2
420*8b438e21SMatthew G. Knepley     requires: double
421*8b438e21SMatthew G. Knepley     args: -dim 3 -simplex 0 -dm_plex_box_faces 1,1,1 -dm_plex_box_lower -1.,-1.,-1. -dm_plex_box_upper 1.,1.,1. -dm_refine 2 -geom_petscspace_degree 2 -petscfe_default_quadrature_order 2 -mesh_transform shell -volume 14.66076571675238 -tol .1
422*8b438e21SMatthew G. Knepley 
423*8b438e21SMatthew G. Knepley   test:
424*8b438e21SMatthew G. Knepley     # Volume: 4/3 \pi (8 - 1)/2 = 14/3 \pi = 14.66076571675238
425*8b438e21SMatthew G. Knepley     suffix: shell_3
426*8b438e21SMatthew G. Knepley     requires: double
427*8b438e21SMatthew G. Knepley     args: -dim 3 -simplex 0 -dm_plex_box_faces 1,1,1 -dm_plex_box_lower -1.,-1.,-1. -dm_plex_box_upper 1.,1.,1. -dm_refine 2 -geom_petscspace_degree 3 -petscfe_default_quadrature_order 3 -mesh_transform shell -volume 14.66076571675238 -tol .02
428*8b438e21SMatthew G. Knepley 
429*8b438e21SMatthew G. Knepley   test:
430*8b438e21SMatthew G. Knepley     # Volume: 4/3 \pi (8 - 1)/2 = 14/3 \pi = 14.66076571675238
431*8b438e21SMatthew G. Knepley     suffix: shell_4
432*8b438e21SMatthew G. Knepley     requires: double
433*8b438e21SMatthew G. Knepley     args: -dim 3 -simplex 0 -dm_plex_box_faces 1,1,1 -dm_plex_box_lower -1.,-1.,-1. -dm_plex_box_upper 1.,1.,1. -dm_refine 2 -geom_petscspace_degree 4 -petscfe_default_quadrature_order 4 -mesh_transform shell -volume 14.66076571675238 -tol .006
434*8b438e21SMatthew G. Knepley 
435*8b438e21SMatthew G. Knepley TEST*/
436