1c4762a1bSJed Brown static char help[] = "Check that a DM can accurately represent and interpolate functions of a given polynomial order\n\n";
2c4762a1bSJed Brown
3c4762a1bSJed Brown #include <petscdmplex.h>
4c4762a1bSJed Brown #include <petscdm.h>
5c4762a1bSJed Brown #include <petscdmda.h>
6c4762a1bSJed Brown #include <petscfe.h>
7c4762a1bSJed Brown #include <petscds.h>
8c4762a1bSJed Brown #include <petscksp.h>
9c4762a1bSJed Brown #include <petscsnes.h>
104446c3cdSksagiyam #include <petscsf.h>
11c4762a1bSJed Brown
12c4762a1bSJed Brown typedef struct {
13c4762a1bSJed Brown /* Domain and mesh definition */
14c4762a1bSJed Brown PetscBool useDA; /* Flag DMDA tensor product mesh */
15c4762a1bSJed Brown PetscBool shearCoords; /* Flag for shear transform */
16c4762a1bSJed Brown PetscBool nonaffineCoords; /* Flag for non-affine transform */
17c4762a1bSJed Brown /* Element definition */
18c4762a1bSJed Brown PetscInt qorder; /* Order of the quadrature */
19c4762a1bSJed Brown PetscInt numComponents; /* Number of field components */
20c4762a1bSJed Brown PetscFE fe; /* The finite element */
21c4762a1bSJed Brown /* Testing space */
22c4762a1bSJed Brown PetscInt porder; /* Order of polynomials to test */
236c2660fcSDarsh Nathawani PetscBool RT; /* Test for Raviart-Thomas elements */
24c4762a1bSJed Brown PetscBool convergence; /* Test for order of convergence */
25c4762a1bSJed Brown PetscBool convRefine; /* Test for convergence using refinement, otherwise use coarsening */
26c4762a1bSJed Brown PetscBool constraints; /* Test local constraints */
27c4762a1bSJed Brown PetscBool tree; /* Test tree routines */
28c4762a1bSJed Brown PetscBool testFEjacobian; /* Test finite element Jacobian assembly */
29c4762a1bSJed Brown PetscBool testFVgrad; /* Test finite difference gradient routine */
30c4762a1bSJed Brown PetscBool testInjector; /* Test finite element injection routines */
31c4762a1bSJed Brown PetscInt treeCell; /* Cell to refine in tree test */
32c4762a1bSJed Brown PetscReal constants[3]; /* Constant values for each dimension */
33c4762a1bSJed Brown } AppCtx;
34c4762a1bSJed Brown
356c2660fcSDarsh Nathawani /*
366c2660fcSDarsh Nathawani Derivatives are set as n_i \partial u_j / \partial x_i
376c2660fcSDarsh Nathawani */
386c2660fcSDarsh Nathawani
39c4762a1bSJed Brown /* u = 1 */
constant(PetscInt dim,PetscReal time,const PetscReal coords[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)40*2a8381b2SBarry Smith PetscErrorCode constant(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
41d71ae5a4SJacob Faibussowitsch {
42c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx;
43c4762a1bSJed Brown PetscInt d;
4430602db0SMatthew G. Knepley for (d = 0; d < dim; ++d) u[d] = user->constants[d];
453ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
46c4762a1bSJed Brown }
constantDer(PetscInt dim,PetscReal time,const PetscReal coords[],const PetscReal n[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)47*2a8381b2SBarry Smith PetscErrorCode constantDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
48d71ae5a4SJacob Faibussowitsch {
49c4762a1bSJed Brown PetscInt d;
5030602db0SMatthew G. Knepley for (d = 0; d < dim; ++d) u[d] = 0.0;
513ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
52c4762a1bSJed Brown }
53c4762a1bSJed Brown
546c2660fcSDarsh Nathawani /* RT_0: u = (1 + x, 1 + y) or (1 + x, 1 + y, 1 + z) */
rt0(PetscInt dim,PetscReal time,const PetscReal coords[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)55*2a8381b2SBarry Smith PetscErrorCode rt0(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
566c2660fcSDarsh Nathawani {
576c2660fcSDarsh Nathawani PetscInt d;
586c2660fcSDarsh Nathawani for (d = 0; d < dim; ++d) u[d] = 1.0 + coords[d];
596c2660fcSDarsh Nathawani return PETSC_SUCCESS;
606c2660fcSDarsh Nathawani }
616c2660fcSDarsh Nathawani
rt0Der(PetscInt dim,PetscReal time,const PetscReal coords[],const PetscReal n[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)62*2a8381b2SBarry Smith PetscErrorCode rt0Der(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
636c2660fcSDarsh Nathawani {
646c2660fcSDarsh Nathawani PetscInt d, e;
656c2660fcSDarsh Nathawani for (d = 0; d < dim; ++d) {
666c2660fcSDarsh Nathawani u[d] = 0.0;
676c2660fcSDarsh Nathawani for (e = 0; e < dim; ++e) u[d] += (d == e ? 1.0 : 0.0) * n[e];
686c2660fcSDarsh Nathawani }
696c2660fcSDarsh Nathawani return PETSC_SUCCESS;
706c2660fcSDarsh Nathawani }
716c2660fcSDarsh Nathawani
726c2660fcSDarsh Nathawani /* u = (x + y, y + x) or (x + z, 2y, z + x) */
linear(PetscInt dim,PetscReal time,const PetscReal coords[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)73*2a8381b2SBarry Smith PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
74d71ae5a4SJacob Faibussowitsch {
75c4762a1bSJed Brown PetscInt d;
766c2660fcSDarsh Nathawani for (d = 0; d < dim; ++d) u[d] = coords[d] + coords[dim - d - 1];
773ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
78c4762a1bSJed Brown }
linearDer(PetscInt dim,PetscReal time,const PetscReal coords[],const PetscReal n[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)79*2a8381b2SBarry Smith PetscErrorCode linearDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
80d71ae5a4SJacob Faibussowitsch {
81c4762a1bSJed Brown PetscInt d, e;
82c4762a1bSJed Brown for (d = 0; d < dim; ++d) {
83c4762a1bSJed Brown u[d] = 0.0;
846c2660fcSDarsh Nathawani for (e = 0; e < dim; ++e) u[d] += ((d == e ? 1. : 0.) + (d == (dim - e - 1) ? 1. : 0.)) * n[e];
856c2660fcSDarsh Nathawani }
866c2660fcSDarsh Nathawani return PETSC_SUCCESS;
876c2660fcSDarsh Nathawani }
886c2660fcSDarsh Nathawani
896c2660fcSDarsh Nathawani /* RT_1: u = (1 + x + y + x^2 + xy, 1 + x + y + xy + y^2) or (1 + x + y + z + x^2 + xy + xz, 1 + x + y + z + xy + y^2 + yz, 1 + x + y + z + xz + yz + z^2) */
rt1(PetscInt dim,PetscReal time,const PetscReal coords[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)90*2a8381b2SBarry Smith PetscErrorCode rt1(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
916c2660fcSDarsh Nathawani {
926c2660fcSDarsh Nathawani if (dim > 2) {
936c2660fcSDarsh Nathawani u[0] = 1.0 + coords[0] + coords[1] + coords[2] + coords[0] * coords[0] + coords[0] * coords[1] + coords[0] * coords[2];
946c2660fcSDarsh Nathawani u[1] = 1.0 + coords[0] + coords[1] + coords[2] + coords[0] * coords[1] + coords[1] * coords[1] + coords[1] * coords[2];
956c2660fcSDarsh Nathawani u[2] = 1.0 + coords[0] + coords[1] + coords[2] + coords[0] * coords[2] + coords[1] * coords[2] + coords[2] * coords[2];
966c2660fcSDarsh Nathawani } else if (dim > 1) {
976c2660fcSDarsh Nathawani u[0] = 1.0 + coords[0] + coords[1] + coords[0] * coords[0] + coords[0] * coords[1];
986c2660fcSDarsh Nathawani u[1] = 1.0 + coords[0] + coords[1] + coords[0] * coords[1] + coords[1] * coords[1];
996c2660fcSDarsh Nathawani }
1006c2660fcSDarsh Nathawani return PETSC_SUCCESS;
1016c2660fcSDarsh Nathawani }
1026c2660fcSDarsh Nathawani
rt1Der(PetscInt dim,PetscReal time,const PetscReal coords[],const PetscReal n[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)103*2a8381b2SBarry Smith PetscErrorCode rt1Der(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
1046c2660fcSDarsh Nathawani {
1056c2660fcSDarsh Nathawani if (dim > 2) {
1066c2660fcSDarsh Nathawani u[0] = (1.0 + 2.0 * coords[0] + coords[1] + coords[2]) * n[0] + (1.0 + coords[0]) * n[1] + (1.0 + coords[0]) * n[2];
1076c2660fcSDarsh Nathawani u[1] = (1.0 + coords[1]) * n[0] + (1.0 + coords[0] + 2.0 * coords[1] + coords[2]) * n[1] + (1.0 + coords[1]) * n[2];
1086c2660fcSDarsh Nathawani u[2] = (1.0 + coords[2]) * n[0] + (1.0 + coords[2]) * n[1] + (1.0 + coords[0] + coords[1] + 2.0 * coords[2]) * n[2];
1096c2660fcSDarsh Nathawani } else if (dim > 1) {
1106c2660fcSDarsh Nathawani u[0] = (1.0 + 2.0 * coords[0] + coords[1]) * n[0] + (1.0 + coords[0]) * n[1];
1116c2660fcSDarsh Nathawani u[1] = (1.0 + coords[1]) * n[0] + (1.0 + coords[0] + 2.0 * coords[1]) * n[1];
112c4762a1bSJed Brown }
1133ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
114c4762a1bSJed Brown }
115c4762a1bSJed Brown
116c4762a1bSJed Brown /* u = x^2 or u = (x^2, xy) or u = (xy, yz, zx) */
quadratic(PetscInt dim,PetscReal time,const PetscReal coords[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)117*2a8381b2SBarry Smith PetscErrorCode quadratic(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
118d71ae5a4SJacob Faibussowitsch {
1199371c9d4SSatish Balay if (dim > 2) {
1209371c9d4SSatish Balay u[0] = coords[0] * coords[1];
1219371c9d4SSatish Balay u[1] = coords[1] * coords[2];
1229371c9d4SSatish Balay u[2] = coords[2] * coords[0];
1239371c9d4SSatish Balay } else if (dim > 1) {
1249371c9d4SSatish Balay u[0] = coords[0] * coords[0];
1259371c9d4SSatish Balay u[1] = coords[0] * coords[1];
1269371c9d4SSatish Balay } else if (dim > 0) {
1279371c9d4SSatish Balay u[0] = coords[0] * coords[0];
1289371c9d4SSatish Balay }
1293ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
130c4762a1bSJed Brown }
quadraticDer(PetscInt dim,PetscReal time,const PetscReal coords[],const PetscReal n[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)131*2a8381b2SBarry Smith PetscErrorCode quadraticDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
132d71ae5a4SJacob Faibussowitsch {
1339371c9d4SSatish Balay if (dim > 2) {
1349371c9d4SSatish Balay u[0] = coords[1] * n[0] + coords[0] * n[1];
1359371c9d4SSatish Balay u[1] = coords[2] * n[1] + coords[1] * n[2];
1369371c9d4SSatish Balay u[2] = coords[2] * n[0] + coords[0] * n[2];
1379371c9d4SSatish Balay } else if (dim > 1) {
1389371c9d4SSatish Balay u[0] = 2.0 * coords[0] * n[0];
1399371c9d4SSatish Balay u[1] = coords[1] * n[0] + coords[0] * n[1];
1409371c9d4SSatish Balay } else if (dim > 0) {
1419371c9d4SSatish Balay u[0] = 2.0 * coords[0] * n[0];
1429371c9d4SSatish Balay }
1433ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
144c4762a1bSJed Brown }
145c4762a1bSJed Brown
146c4762a1bSJed Brown /* u = x^3 or u = (x^3, x^2y) or u = (x^2y, y^2z, z^2x) */
cubic(PetscInt dim,PetscReal time,const PetscReal coords[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)147*2a8381b2SBarry Smith PetscErrorCode cubic(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
148d71ae5a4SJacob Faibussowitsch {
1499371c9d4SSatish Balay if (dim > 2) {
1509371c9d4SSatish Balay u[0] = coords[0] * coords[0] * coords[1];
1519371c9d4SSatish Balay u[1] = coords[1] * coords[1] * coords[2];
1529371c9d4SSatish Balay u[2] = coords[2] * coords[2] * coords[0];
1539371c9d4SSatish Balay } else if (dim > 1) {
1549371c9d4SSatish Balay u[0] = coords[0] * coords[0] * coords[0];
1559371c9d4SSatish Balay u[1] = coords[0] * coords[0] * coords[1];
1569371c9d4SSatish Balay } else if (dim > 0) {
1579371c9d4SSatish Balay u[0] = coords[0] * coords[0] * coords[0];
1589371c9d4SSatish Balay }
1593ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
160c4762a1bSJed Brown }
cubicDer(PetscInt dim,PetscReal time,const PetscReal coords[],const PetscReal n[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)161*2a8381b2SBarry Smith PetscErrorCode cubicDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
162d71ae5a4SJacob Faibussowitsch {
1639371c9d4SSatish Balay if (dim > 2) {
1649371c9d4SSatish Balay u[0] = 2.0 * coords[0] * coords[1] * n[0] + coords[0] * coords[0] * n[1];
1659371c9d4SSatish Balay u[1] = 2.0 * coords[1] * coords[2] * n[1] + coords[1] * coords[1] * n[2];
1669371c9d4SSatish Balay u[2] = 2.0 * coords[2] * coords[0] * n[2] + coords[2] * coords[2] * n[0];
1679371c9d4SSatish Balay } else if (dim > 1) {
1689371c9d4SSatish Balay u[0] = 3.0 * coords[0] * coords[0] * n[0];
1699371c9d4SSatish Balay u[1] = 2.0 * coords[0] * coords[1] * n[0] + coords[0] * coords[0] * n[1];
1709371c9d4SSatish Balay } else if (dim > 0) {
1719371c9d4SSatish Balay u[0] = 3.0 * coords[0] * coords[0] * n[0];
1729371c9d4SSatish Balay }
1733ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
174c4762a1bSJed Brown }
175c4762a1bSJed Brown
176c4762a1bSJed Brown /* u = tanh(x) */
trig(PetscInt dim,PetscReal time,const PetscReal coords[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)177*2a8381b2SBarry Smith PetscErrorCode trig(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
178d71ae5a4SJacob Faibussowitsch {
179c4762a1bSJed Brown PetscInt d;
18030602db0SMatthew G. Knepley for (d = 0; d < dim; ++d) u[d] = PetscTanhReal(coords[d] - 0.5);
1813ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
182c4762a1bSJed Brown }
trigDer(PetscInt dim,PetscReal time,const PetscReal coords[],const PetscReal n[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)183*2a8381b2SBarry Smith PetscErrorCode trigDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
184d71ae5a4SJacob Faibussowitsch {
185c4762a1bSJed Brown PetscInt d;
18630602db0SMatthew G. Knepley for (d = 0; d < dim; ++d) u[d] = 1.0 / PetscSqr(PetscCoshReal(coords[d] - 0.5)) * n[d];
1873ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
188c4762a1bSJed Brown }
189c4762a1bSJed Brown
ProcessOptions(MPI_Comm comm,AppCtx * options)190d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
191d71ae5a4SJacob Faibussowitsch {
192c4762a1bSJed Brown PetscInt n = 3;
193c4762a1bSJed Brown
194c4762a1bSJed Brown PetscFunctionBeginUser;
19530602db0SMatthew G. Knepley options->useDA = PETSC_FALSE;
196c4762a1bSJed Brown options->shearCoords = PETSC_FALSE;
197c4762a1bSJed Brown options->nonaffineCoords = PETSC_FALSE;
198c4762a1bSJed Brown options->qorder = 0;
199c4762a1bSJed Brown options->numComponents = PETSC_DEFAULT;
200c4762a1bSJed Brown options->porder = 0;
2016c2660fcSDarsh Nathawani options->RT = PETSC_FALSE;
202c4762a1bSJed Brown options->convergence = PETSC_FALSE;
203c4762a1bSJed Brown options->convRefine = PETSC_TRUE;
204c4762a1bSJed Brown options->constraints = PETSC_FALSE;
205c4762a1bSJed Brown options->tree = PETSC_FALSE;
206c4762a1bSJed Brown options->treeCell = 0;
207c4762a1bSJed Brown options->testFEjacobian = PETSC_FALSE;
208c4762a1bSJed Brown options->testFVgrad = PETSC_FALSE;
209c4762a1bSJed Brown options->testInjector = PETSC_FALSE;
210c4762a1bSJed Brown options->constants[0] = 1.0;
211c4762a1bSJed Brown options->constants[1] = 2.0;
212c4762a1bSJed Brown options->constants[2] = 3.0;
213c4762a1bSJed Brown
214d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Projection Test Options", "DMPlex");
2159566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-use_da", "Flag for DMDA mesh", "ex3.c", options->useDA, &options->useDA, NULL));
2169566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-shear_coords", "Transform coordinates with a shear", "ex3.c", options->shearCoords, &options->shearCoords, NULL));
2179566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-non_affine_coords", "Transform coordinates with a non-affine transform", "ex3.c", options->nonaffineCoords, &options->nonaffineCoords, NULL));
2189566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-qorder", "The quadrature order", "ex3.c", options->qorder, &options->qorder, NULL, 0));
2199566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-num_comp", "The number of field components", "ex3.c", options->numComponents, &options->numComponents, NULL, PETSC_DEFAULT));
2209566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-porder", "The order of polynomials to test", "ex3.c", options->porder, &options->porder, NULL, 0));
2216c2660fcSDarsh Nathawani PetscCall(PetscOptionsBool("-RT", "Use the Raviart-Thomas elements", "ex3.c", options->RT, &options->RT, NULL));
2229566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-convergence", "Check the convergence rate", "ex3.c", options->convergence, &options->convergence, NULL));
2239566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-conv_refine", "Use refinement for the convergence rate", "ex3.c", options->convRefine, &options->convRefine, NULL));
2249566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-constraints", "Test local constraints (serial only)", "ex3.c", options->constraints, &options->constraints, NULL));
2259566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-tree", "Test tree routines", "ex3.c", options->tree, &options->tree, NULL));
2269566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-tree_cell", "cell to refine in tree test", "ex3.c", options->treeCell, &options->treeCell, NULL, 0));
2279566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_fe_jacobian", "Test finite element Jacobian assembly", "ex3.c", options->testFEjacobian, &options->testFEjacobian, NULL));
2289566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_fv_grad", "Test finite volume gradient reconstruction", "ex3.c", options->testFVgrad, &options->testFVgrad, NULL));
2299566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_injector", "Test finite element injection", "ex3.c", options->testInjector, &options->testInjector, NULL));
2309566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-constants", "Set the constant values", "ex3.c", options->constants, &n, NULL));
231d0609cedSBarry Smith PetscOptionsEnd();
2323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
233c4762a1bSJed Brown }
234c4762a1bSJed Brown
TransformCoordinates(DM dm,AppCtx * user)235d71ae5a4SJacob Faibussowitsch static PetscErrorCode TransformCoordinates(DM dm, AppCtx *user)
236d71ae5a4SJacob Faibussowitsch {
237c4762a1bSJed Brown PetscSection coordSection;
238c4762a1bSJed Brown Vec coordinates;
239c4762a1bSJed Brown PetscScalar *coords;
240c4762a1bSJed Brown PetscInt vStart, vEnd, v;
241c4762a1bSJed Brown
242c4762a1bSJed Brown PetscFunctionBeginUser;
243c4762a1bSJed Brown if (user->nonaffineCoords) {
244c4762a1bSJed Brown /* x' = r^(1/p) (x/r), y' = r^(1/p) (y/r), z' = z */
2459566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm, &coordSection));
2469566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
2479566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
2489566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordinates, &coords));
249c4762a1bSJed Brown for (v = vStart; v < vEnd; ++v) {
250c4762a1bSJed Brown PetscInt dof, off;
251c4762a1bSJed Brown PetscReal p = 4.0, r;
252c4762a1bSJed Brown
2539566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(coordSection, v, &dof));
2549566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(coordSection, v, &off));
255c4762a1bSJed Brown switch (dof) {
256c4762a1bSJed Brown case 2:
257c4762a1bSJed Brown r = PetscSqr(PetscRealPart(coords[off + 0])) + PetscSqr(PetscRealPart(coords[off + 1]));
258c4762a1bSJed Brown coords[off + 0] = r == 0.0 ? 0.0 : PetscPowReal(r, (1 - p) / (2 * p)) * coords[off + 0];
259c4762a1bSJed Brown coords[off + 1] = r == 0.0 ? 0.0 : PetscPowReal(r, (1 - p) / (2 * p)) * coords[off + 1];
260c4762a1bSJed Brown break;
261c4762a1bSJed Brown case 3:
262c4762a1bSJed Brown r = PetscSqr(PetscRealPart(coords[off + 0])) + PetscSqr(PetscRealPart(coords[off + 1]));
263c4762a1bSJed Brown coords[off + 0] = r == 0.0 ? 0.0 : PetscPowReal(r, (1 - p) / (2 * p)) * coords[off + 0];
264c4762a1bSJed Brown coords[off + 1] = r == 0.0 ? 0.0 : PetscPowReal(r, (1 - p) / (2 * p)) * coords[off + 1];
265c4762a1bSJed Brown coords[off + 2] = coords[off + 2];
266c4762a1bSJed Brown break;
267c4762a1bSJed Brown }
268c4762a1bSJed Brown }
2699566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coordinates, &coords));
270c4762a1bSJed Brown }
271c4762a1bSJed Brown if (user->shearCoords) {
272c4762a1bSJed Brown /* x' = x + m y + m z, y' = y + m z, z' = z */
2739566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm, &coordSection));
2749566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
2759566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
2769566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordinates, &coords));
277c4762a1bSJed Brown for (v = vStart; v < vEnd; ++v) {
278c4762a1bSJed Brown PetscInt dof, off;
279c4762a1bSJed Brown PetscReal m = 1.0;
280c4762a1bSJed Brown
2819566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(coordSection, v, &dof));
2829566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(coordSection, v, &off));
283c4762a1bSJed Brown switch (dof) {
284c4762a1bSJed Brown case 2:
285c4762a1bSJed Brown coords[off + 0] = coords[off + 0] + m * coords[off + 1];
286c4762a1bSJed Brown coords[off + 1] = coords[off + 1];
287c4762a1bSJed Brown break;
288c4762a1bSJed Brown case 3:
289c4762a1bSJed Brown coords[off + 0] = coords[off + 0] + m * coords[off + 1] + m * coords[off + 2];
290c4762a1bSJed Brown coords[off + 1] = coords[off + 1] + m * coords[off + 2];
291c4762a1bSJed Brown coords[off + 2] = coords[off + 2];
292c4762a1bSJed Brown break;
293c4762a1bSJed Brown }
294c4762a1bSJed Brown }
2959566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coordinates, &coords));
296c4762a1bSJed Brown }
2973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
298c4762a1bSJed Brown }
299c4762a1bSJed Brown
CreateMesh(MPI_Comm comm,AppCtx * user,DM * dm)300d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
301d71ae5a4SJacob Faibussowitsch {
30230602db0SMatthew G. Knepley PetscInt dim = 2;
30330602db0SMatthew G. Knepley PetscBool simplex;
304c4762a1bSJed Brown
305c4762a1bSJed Brown PetscFunctionBeginUser;
30630602db0SMatthew G. Knepley if (user->useDA) {
30730602db0SMatthew G. Knepley switch (dim) {
308c4762a1bSJed Brown case 2:
3099566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(comm, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, 2, 2, PETSC_DETERMINE, PETSC_DETERMINE, 1, 1, NULL, NULL, dm));
3109566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm));
3119566063dSJacob Faibussowitsch PetscCall(DMSetUp(*dm));
3129566063dSJacob Faibussowitsch PetscCall(DMDASetVertexCoordinates(*dm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
313c4762a1bSJed Brown break;
314d71ae5a4SJacob Faibussowitsch default:
315d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot create structured mesh of dimension %" PetscInt_FMT, dim);
316c4762a1bSJed Brown }
3179566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*dm, "Hexahedral Mesh"));
31830602db0SMatthew G. Knepley } else {
3199566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm));
3209566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX));
3219566063dSJacob Faibussowitsch PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
3229566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm));
323c4762a1bSJed Brown
3249566063dSJacob Faibussowitsch PetscCall(DMGetDimension(*dm, &dim));
3259566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(*dm, &simplex));
3265440e5dcSBarry Smith PetscCallMPI(MPI_Bcast(&simplex, 1, MPI_C_BOOL, 0, comm));
327c4762a1bSJed Brown if (user->tree) {
32830602db0SMatthew G. Knepley DM refTree, ncdm = NULL;
329c4762a1bSJed Brown
3309566063dSJacob Faibussowitsch PetscCall(DMPlexCreateDefaultReferenceTree(comm, dim, simplex, &refTree));
3319566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(refTree, NULL, "-reftree_dm_view"));
3329566063dSJacob Faibussowitsch PetscCall(DMPlexSetReferenceTree(*dm, refTree));
3339566063dSJacob Faibussowitsch PetscCall(DMDestroy(&refTree));
3349566063dSJacob Faibussowitsch PetscCall(DMPlexTreeRefineCell(*dm, user->treeCell, &ncdm));
335c4762a1bSJed Brown if (ncdm) {
3369566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm));
337c4762a1bSJed Brown *dm = ncdm;
3389566063dSJacob Faibussowitsch PetscCall(DMPlexSetRefinementUniform(*dm, PETSC_FALSE));
339c4762a1bSJed Brown }
3409566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*dm, "tree_"));
3419566063dSJacob Faibussowitsch PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
3429566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm));
3439566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
344c4762a1bSJed Brown } else {
3459566063dSJacob Faibussowitsch PetscCall(DMPlexSetRefinementUniform(*dm, PETSC_TRUE));
346c4762a1bSJed Brown }
3479566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*dm, "dist_"));
3489566063dSJacob Faibussowitsch PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
3499566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm));
3509566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*dm, NULL));
3519566063dSJacob Faibussowitsch if (simplex) PetscCall(PetscObjectSetName((PetscObject)*dm, "Simplicial Mesh"));
3529566063dSJacob Faibussowitsch else PetscCall(PetscObjectSetName((PetscObject)*dm, "Hexahedral Mesh"));
353c4762a1bSJed Brown }
3549566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm));
3559566063dSJacob Faibussowitsch PetscCall(TransformCoordinates(*dm, user));
3569566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
3573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
358c4762a1bSJed Brown }
359c4762a1bSJed Brown
simple_mass(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,PetscReal u_tShift,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar g0[])360d71ae5a4SJacob Faibussowitsch static void simple_mass(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
361d71ae5a4SJacob Faibussowitsch {
362c4762a1bSJed Brown PetscInt d, e;
363ad540459SPierre Jolivet for (d = 0, e = 0; d < dim; d++, e += dim + 1) g0[e] = 1.;
364c4762a1bSJed Brown }
365c4762a1bSJed Brown
366c4762a1bSJed Brown /* < \nabla v, 1/2(\nabla u + {\nabla u}^T) > */
symmetric_gradient_inner_product(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,PetscReal u_tShift,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar C[])367d71ae5a4SJacob Faibussowitsch static void symmetric_gradient_inner_product(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar C[])
368d71ae5a4SJacob Faibussowitsch {
369c4762a1bSJed Brown PetscInt compI, compJ, d, e;
370c4762a1bSJed Brown
371c4762a1bSJed Brown for (compI = 0; compI < dim; ++compI) {
372c4762a1bSJed Brown for (compJ = 0; compJ < dim; ++compJ) {
373c4762a1bSJed Brown for (d = 0; d < dim; ++d) {
374c4762a1bSJed Brown for (e = 0; e < dim; e++) {
375c4762a1bSJed Brown if (d == e && d == compI && d == compJ) {
376c4762a1bSJed Brown C[((compI * dim + compJ) * dim + d) * dim + e] = 1.0;
377c4762a1bSJed Brown } else if ((d == compJ && e == compI) || (d == e && compI == compJ)) {
378c4762a1bSJed Brown C[((compI * dim + compJ) * dim + d) * dim + e] = 0.5;
379c4762a1bSJed Brown } else {
380c4762a1bSJed Brown C[((compI * dim + compJ) * dim + d) * dim + e] = 0.0;
381c4762a1bSJed Brown }
382c4762a1bSJed Brown }
383c4762a1bSJed Brown }
384c4762a1bSJed Brown }
385c4762a1bSJed Brown }
386c4762a1bSJed Brown }
387c4762a1bSJed Brown
SetupSection(DM dm,AppCtx * user)388d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupSection(DM dm, AppCtx *user)
389d71ae5a4SJacob Faibussowitsch {
390c4762a1bSJed Brown PetscFunctionBeginUser;
39130602db0SMatthew G. Knepley if (user->constraints) {
392c4762a1bSJed Brown /* test local constraints */
393c4762a1bSJed Brown DM coordDM;
394c4762a1bSJed Brown PetscInt fStart, fEnd, f, vStart, vEnd, v;
395c4762a1bSJed Brown PetscInt edgesx = 2, vertsx;
396c4762a1bSJed Brown PetscInt edgesy = 2, vertsy;
397c4762a1bSJed Brown PetscMPIInt size;
398c4762a1bSJed Brown PetscInt numConst;
399c4762a1bSJed Brown PetscSection aSec;
400c4762a1bSJed Brown PetscInt *anchors;
401c4762a1bSJed Brown PetscInt offset;
402c4762a1bSJed Brown IS aIS;
403c4762a1bSJed Brown MPI_Comm comm = PetscObjectComm((PetscObject)dm);
404c4762a1bSJed Brown
4059566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size));
40608401ef6SPierre Jolivet PetscCheck(size <= 1, comm, PETSC_ERR_SUP, "Local constraint test can only be performed in serial");
407c4762a1bSJed Brown
408c4762a1bSJed Brown /* we are going to test constraints by using them to enforce periodicity
409c4762a1bSJed Brown * in one direction, and comparing to the existing method of enforcing
410c4762a1bSJed Brown * periodicity */
411c4762a1bSJed Brown
412c4762a1bSJed Brown /* first create the coordinate section so that it does not clone the
413c4762a1bSJed Brown * constraints */
4149566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &coordDM));
415c4762a1bSJed Brown
416c4762a1bSJed Brown /* create the constrained-to-anchor section */
4179566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
4189566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 1, &fStart, &fEnd));
4199566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(PETSC_COMM_SELF, &aSec));
4209566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(aSec, PetscMin(fStart, vStart), PetscMax(fEnd, vEnd)));
421c4762a1bSJed Brown
422c4762a1bSJed Brown /* define the constraints */
4239566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-da_grid_x", &edgesx, NULL));
4249566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-da_grid_y", &edgesy, NULL));
425c4762a1bSJed Brown vertsx = edgesx + 1;
426c4762a1bSJed Brown vertsy = edgesy + 1;
427c4762a1bSJed Brown numConst = vertsy + edgesy;
4289566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numConst, &anchors));
429c4762a1bSJed Brown offset = 0;
430c4762a1bSJed Brown for (v = vStart + edgesx; v < vEnd; v += vertsx) {
4319566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(aSec, v, 1));
432c4762a1bSJed Brown anchors[offset++] = v - edgesx;
433c4762a1bSJed Brown }
434c4762a1bSJed Brown for (f = fStart + edgesx * vertsy + edgesx * edgesy; f < fEnd; f++) {
4359566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(aSec, f, 1));
436c4762a1bSJed Brown anchors[offset++] = f - edgesx * edgesy;
437c4762a1bSJed Brown }
4389566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(aSec));
4399566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, numConst, anchors, PETSC_OWN_POINTER, &aIS));
440c4762a1bSJed Brown
4419566063dSJacob Faibussowitsch PetscCall(DMPlexSetAnchors(dm, aSec, aIS));
4429566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&aSec));
4439566063dSJacob Faibussowitsch PetscCall(ISDestroy(&aIS));
444c4762a1bSJed Brown }
4459566063dSJacob Faibussowitsch PetscCall(DMSetNumFields(dm, 1));
4469566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)user->fe));
4479566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm));
44830602db0SMatthew G. Knepley if (user->constraints) {
449c4762a1bSJed Brown /* test getting local constraint matrix that matches section */
450c4762a1bSJed Brown PetscSection aSec;
451c4762a1bSJed Brown IS aIS;
452c4762a1bSJed Brown
4539566063dSJacob Faibussowitsch PetscCall(DMPlexGetAnchors(dm, &aSec, &aIS));
454c4762a1bSJed Brown if (aSec) {
455c4762a1bSJed Brown PetscDS ds;
456c4762a1bSJed Brown PetscSection cSec, section;
457c4762a1bSJed Brown PetscInt cStart, cEnd, c, numComp;
458c4762a1bSJed Brown Mat cMat, mass;
459c4762a1bSJed Brown Vec local;
460c4762a1bSJed Brown const PetscInt *anchors;
461c4762a1bSJed Brown
4629566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, §ion));
4630b4b7b1cSBarry Smith /* this creates the matrix and preallocates the matrix nonzero structure: we
464c4762a1bSJed Brown * just have to fill in the values */
4659566063dSJacob Faibussowitsch PetscCall(DMGetDefaultConstraints(dm, &cSec, &cMat, NULL));
4669566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(cSec, &cStart, &cEnd));
4679566063dSJacob Faibussowitsch PetscCall(ISGetIndices(aIS, &anchors));
4689566063dSJacob Faibussowitsch PetscCall(PetscFEGetNumComponents(user->fe, &numComp));
469c4762a1bSJed Brown for (c = cStart; c < cEnd; c++) {
470c4762a1bSJed Brown PetscInt cDof;
471c4762a1bSJed Brown
472c4762a1bSJed Brown /* is this point constrained? (does it have an anchor?) */
4739566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(aSec, c, &cDof));
474c4762a1bSJed Brown if (cDof) {
475c4762a1bSJed Brown PetscInt cOff, a, aDof, aOff, j;
47663a3b9bcSJacob Faibussowitsch PetscCheck(cDof == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Found %" PetscInt_FMT " anchor points: should be just one", cDof);
477c4762a1bSJed Brown
478c4762a1bSJed Brown /* find the anchor point */
4799566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(aSec, c, &cOff));
480c4762a1bSJed Brown a = anchors[cOff];
481c4762a1bSJed Brown
482c4762a1bSJed Brown /* find the constrained dofs (row in constraint matrix) */
4839566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(cSec, c, &cDof));
4849566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(cSec, c, &cOff));
485c4762a1bSJed Brown
486c4762a1bSJed Brown /* find the anchor dofs (column in constraint matrix) */
4879566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, a, &aDof));
4889566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, a, &aOff));
489c4762a1bSJed Brown
49063a3b9bcSJacob Faibussowitsch PetscCheck(cDof == aDof, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point and anchor have different number of dofs: %" PetscInt_FMT ", %" PetscInt_FMT, cDof, aDof);
49163a3b9bcSJacob Faibussowitsch PetscCheck(cDof % numComp == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point dofs not divisible by field components: %" PetscInt_FMT ", %" PetscInt_FMT, cDof, numComp);
492c4762a1bSJed Brown
493c4762a1bSJed Brown /* put in a simple equality constraint */
49448a46eb9SPierre Jolivet for (j = 0; j < cDof; j++) PetscCall(MatSetValue(cMat, cOff + j, aOff + j, 1., INSERT_VALUES));
495c4762a1bSJed Brown }
496c4762a1bSJed Brown }
4979566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(cMat, MAT_FINAL_ASSEMBLY));
4989566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(cMat, MAT_FINAL_ASSEMBLY));
4999566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(aIS, &anchors));
500c4762a1bSJed Brown
501c4762a1bSJed Brown /* Now that we have constructed the constraint matrix, any FE matrix
502c4762a1bSJed Brown * that we construct will apply the constraints during construction */
503c4762a1bSJed Brown
5049566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(dm, &mass));
505*2a8381b2SBarry Smith /* get a local variable to serve as the solution */
5069566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &local));
5079566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds));
508c4762a1bSJed Brown /* set the jacobian to be the mass matrix */
5099566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, simple_mass, NULL, NULL, NULL));
510c4762a1bSJed Brown /* build the mass matrix */
5119566063dSJacob Faibussowitsch PetscCall(DMPlexSNESComputeJacobianFEM(dm, local, mass, mass, NULL));
5129566063dSJacob Faibussowitsch PetscCall(MatView(mass, PETSC_VIEWER_STDOUT_WORLD));
5139566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mass));
5149566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &local));
515c4762a1bSJed Brown }
516c4762a1bSJed Brown }
5173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
518c4762a1bSJed Brown }
519c4762a1bSJed Brown
TestFEJacobian(DM dm,AppCtx * user)520d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestFEJacobian(DM dm, AppCtx *user)
521d71ae5a4SJacob Faibussowitsch {
522c4762a1bSJed Brown PetscFunctionBeginUser;
52330602db0SMatthew G. Knepley if (!user->useDA) {
524c4762a1bSJed Brown Vec local;
525c4762a1bSJed Brown const Vec *vecs;
526c4762a1bSJed Brown Mat E;
527c4762a1bSJed Brown MatNullSpace sp;
528c4762a1bSJed Brown PetscBool isNullSpace, hasConst;
52930602db0SMatthew G. Knepley PetscInt dim, n, i;
530c4762a1bSJed Brown Vec res = NULL, localX, localRes;
531c4762a1bSJed Brown PetscDS ds;
532c4762a1bSJed Brown
5339566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim));
53463a3b9bcSJacob Faibussowitsch PetscCheck(user->numComponents == dim, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "The number of components %" PetscInt_FMT " must be equal to the dimension %" PetscInt_FMT " for this test", user->numComponents, dim);
5359566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds));
5369566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, symmetric_gradient_inner_product));
5379566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(dm, &E));
5389566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &local));
5399566063dSJacob Faibussowitsch PetscCall(DMPlexSNESComputeJacobianFEM(dm, local, E, E, NULL));
5409566063dSJacob Faibussowitsch PetscCall(DMPlexCreateRigidBody(dm, 0, &sp));
5419566063dSJacob Faibussowitsch PetscCall(MatNullSpaceGetVecs(sp, &hasConst, &n, &vecs));
5429566063dSJacob Faibussowitsch if (n) PetscCall(VecDuplicate(vecs[0], &res));
5439566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(dm, &localX));
5449566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(dm, &localRes));
545c4762a1bSJed Brown for (i = 0; i < n; i++) { /* also test via matrix-free Jacobian application */
546c4762a1bSJed Brown PetscReal resNorm;
547c4762a1bSJed Brown
5489566063dSJacob Faibussowitsch PetscCall(VecSet(localRes, 0.));
5499566063dSJacob Faibussowitsch PetscCall(VecSet(localX, 0.));
5509566063dSJacob Faibussowitsch PetscCall(VecSet(local, 0.));
5519566063dSJacob Faibussowitsch PetscCall(VecSet(res, 0.));
5529566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm, vecs[i], INSERT_VALUES, localX));
5539566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm, vecs[i], INSERT_VALUES, localX));
5549566063dSJacob Faibussowitsch PetscCall(DMSNESComputeJacobianAction(dm, local, localX, localRes, NULL));
5559566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(dm, localRes, ADD_VALUES, res));
5569566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(dm, localRes, ADD_VALUES, res));
5579566063dSJacob Faibussowitsch PetscCall(VecNorm(res, NORM_2, &resNorm));
55848a46eb9SPierre Jolivet if (resNorm > PETSC_SMALL) PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Symmetric gradient action null space vector %" PetscInt_FMT " residual: %E\n", i, (double)resNorm));
559c4762a1bSJed Brown }
5609566063dSJacob Faibussowitsch PetscCall(VecDestroy(&localRes));
5619566063dSJacob Faibussowitsch PetscCall(VecDestroy(&localX));
5629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&res));
5639566063dSJacob Faibussowitsch PetscCall(MatNullSpaceTest(sp, E, &isNullSpace));
564c4762a1bSJed Brown if (isNullSpace) {
5659566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Symmetric gradient null space: PASS\n"));
566c4762a1bSJed Brown } else {
5679566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Symmetric gradient null space: FAIL\n"));
568c4762a1bSJed Brown }
5699566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&sp));
5709566063dSJacob Faibussowitsch PetscCall(MatDestroy(&E));
5719566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &local));
572c4762a1bSJed Brown }
5733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
574c4762a1bSJed Brown }
575c4762a1bSJed Brown
TestInjector(DM dm,AppCtx * user)576d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestInjector(DM dm, AppCtx *user)
577d71ae5a4SJacob Faibussowitsch {
578c4762a1bSJed Brown DM refTree;
579c4762a1bSJed Brown PetscMPIInt rank;
580c4762a1bSJed Brown
581c4762a1bSJed Brown PetscFunctionBegin;
5829566063dSJacob Faibussowitsch PetscCall(DMPlexGetReferenceTree(dm, &refTree));
583c4762a1bSJed Brown if (refTree) {
584c4762a1bSJed Brown Mat inj;
585c4762a1bSJed Brown
5869566063dSJacob Faibussowitsch PetscCall(DMPlexComputeInjectorReferenceTree(refTree, &inj));
5879566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)inj, "Reference Tree Injector"));
5889566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
58948a46eb9SPierre Jolivet if (rank == 0) PetscCall(MatView(inj, PETSC_VIEWER_STDOUT_SELF));
5909566063dSJacob Faibussowitsch PetscCall(MatDestroy(&inj));
591c4762a1bSJed Brown }
5923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
593c4762a1bSJed Brown }
594c4762a1bSJed Brown
TestFVGrad(DM dm,AppCtx * user)595d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestFVGrad(DM dm, AppCtx *user)
596d71ae5a4SJacob Faibussowitsch {
597c4762a1bSJed Brown MPI_Comm comm;
598c4762a1bSJed Brown DM dmRedist, dmfv, dmgrad, dmCell, refTree;
599c4762a1bSJed Brown PetscFV fv;
60030602db0SMatthew G. Knepley PetscInt dim, nvecs, v, cStart, cEnd, cEndInterior;
601c4762a1bSJed Brown PetscMPIInt size;
602c4762a1bSJed Brown Vec cellgeom, grad, locGrad;
603c4762a1bSJed Brown const PetscScalar *cgeom;
604c4762a1bSJed Brown PetscReal allVecMaxDiff = 0., fvTol = 100. * PETSC_MACHINE_EPSILON;
605c4762a1bSJed Brown
606c4762a1bSJed Brown PetscFunctionBeginUser;
607c4762a1bSJed Brown comm = PetscObjectComm((PetscObject)dm);
6089566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim));
609c4762a1bSJed Brown /* duplicate DM, give dup. a FV discretization */
6109566063dSJacob Faibussowitsch PetscCall(DMSetBasicAdjacency(dm, PETSC_TRUE, PETSC_FALSE));
6119566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size));
612c4762a1bSJed Brown dmRedist = NULL;
61348a46eb9SPierre Jolivet if (size > 1) PetscCall(DMPlexDistributeOverlap(dm, 1, NULL, &dmRedist));
614c4762a1bSJed Brown if (!dmRedist) {
615c4762a1bSJed Brown dmRedist = dm;
6169566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)dmRedist));
617c4762a1bSJed Brown }
6189566063dSJacob Faibussowitsch PetscCall(PetscFVCreate(comm, &fv));
6199566063dSJacob Faibussowitsch PetscCall(PetscFVSetType(fv, PETSCFVLEASTSQUARES));
6209566063dSJacob Faibussowitsch PetscCall(PetscFVSetNumComponents(fv, user->numComponents));
6219566063dSJacob Faibussowitsch PetscCall(PetscFVSetSpatialDimension(fv, dim));
6229566063dSJacob Faibussowitsch PetscCall(PetscFVSetFromOptions(fv));
6239566063dSJacob Faibussowitsch PetscCall(PetscFVSetUp(fv));
6244446c3cdSksagiyam {
6254446c3cdSksagiyam PetscSF pointSF;
6264446c3cdSksagiyam DMLabel label;
6274446c3cdSksagiyam
6284446c3cdSksagiyam PetscCall(DMCreateLabel(dmRedist, "Face Sets"));
6294446c3cdSksagiyam PetscCall(DMGetLabel(dmRedist, "Face Sets", &label));
6304446c3cdSksagiyam PetscCall(DMGetPointSF(dmRedist, &pointSF));
6314446c3cdSksagiyam PetscCall(PetscObjectReference((PetscObject)pointSF));
6324446c3cdSksagiyam PetscCall(DMSetPointSF(dmRedist, NULL));
6334446c3cdSksagiyam PetscCall(DMPlexMarkBoundaryFaces(dmRedist, 1, label));
6344446c3cdSksagiyam PetscCall(DMSetPointSF(dmRedist, pointSF));
6354446c3cdSksagiyam PetscCall(PetscSFDestroy(&pointSF));
6364446c3cdSksagiyam }
6379566063dSJacob Faibussowitsch PetscCall(DMPlexConstructGhostCells(dmRedist, NULL, NULL, &dmfv));
6389566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmRedist));
6399566063dSJacob Faibussowitsch PetscCall(DMSetNumFields(dmfv, 1));
6409566063dSJacob Faibussowitsch PetscCall(DMSetField(dmfv, 0, NULL, (PetscObject)fv));
6419566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dmfv));
6429566063dSJacob Faibussowitsch PetscCall(DMPlexGetReferenceTree(dm, &refTree));
6439566063dSJacob Faibussowitsch if (refTree) PetscCall(DMCopyDisc(dmfv, refTree));
6449566063dSJacob Faibussowitsch PetscCall(DMPlexGetGradientDM(dmfv, fv, &dmgrad));
6459566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dmfv, 0, &cStart, &cEnd));
64630602db0SMatthew G. Knepley nvecs = dim * (dim + 1) / 2;
6479566063dSJacob Faibussowitsch PetscCall(DMPlexGetGeometryFVM(dmfv, NULL, &cellgeom, NULL));
6489566063dSJacob Faibussowitsch PetscCall(VecGetDM(cellgeom, &dmCell));
6499566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(cellgeom, &cgeom));
6509566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmgrad, &grad));
6519566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dmgrad, &locGrad));
6522827ebadSStefano Zampini PetscCall(DMPlexGetCellTypeStratum(dmgrad, DM_POLYTOPE_FV_GHOST, &cEndInterior, NULL));
653c4762a1bSJed Brown cEndInterior = (cEndInterior < 0) ? cEnd : cEndInterior;
654c4762a1bSJed Brown for (v = 0; v < nvecs; v++) {
655c4762a1bSJed Brown Vec locX;
656c4762a1bSJed Brown PetscInt c;
657c4762a1bSJed Brown PetscScalar trueGrad[3][3] = {{0.}};
658c4762a1bSJed Brown const PetscScalar *gradArray;
659c4762a1bSJed Brown PetscReal maxDiff, maxDiffGlob;
660c4762a1bSJed Brown
6619566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dmfv, &locX));
662c4762a1bSJed Brown /* get the local projection of the rigid body mode */
663c4762a1bSJed Brown for (c = cStart; c < cEnd; c++) {
664c4762a1bSJed Brown PetscFVCellGeom *cg;
665c4762a1bSJed Brown PetscScalar cx[3] = {0., 0., 0.};
666c4762a1bSJed Brown
6679566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRead(dmCell, c, cgeom, &cg));
66830602db0SMatthew G. Knepley if (v < dim) {
669c4762a1bSJed Brown cx[v] = 1.;
670c4762a1bSJed Brown } else {
67130602db0SMatthew G. Knepley PetscInt w = v - dim;
672c4762a1bSJed Brown
67330602db0SMatthew G. Knepley cx[(w + 1) % dim] = cg->centroid[(w + 2) % dim];
67430602db0SMatthew G. Knepley cx[(w + 2) % dim] = -cg->centroid[(w + 1) % dim];
675c4762a1bSJed Brown }
6769566063dSJacob Faibussowitsch PetscCall(DMPlexVecSetClosure(dmfv, NULL, locX, c, cx, INSERT_ALL_VALUES));
677c4762a1bSJed Brown }
678c4762a1bSJed Brown /* TODO: this isn't in any header */
6799566063dSJacob Faibussowitsch PetscCall(DMPlexReconstructGradientsFVM(dmfv, locX, grad));
6809566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dmgrad, grad, INSERT_VALUES, locGrad));
6819566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dmgrad, grad, INSERT_VALUES, locGrad));
6829566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(locGrad, &gradArray));
683c4762a1bSJed Brown /* compare computed gradient to exact gradient */
68430602db0SMatthew G. Knepley if (v >= dim) {
68530602db0SMatthew G. Knepley PetscInt w = v - dim;
686c4762a1bSJed Brown
68730602db0SMatthew G. Knepley trueGrad[(w + 1) % dim][(w + 2) % dim] = 1.;
68830602db0SMatthew G. Knepley trueGrad[(w + 2) % dim][(w + 1) % dim] = -1.;
689c4762a1bSJed Brown }
690c4762a1bSJed Brown maxDiff = 0.;
691c4762a1bSJed Brown for (c = cStart; c < cEndInterior; c++) {
692c4762a1bSJed Brown PetscScalar *compGrad;
693c4762a1bSJed Brown PetscInt i, j, k;
694c4762a1bSJed Brown PetscReal FrobDiff = 0.;
695c4762a1bSJed Brown
6969566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRead(dmgrad, c, gradArray, &compGrad));
697c4762a1bSJed Brown
69830602db0SMatthew G. Knepley for (i = 0, k = 0; i < dim; i++) {
69930602db0SMatthew G. Knepley for (j = 0; j < dim; j++, k++) {
700c4762a1bSJed Brown PetscScalar diff = compGrad[k] - trueGrad[i][j];
701c4762a1bSJed Brown FrobDiff += PetscRealPart(diff * PetscConj(diff));
702c4762a1bSJed Brown }
703c4762a1bSJed Brown }
704c4762a1bSJed Brown FrobDiff = PetscSqrtReal(FrobDiff);
705c4762a1bSJed Brown maxDiff = PetscMax(maxDiff, FrobDiff);
706c4762a1bSJed Brown }
707462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&maxDiff, &maxDiffGlob, 1, MPIU_REAL, MPIU_MAX, comm));
708c4762a1bSJed Brown allVecMaxDiff = PetscMax(allVecMaxDiff, maxDiffGlob);
7099566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(locGrad, &gradArray));
7109566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dmfv, &locX));
711c4762a1bSJed Brown }
712c4762a1bSJed Brown if (allVecMaxDiff < fvTol) {
7139566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Finite volume gradient reconstruction: PASS\n"));
714c4762a1bSJed Brown } else {
71563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Finite volume gradient reconstruction: FAIL at tolerance %g with max difference %g\n", (double)fvTol, (double)allVecMaxDiff));
716c4762a1bSJed Brown }
7179566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dmgrad, &locGrad));
7189566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmgrad, &grad));
7199566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(cellgeom, &cgeom));
7209566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmfv));
7219566063dSJacob Faibussowitsch PetscCall(PetscFVDestroy(&fv));
7223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
723c4762a1bSJed Brown }
724c4762a1bSJed Brown
ComputeError(DM dm,PetscErrorCode (** exactFuncs)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),PetscErrorCode (** exactFuncDers)(PetscInt,PetscReal,const PetscReal[],const PetscReal[],PetscInt,PetscScalar *,void *),void ** exactCtxs,PetscReal * error,PetscReal * errorDer,AppCtx * user)725d71ae5a4SJacob Faibussowitsch static PetscErrorCode ComputeError(DM dm, PetscErrorCode (**exactFuncs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), PetscErrorCode (**exactFuncDers)(PetscInt, PetscReal, const PetscReal[], const PetscReal[], PetscInt, PetscScalar *, void *), void **exactCtxs, PetscReal *error, PetscReal *errorDer, AppCtx *user)
726d71ae5a4SJacob Faibussowitsch {
727c4762a1bSJed Brown Vec u;
728c4762a1bSJed Brown PetscReal n[3] = {1.0, 1.0, 1.0};
729c4762a1bSJed Brown
730c4762a1bSJed Brown PetscFunctionBeginUser;
7319566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &u));
732c4762a1bSJed Brown /* Project function into FE function space */
7339566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, 0.0, exactFuncs, exactCtxs, INSERT_ALL_VALUES, u));
7349566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-projection_view"));
735c4762a1bSJed Brown /* Compare approximation to exact in L_2 */
7369566063dSJacob Faibussowitsch PetscCall(DMComputeL2Diff(dm, 0.0, exactFuncs, exactCtxs, u, error));
7379566063dSJacob Faibussowitsch PetscCall(DMComputeL2GradientDiff(dm, 0.0, exactFuncDers, exactCtxs, u, n, errorDer));
7389566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &u));
7393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
740c4762a1bSJed Brown }
741c4762a1bSJed Brown
CheckFunctions(DM dm,PetscInt order,AppCtx * user)742d71ae5a4SJacob Faibussowitsch static PetscErrorCode CheckFunctions(DM dm, PetscInt order, AppCtx *user)
743d71ae5a4SJacob Faibussowitsch {
744*2a8381b2SBarry Smith PetscErrorCode (*exactFuncs[1])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, PetscCtx ctx);
745*2a8381b2SBarry Smith PetscErrorCode (*exactFuncDers[1])(PetscInt dim, PetscReal time, const PetscReal x[], const PetscReal n[], PetscInt Nf, PetscScalar *u, PetscCtx ctx);
746c4762a1bSJed Brown void *exactCtxs[3];
747c4762a1bSJed Brown MPI_Comm comm;
748c4762a1bSJed Brown PetscReal error, errorDer, tol = PETSC_SMALL;
749c4762a1bSJed Brown
750c4762a1bSJed Brown PetscFunctionBeginUser;
751c4762a1bSJed Brown exactCtxs[0] = user;
752c4762a1bSJed Brown exactCtxs[1] = user;
753c4762a1bSJed Brown exactCtxs[2] = user;
7549566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
755c4762a1bSJed Brown /* Setup functions to approximate */
756c4762a1bSJed Brown switch (order) {
757c4762a1bSJed Brown case 0:
7586c2660fcSDarsh Nathawani if (user->RT) {
7596c2660fcSDarsh Nathawani exactFuncs[0] = rt0;
7606c2660fcSDarsh Nathawani exactFuncDers[0] = rt0Der;
7616c2660fcSDarsh Nathawani } else {
762c4762a1bSJed Brown exactFuncs[0] = constant;
763c4762a1bSJed Brown exactFuncDers[0] = constantDer;
7646c2660fcSDarsh Nathawani }
765c4762a1bSJed Brown break;
766c4762a1bSJed Brown case 1:
7676c2660fcSDarsh Nathawani if (user->RT) {
7686c2660fcSDarsh Nathawani exactFuncs[0] = rt1;
7696c2660fcSDarsh Nathawani exactFuncDers[0] = rt1Der;
7706c2660fcSDarsh Nathawani } else {
771c4762a1bSJed Brown exactFuncs[0] = linear;
772c4762a1bSJed Brown exactFuncDers[0] = linearDer;
7736c2660fcSDarsh Nathawani }
774c4762a1bSJed Brown break;
775c4762a1bSJed Brown case 2:
776c4762a1bSJed Brown exactFuncs[0] = quadratic;
777c4762a1bSJed Brown exactFuncDers[0] = quadraticDer;
778c4762a1bSJed Brown break;
779c4762a1bSJed Brown case 3:
780c4762a1bSJed Brown exactFuncs[0] = cubic;
781c4762a1bSJed Brown exactFuncDers[0] = cubicDer;
782c4762a1bSJed Brown break;
783d71ae5a4SJacob Faibussowitsch default:
784d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Could not determine functions to test for order %" PetscInt_FMT, order);
785c4762a1bSJed Brown }
7869566063dSJacob Faibussowitsch PetscCall(ComputeError(dm, exactFuncs, exactFuncDers, exactCtxs, &error, &errorDer, user));
787c4762a1bSJed Brown /* Report result */
78863a3b9bcSJacob Faibussowitsch if (error > tol) PetscCall(PetscPrintf(comm, "Function tests FAIL for order %" PetscInt_FMT " at tolerance %g error %g\n", order, (double)tol, (double)error));
78963a3b9bcSJacob Faibussowitsch else PetscCall(PetscPrintf(comm, "Function tests pass for order %" PetscInt_FMT " at tolerance %g\n", order, (double)tol));
79063a3b9bcSJacob Faibussowitsch if (errorDer > tol) PetscCall(PetscPrintf(comm, "Function tests FAIL for order %" PetscInt_FMT " derivatives at tolerance %g error %g\n", order, (double)tol, (double)errorDer));
79163a3b9bcSJacob Faibussowitsch else PetscCall(PetscPrintf(comm, "Function tests pass for order %" PetscInt_FMT " derivatives at tolerance %g\n", order, (double)tol));
7923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
793c4762a1bSJed Brown }
794c4762a1bSJed Brown
CheckInterpolation(DM dm,PetscBool checkRestrict,PetscInt order,AppCtx * user)795d71ae5a4SJacob Faibussowitsch static PetscErrorCode CheckInterpolation(DM dm, PetscBool checkRestrict, PetscInt order, AppCtx *user)
796d71ae5a4SJacob Faibussowitsch {
797*2a8381b2SBarry Smith PetscErrorCode (*exactFuncs[1])(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, PetscCtx ctx);
798*2a8381b2SBarry Smith PetscErrorCode (*exactFuncDers[1])(PetscInt, PetscReal, const PetscReal x[], const PetscReal n[], PetscInt, PetscScalar *u, PetscCtx ctx);
799c4762a1bSJed Brown PetscReal n[3] = {1.0, 1.0, 1.0};
800c4762a1bSJed Brown void *exactCtxs[3];
801c4762a1bSJed Brown DM rdm, idm, fdm;
802c4762a1bSJed Brown Mat Interp;
803c4762a1bSJed Brown Vec iu, fu, scaling;
804c4762a1bSJed Brown MPI_Comm comm;
80530602db0SMatthew G. Knepley PetscInt dim;
806c4762a1bSJed Brown PetscReal error, errorDer, tol = PETSC_SMALL;
807c4762a1bSJed Brown
808c4762a1bSJed Brown PetscFunctionBeginUser;
809c4762a1bSJed Brown exactCtxs[0] = user;
810c4762a1bSJed Brown exactCtxs[1] = user;
811c4762a1bSJed Brown exactCtxs[2] = user;
8129566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
8139566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim));
8149566063dSJacob Faibussowitsch PetscCall(DMRefine(dm, comm, &rdm));
8159566063dSJacob Faibussowitsch PetscCall(DMSetCoarseDM(rdm, dm));
8169b750a2fSMatthew G. Knepley PetscCall(DMGetCoordinatesLocalSetUp(rdm));
8179566063dSJacob Faibussowitsch PetscCall(DMPlexSetRegularRefinement(rdm, user->convRefine));
81830602db0SMatthew G. Knepley if (user->tree) {
819c4762a1bSJed Brown DM refTree;
8209566063dSJacob Faibussowitsch PetscCall(DMPlexGetReferenceTree(dm, &refTree));
8219566063dSJacob Faibussowitsch PetscCall(DMPlexSetReferenceTree(rdm, refTree));
822c4762a1bSJed Brown }
8239566063dSJacob Faibussowitsch if (user->useDA) PetscCall(DMDASetVertexCoordinates(rdm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
8249566063dSJacob Faibussowitsch PetscCall(SetupSection(rdm, user));
825c4762a1bSJed Brown /* Setup functions to approximate */
826c4762a1bSJed Brown switch (order) {
827c4762a1bSJed Brown case 0:
828c4762a1bSJed Brown exactFuncs[0] = constant;
829c4762a1bSJed Brown exactFuncDers[0] = constantDer;
830c4762a1bSJed Brown break;
831c4762a1bSJed Brown case 1:
832c4762a1bSJed Brown exactFuncs[0] = linear;
833c4762a1bSJed Brown exactFuncDers[0] = linearDer;
834c4762a1bSJed Brown break;
835c4762a1bSJed Brown case 2:
836c4762a1bSJed Brown exactFuncs[0] = quadratic;
837c4762a1bSJed Brown exactFuncDers[0] = quadraticDer;
838c4762a1bSJed Brown break;
839c4762a1bSJed Brown case 3:
840c4762a1bSJed Brown exactFuncs[0] = cubic;
841c4762a1bSJed Brown exactFuncDers[0] = cubicDer;
842c4762a1bSJed Brown break;
843d71ae5a4SJacob Faibussowitsch default:
844d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Could not determine functions to test for dimension %" PetscInt_FMT " order %" PetscInt_FMT, dim, order);
845c4762a1bSJed Brown }
846c4762a1bSJed Brown idm = checkRestrict ? rdm : dm;
847c4762a1bSJed Brown fdm = checkRestrict ? dm : rdm;
8489566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(idm, &iu));
8499566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(fdm, &fu));
8509566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(dm, user));
8519566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(rdm, user));
8529566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(dm, rdm, &Interp, &scaling));
853c4762a1bSJed Brown /* Project function into initial FE function space */
8549566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(idm, 0.0, exactFuncs, exactCtxs, INSERT_ALL_VALUES, iu));
855c4762a1bSJed Brown /* Interpolate function into final FE function space */
8569371c9d4SSatish Balay if (checkRestrict) {
8579371c9d4SSatish Balay PetscCall(MatRestrict(Interp, iu, fu));
8589371c9d4SSatish Balay PetscCall(VecPointwiseMult(fu, scaling, fu));
8599371c9d4SSatish Balay } else PetscCall(MatInterpolate(Interp, iu, fu));
860c4762a1bSJed Brown /* Compare approximation to exact in L_2 */
8619b750a2fSMatthew G. Knepley PetscCall(DMGetCoordinatesLocalSetUp(fdm));
8629566063dSJacob Faibussowitsch PetscCall(DMComputeL2Diff(fdm, 0.0, exactFuncs, exactCtxs, fu, &error));
8639566063dSJacob Faibussowitsch PetscCall(DMComputeL2GradientDiff(fdm, 0.0, exactFuncDers, exactCtxs, fu, n, &errorDer));
864c4762a1bSJed Brown /* Report result */
86563a3b9bcSJacob Faibussowitsch if (error > tol) PetscCall(PetscPrintf(comm, "Interpolation tests FAIL for order %" PetscInt_FMT " at tolerance %g error %g\n", order, (double)tol, (double)error));
86663a3b9bcSJacob Faibussowitsch else PetscCall(PetscPrintf(comm, "Interpolation tests pass for order %" PetscInt_FMT " at tolerance %g\n", order, (double)tol));
86763a3b9bcSJacob Faibussowitsch if (errorDer > tol) PetscCall(PetscPrintf(comm, "Interpolation tests FAIL for order %" PetscInt_FMT " derivatives at tolerance %g error %g\n", order, (double)tol, (double)errorDer));
86863a3b9bcSJacob Faibussowitsch else PetscCall(PetscPrintf(comm, "Interpolation tests pass for order %" PetscInt_FMT " derivatives at tolerance %g\n", order, (double)tol));
8699566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(idm, &iu));
8709566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(fdm, &fu));
8719566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Interp));
8729566063dSJacob Faibussowitsch PetscCall(VecDestroy(&scaling));
8739566063dSJacob Faibussowitsch PetscCall(DMDestroy(&rdm));
8743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
875c4762a1bSJed Brown }
876c4762a1bSJed Brown
CheckConvergence(DM dm,PetscInt Nr,AppCtx * user)877d71ae5a4SJacob Faibussowitsch static PetscErrorCode CheckConvergence(DM dm, PetscInt Nr, AppCtx *user)
878d71ae5a4SJacob Faibussowitsch {
879c4762a1bSJed Brown DM odm = dm, rdm = NULL, cdm = NULL;
880*2a8381b2SBarry Smith PetscErrorCode (*exactFuncs[1])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, PetscCtx ctx) = {trig};
881*2a8381b2SBarry Smith PetscErrorCode (*exactFuncDers[1])(PetscInt dim, PetscReal time, const PetscReal x[], const PetscReal n[], PetscInt Nf, PetscScalar *u, PetscCtx ctx) = {trigDer};
882c4762a1bSJed Brown void *exactCtxs[3];
883c4762a1bSJed Brown PetscInt r, c, cStart, cEnd;
884c4762a1bSJed Brown PetscReal errorOld, errorDerOld, error, errorDer, rel, len, lenOld;
885c4762a1bSJed Brown double p;
886c4762a1bSJed Brown
887c4762a1bSJed Brown PetscFunctionBeginUser;
8883ba16761SJacob Faibussowitsch if (!user->convergence) PetscFunctionReturn(PETSC_SUCCESS);
889c4762a1bSJed Brown exactCtxs[0] = user;
890c4762a1bSJed Brown exactCtxs[1] = user;
891c4762a1bSJed Brown exactCtxs[2] = user;
8929566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)odm));
893c4762a1bSJed Brown if (!user->convRefine) {
894c4762a1bSJed Brown for (r = 0; r < Nr; ++r) {
8959566063dSJacob Faibussowitsch PetscCall(DMRefine(odm, PetscObjectComm((PetscObject)dm), &rdm));
8969566063dSJacob Faibussowitsch PetscCall(DMDestroy(&odm));
897c4762a1bSJed Brown odm = rdm;
898c4762a1bSJed Brown }
8999566063dSJacob Faibussowitsch PetscCall(SetupSection(odm, user));
900c4762a1bSJed Brown }
9019566063dSJacob Faibussowitsch PetscCall(ComputeError(odm, exactFuncs, exactFuncDers, exactCtxs, &errorOld, &errorDerOld, user));
902c4762a1bSJed Brown if (user->convRefine) {
903c4762a1bSJed Brown for (r = 0; r < Nr; ++r) {
9049566063dSJacob Faibussowitsch PetscCall(DMRefine(odm, PetscObjectComm((PetscObject)dm), &rdm));
9059566063dSJacob Faibussowitsch if (user->useDA) PetscCall(DMDASetVertexCoordinates(rdm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
9069566063dSJacob Faibussowitsch PetscCall(SetupSection(rdm, user));
9079566063dSJacob Faibussowitsch PetscCall(ComputeError(rdm, exactFuncs, exactFuncDers, exactCtxs, &error, &errorDer, user));
908c4762a1bSJed Brown p = PetscLog2Real(errorOld / error);
909300f1712SStefano Zampini PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Function convergence rate at refinement %" PetscInt_FMT ": %.2f\n", r, p));
910c4762a1bSJed Brown p = PetscLog2Real(errorDerOld / errorDer);
911300f1712SStefano Zampini PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Derivative convergence rate at refinement %" PetscInt_FMT ": %.2f\n", r, p));
9129566063dSJacob Faibussowitsch PetscCall(DMDestroy(&odm));
913c4762a1bSJed Brown odm = rdm;
914c4762a1bSJed Brown errorOld = error;
915c4762a1bSJed Brown errorDerOld = errorDer;
916c4762a1bSJed Brown }
917c4762a1bSJed Brown } else {
9189566063dSJacob Faibussowitsch /* PetscCall(ComputeLongestEdge(dm, &lenOld)); */
9199566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(odm, 0, &cStart, &cEnd));
920c4762a1bSJed Brown lenOld = cEnd - cStart;
921c4762a1bSJed Brown for (c = 0; c < Nr; ++c) {
9229566063dSJacob Faibussowitsch PetscCall(DMCoarsen(odm, PetscObjectComm((PetscObject)dm), &cdm));
9239566063dSJacob Faibussowitsch if (user->useDA) PetscCall(DMDASetVertexCoordinates(cdm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
9249566063dSJacob Faibussowitsch PetscCall(SetupSection(cdm, user));
9259566063dSJacob Faibussowitsch PetscCall(ComputeError(cdm, exactFuncs, exactFuncDers, exactCtxs, &error, &errorDer, user));
9269566063dSJacob Faibussowitsch /* PetscCall(ComputeLongestEdge(cdm, &len)); */
9279566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(cdm, 0, &cStart, &cEnd));
928c4762a1bSJed Brown len = cEnd - cStart;
929c4762a1bSJed Brown rel = error / errorOld;
930c4762a1bSJed Brown p = PetscLogReal(rel) / PetscLogReal(lenOld / len);
931300f1712SStefano Zampini PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Function convergence rate at coarsening %" PetscInt_FMT ": %.2f\n", c, p));
932c4762a1bSJed Brown rel = errorDer / errorDerOld;
933c4762a1bSJed Brown p = PetscLogReal(rel) / PetscLogReal(lenOld / len);
934300f1712SStefano Zampini PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Derivative convergence rate at coarsening %" PetscInt_FMT ": %.2f\n", c, p));
9359566063dSJacob Faibussowitsch PetscCall(DMDestroy(&odm));
936c4762a1bSJed Brown odm = cdm;
937c4762a1bSJed Brown errorOld = error;
938c4762a1bSJed Brown errorDerOld = errorDer;
939c4762a1bSJed Brown lenOld = len;
940c4762a1bSJed Brown }
941c4762a1bSJed Brown }
9429566063dSJacob Faibussowitsch PetscCall(DMDestroy(&odm));
9433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
944c4762a1bSJed Brown }
945c4762a1bSJed Brown
main(int argc,char ** argv)946d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
947d71ae5a4SJacob Faibussowitsch {
948c4762a1bSJed Brown DM dm;
949c4762a1bSJed Brown AppCtx user; /* user-defined work context */
95030602db0SMatthew G. Knepley PetscInt dim = 2;
95130602db0SMatthew G. Knepley PetscBool simplex = PETSC_FALSE;
952c4762a1bSJed Brown
953327415f7SBarry Smith PetscFunctionBeginUser;
9549566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
9559566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
9569566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
95730602db0SMatthew G. Knepley if (!user.useDA) {
9589566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim));
9599566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex));
96030602db0SMatthew G. Knepley }
9619566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetFromOptions(dm));
96230602db0SMatthew G. Knepley user.numComponents = user.numComponents < 0 ? dim : user.numComponents;
9639566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_WORLD, dim, user.numComponents, simplex, NULL, user.qorder, &user.fe));
9649566063dSJacob Faibussowitsch PetscCall(SetupSection(dm, &user));
9659566063dSJacob Faibussowitsch if (user.testFEjacobian) PetscCall(TestFEJacobian(dm, &user));
9669566063dSJacob Faibussowitsch if (user.testFVgrad) PetscCall(TestFVGrad(dm, &user));
9679566063dSJacob Faibussowitsch if (user.testInjector) PetscCall(TestInjector(dm, &user));
9689566063dSJacob Faibussowitsch PetscCall(CheckFunctions(dm, user.porder, &user));
969c4762a1bSJed Brown {
970c4762a1bSJed Brown PetscDualSpace dsp;
971c4762a1bSJed Brown PetscInt k;
972c4762a1bSJed Brown
9739566063dSJacob Faibussowitsch PetscCall(PetscFEGetDualSpace(user.fe, &dsp));
9749566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDeRahm(dsp, &k));
975e3154614SMatthew G. Knepley if (dim == 2 && user.constraints == PETSC_FALSE && user.tree == PETSC_FALSE && k == 0) {
9769566063dSJacob Faibussowitsch PetscCall(CheckInterpolation(dm, PETSC_FALSE, user.porder, &user));
9779566063dSJacob Faibussowitsch PetscCall(CheckInterpolation(dm, PETSC_TRUE, user.porder, &user));
978c4762a1bSJed Brown }
979c4762a1bSJed Brown }
9809566063dSJacob Faibussowitsch PetscCall(CheckConvergence(dm, 3, &user));
9819566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&user.fe));
9829566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm));
9839566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
984b122ec5aSJacob Faibussowitsch return 0;
985c4762a1bSJed Brown }
986c4762a1bSJed Brown
987c4762a1bSJed Brown /*TEST
988c4762a1bSJed Brown
989c4762a1bSJed Brown test:
990c4762a1bSJed Brown suffix: 1
991c4762a1bSJed Brown requires: triangle
992c4762a1bSJed Brown
993e3154614SMatthew G. Knepley # 2D P_0 on a triangle
994e3154614SMatthew G. Knepley test:
995e3154614SMatthew G. Knepley suffix: p0_2d_0
996e3154614SMatthew G. Knepley requires: triangle
997e3154614SMatthew G. Knepley args: -petscspace_degree 1 -qorder 1 -convergence
998e3154614SMatthew G. Knepley test:
999e3154614SMatthew G. Knepley suffix: p0_2d_1
1000e3154614SMatthew G. Knepley requires: triangle
1001e3154614SMatthew G. Knepley args: -petscspace_degree 1 -qorder 1 -porder 1
1002e3154614SMatthew G. Knepley
1003c4762a1bSJed Brown # 2D P_1 on a triangle
1004c4762a1bSJed Brown test:
1005c4762a1bSJed Brown suffix: p1_2d_0
1006c4762a1bSJed Brown requires: triangle
1007c4762a1bSJed Brown args: -petscspace_degree 1 -qorder 1 -convergence
1008c4762a1bSJed Brown test:
1009c4762a1bSJed Brown suffix: p1_2d_1
1010c4762a1bSJed Brown requires: triangle
1011c4762a1bSJed Brown args: -petscspace_degree 1 -qorder 1 -porder 1
1012c4762a1bSJed Brown test:
1013c4762a1bSJed Brown suffix: p1_2d_2
1014c4762a1bSJed Brown requires: triangle
1015c4762a1bSJed Brown args: -petscspace_degree 1 -qorder 1 -porder 2
1016c4762a1bSJed Brown test:
1017c4762a1bSJed Brown suffix: p1_2d_3
1018e788613bSJoe Wallwork requires: triangle mmg
1019c4762a1bSJed Brown args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -convergence -conv_refine 0
1020c4762a1bSJed Brown test:
1021c4762a1bSJed Brown suffix: p1_2d_4
1022e788613bSJoe Wallwork requires: triangle mmg
1023c4762a1bSJed Brown args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 1 -conv_refine 0
1024c4762a1bSJed Brown test:
1025c4762a1bSJed Brown suffix: p1_2d_5
1026e788613bSJoe Wallwork requires: triangle mmg
1027c4762a1bSJed Brown args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 2 -conv_refine 0
1028c4762a1bSJed Brown
1029c4762a1bSJed Brown # 3D P_1 on a tetrahedron
1030c4762a1bSJed Brown test:
1031c4762a1bSJed Brown suffix: p1_3d_0
1032c4762a1bSJed Brown requires: ctetgen
103330602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -convergence
1034c4762a1bSJed Brown test:
1035c4762a1bSJed Brown suffix: p1_3d_1
1036c4762a1bSJed Brown requires: ctetgen
103730602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -porder 1
1038c4762a1bSJed Brown test:
1039c4762a1bSJed Brown suffix: p1_3d_2
1040c4762a1bSJed Brown requires: ctetgen
104130602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -porder 2
1042c4762a1bSJed Brown test:
1043c4762a1bSJed Brown suffix: p1_3d_3
1044e788613bSJoe Wallwork requires: ctetgen mmg
104530602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -convergence -conv_refine 0
1046c4762a1bSJed Brown test:
1047c4762a1bSJed Brown suffix: p1_3d_4
1048e788613bSJoe Wallwork requires: ctetgen mmg
104930602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 1 -conv_refine 0
1050c4762a1bSJed Brown test:
1051c4762a1bSJed Brown suffix: p1_3d_5
1052e788613bSJoe Wallwork requires: ctetgen mmg
105330602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 2 -conv_refine 0
1054c4762a1bSJed Brown
1055c4762a1bSJed Brown # 2D P_2 on a triangle
1056c4762a1bSJed Brown test:
1057c4762a1bSJed Brown suffix: p2_2d_0
1058c4762a1bSJed Brown requires: triangle
1059c4762a1bSJed Brown args: -petscspace_degree 2 -qorder 2 -convergence
1060c4762a1bSJed Brown test:
1061c4762a1bSJed Brown suffix: p2_2d_1
1062c4762a1bSJed Brown requires: triangle
1063c4762a1bSJed Brown args: -petscspace_degree 2 -qorder 2 -porder 1
1064c4762a1bSJed Brown test:
1065c4762a1bSJed Brown suffix: p2_2d_2
1066c4762a1bSJed Brown requires: triangle
1067c4762a1bSJed Brown args: -petscspace_degree 2 -qorder 2 -porder 2
1068c4762a1bSJed Brown test:
1069c4762a1bSJed Brown suffix: p2_2d_3
1070e788613bSJoe Wallwork requires: triangle mmg
1071c4762a1bSJed Brown args: -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -convergence -conv_refine 0
1072c4762a1bSJed Brown test:
1073c4762a1bSJed Brown suffix: p2_2d_4
1074e788613bSJoe Wallwork requires: triangle mmg
1075c4762a1bSJed Brown args: -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder 1 -conv_refine 0
1076c4762a1bSJed Brown test:
1077c4762a1bSJed Brown suffix: p2_2d_5
1078e788613bSJoe Wallwork requires: triangle mmg
1079c4762a1bSJed Brown args: -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder 2 -conv_refine 0
1080c4762a1bSJed Brown
1081c4762a1bSJed Brown # 3D P_2 on a tetrahedron
1082c4762a1bSJed Brown test:
1083c4762a1bSJed Brown suffix: p2_3d_0
1084c4762a1bSJed Brown requires: ctetgen
108530602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -convergence
1086c4762a1bSJed Brown test:
1087c4762a1bSJed Brown suffix: p2_3d_1
1088c4762a1bSJed Brown requires: ctetgen
108930602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -porder 1
1090c4762a1bSJed Brown test:
1091c4762a1bSJed Brown suffix: p2_3d_2
1092c4762a1bSJed Brown requires: ctetgen
109330602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -porder 2
1094c4762a1bSJed Brown test:
1095c4762a1bSJed Brown suffix: p2_3d_3
1096e788613bSJoe Wallwork requires: ctetgen mmg
109730602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -convergence -conv_refine 0
1098c4762a1bSJed Brown test:
1099c4762a1bSJed Brown suffix: p2_3d_4
1100e788613bSJoe Wallwork requires: ctetgen mmg
110130602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder 1 -conv_refine 0
1102c4762a1bSJed Brown test:
1103c4762a1bSJed Brown suffix: p2_3d_5
1104e788613bSJoe Wallwork requires: ctetgen mmg
110530602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder 2 -conv_refine 0
1106c4762a1bSJed Brown
1107c4762a1bSJed Brown # 2D Q_1 on a quadrilaterial DA
1108c4762a1bSJed Brown test:
1109c4762a1bSJed Brown suffix: q1_2d_da_0
11100338c944SBarry Smith TODO: broken
111130602db0SMatthew G. Knepley args: -use_da 1 -petscspace_degree 1 -qorder 1 -convergence
1112c4762a1bSJed Brown test:
1113c4762a1bSJed Brown suffix: q1_2d_da_1
11140338c944SBarry Smith TODO: broken
111530602db0SMatthew G. Knepley args: -use_da 1 -petscspace_degree 1 -qorder 1 -porder 1
1116c4762a1bSJed Brown test:
1117c4762a1bSJed Brown suffix: q1_2d_da_2
11180338c944SBarry Smith TODO: broken
111930602db0SMatthew G. Knepley args: -use_da 1 -petscspace_degree 1 -qorder 1 -porder 2
1120c4762a1bSJed Brown
1121e3154614SMatthew G. Knepley # 2D P_0 on a quadrilaterial Plex
1122e3154614SMatthew G. Knepley test:
1123e3154614SMatthew G. Knepley suffix: p0_2d_plex_0
1124e3154614SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 0 -qorder 1 -convergence
1125e3154614SMatthew G. Knepley test:
1126e3154614SMatthew G. Knepley suffix: p0_2d_plex_1
1127e3154614SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 0 -qorder 1 -porder 1
1128e3154614SMatthew G. Knepley
1129c4762a1bSJed Brown # 2D Q_1 on a quadrilaterial Plex
1130c4762a1bSJed Brown test:
1131c4762a1bSJed Brown suffix: q1_2d_plex_0
113230602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -convergence
1133c4762a1bSJed Brown test:
1134c4762a1bSJed Brown suffix: q1_2d_plex_1
113530602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -porder 1
1136c4762a1bSJed Brown test:
1137c4762a1bSJed Brown suffix: q1_2d_plex_2
113830602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -porder 2
1139c4762a1bSJed Brown test:
1140c4762a1bSJed Brown suffix: q1_2d_plex_3
114130602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -porder 1 -shear_coords
1142c4762a1bSJed Brown test:
1143c4762a1bSJed Brown suffix: q1_2d_plex_4
114430602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -porder 2 -shear_coords
1145c4762a1bSJed Brown test:
1146c4762a1bSJed Brown suffix: q1_2d_plex_5
114730602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -petscspace_type tensor -qorder 1 -porder 0 -non_affine_coords -convergence
1148c4762a1bSJed Brown test:
1149c4762a1bSJed Brown suffix: q1_2d_plex_6
115030602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -petscspace_type tensor -qorder 1 -porder 1 -non_affine_coords -convergence
1151c4762a1bSJed Brown test:
1152c4762a1bSJed Brown suffix: q1_2d_plex_7
115330602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -petscspace_type tensor -qorder 1 -porder 2 -non_affine_coords -convergence
11549b750a2fSMatthew G. Knepley test:
11559b750a2fSMatthew G. Knepley suffix: q1_2d_plex_8
11569b750a2fSMatthew G. Knepley requires: triangle
11579b750a2fSMatthew G. Knepley args: -dist_dm_refine 1 -dist_dm_plex_transform_type refine_tobox -petscspace_degree 1 -qorder 1 -convergence
1158c4762a1bSJed Brown
1159c4762a1bSJed Brown # 2D Q_2 on a quadrilaterial
1160e3154614SMatthew G. Knepley # The derivative interpolation test fails in single because we lose precision
1161c4762a1bSJed Brown test:
1162c4762a1bSJed Brown suffix: q2_2d_plex_0
1163e3154614SMatthew G. Knepley requires: !single
116430602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -convergence
1165c4762a1bSJed Brown test:
1166c4762a1bSJed Brown suffix: q2_2d_plex_1
116730602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 1
1168c4762a1bSJed Brown test:
1169c4762a1bSJed Brown suffix: q2_2d_plex_2
117030602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 2
1171c4762a1bSJed Brown test:
1172c4762a1bSJed Brown suffix: q2_2d_plex_3
117330602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 1 -shear_coords
1174c4762a1bSJed Brown test:
1175c4762a1bSJed Brown suffix: q2_2d_plex_4
117630602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 2 -shear_coords
1177e3154614SMatthew G. Knepley # The derivative interpolation test fails in single because we lose precision
1178c4762a1bSJed Brown test:
1179c4762a1bSJed Brown suffix: q2_2d_plex_5
1180e3154614SMatthew G. Knepley requires: !single
118130602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_type tensor -qorder 2 -porder 0 -non_affine_coords -convergence
1182e3154614SMatthew G. Knepley # The derivative interpolation test fails in single because we lose precision
1183c4762a1bSJed Brown test:
1184c4762a1bSJed Brown suffix: q2_2d_plex_6
1185e3154614SMatthew G. Knepley requires: !single
118630602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_type tensor -qorder 2 -porder 1 -non_affine_coords -convergence
1187c4762a1bSJed Brown test:
1188c4762a1bSJed Brown suffix: q2_2d_plex_7
118930602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_type tensor -qorder 2 -porder 2 -non_affine_coords -convergence
1190c4762a1bSJed Brown
1191c4762a1bSJed Brown # 2D P_3 on a triangle
1192c4762a1bSJed Brown test:
1193c4762a1bSJed Brown suffix: p3_2d_0
1194c4762a1bSJed Brown requires: triangle !single
1195c4762a1bSJed Brown args: -petscspace_degree 3 -qorder 3 -convergence
1196c4762a1bSJed Brown test:
1197c4762a1bSJed Brown suffix: p3_2d_1
1198c4762a1bSJed Brown requires: triangle !single
1199c4762a1bSJed Brown args: -petscspace_degree 3 -qorder 3 -porder 1
1200c4762a1bSJed Brown test:
1201c4762a1bSJed Brown suffix: p3_2d_2
1202c4762a1bSJed Brown requires: triangle !single
1203c4762a1bSJed Brown args: -petscspace_degree 3 -qorder 3 -porder 2
1204c4762a1bSJed Brown test:
1205c4762a1bSJed Brown suffix: p3_2d_3
1206c4762a1bSJed Brown requires: triangle !single
1207c4762a1bSJed Brown args: -petscspace_degree 3 -qorder 3 -porder 3
1208c4762a1bSJed Brown test:
1209c4762a1bSJed Brown suffix: p3_2d_4
1210e788613bSJoe Wallwork requires: triangle mmg
1211c4762a1bSJed Brown args: -petscspace_degree 3 -qorder 3 -dm_plex_hash_location -convergence -conv_refine 0
1212c4762a1bSJed Brown test:
1213c4762a1bSJed Brown suffix: p3_2d_5
1214e788613bSJoe Wallwork requires: triangle mmg
1215c4762a1bSJed Brown args: -petscspace_degree 3 -qorder 3 -dm_plex_hash_location -porder 1 -conv_refine 0
1216c4762a1bSJed Brown test:
1217c4762a1bSJed Brown suffix: p3_2d_6
1218e788613bSJoe Wallwork requires: triangle mmg
1219c4762a1bSJed Brown args: -petscspace_degree 3 -qorder 3 -dm_plex_hash_location -porder 3 -conv_refine 0
1220c4762a1bSJed Brown
1221c4762a1bSJed Brown # 2D Q_3 on a quadrilaterial
1222c4762a1bSJed Brown test:
1223c4762a1bSJed Brown suffix: q3_2d_0
122499a192c5SJunchao Zhang requires: !single
122530602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 3 -qorder 3 -convergence
1226c4762a1bSJed Brown test:
1227c4762a1bSJed Brown suffix: q3_2d_1
122899a192c5SJunchao Zhang requires: !single
122930602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 3 -qorder 3 -porder 1
1230c4762a1bSJed Brown test:
1231c4762a1bSJed Brown suffix: q3_2d_2
123299a192c5SJunchao Zhang requires: !single
123330602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 3 -qorder 3 -porder 2
1234c4762a1bSJed Brown test:
1235c4762a1bSJed Brown suffix: q3_2d_3
123699a192c5SJunchao Zhang requires: !single
123730602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 3 -qorder 3 -porder 3
1238c4762a1bSJed Brown
1239c4762a1bSJed Brown # 2D P_1disc on a triangle/quadrilateral
1240c4762a1bSJed Brown test:
1241c4762a1bSJed Brown suffix: p1d_2d_0
1242c4762a1bSJed Brown requires: triangle
1243c4762a1bSJed Brown args: -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -convergence
1244c4762a1bSJed Brown test:
1245c4762a1bSJed Brown suffix: p1d_2d_1
1246c4762a1bSJed Brown requires: triangle
1247c4762a1bSJed Brown args: -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder 1
1248c4762a1bSJed Brown test:
1249c4762a1bSJed Brown suffix: p1d_2d_2
1250c4762a1bSJed Brown requires: triangle
1251c4762a1bSJed Brown args: -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder 2
1252c4762a1bSJed Brown test:
1253c4762a1bSJed Brown suffix: p1d_2d_3
1254c4762a1bSJed Brown requires: triangle
125530602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -convergence
1256c4762a1bSJed Brown filter: sed -e "s/convergence rate at refinement 0: 2/convergence rate at refinement 0: 1.9/g"
1257c4762a1bSJed Brown test:
1258c4762a1bSJed Brown suffix: p1d_2d_4
1259c4762a1bSJed Brown requires: triangle
126030602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder 1
1261c4762a1bSJed Brown test:
1262c4762a1bSJed Brown suffix: p1d_2d_5
1263c4762a1bSJed Brown requires: triangle
126430602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder 2
1265c4762a1bSJed Brown
1266c4762a1bSJed Brown # 2D BDM_1 on a triangle
1267c4762a1bSJed Brown test:
1268c4762a1bSJed Brown suffix: bdm1_2d_0
1269c4762a1bSJed Brown requires: triangle
1270c4762a1bSJed Brown args: -petscspace_degree 1 -petscdualspace_type bdm \
1271c4762a1bSJed Brown -num_comp 2 -qorder 1 -convergence
1272c4762a1bSJed Brown test:
1273c4762a1bSJed Brown suffix: bdm1_2d_1
1274c4762a1bSJed Brown requires: triangle
1275c4762a1bSJed Brown args: -petscspace_degree 1 -petscdualspace_type bdm \
1276c4762a1bSJed Brown -num_comp 2 -qorder 1 -porder 1
1277c4762a1bSJed Brown test:
1278c4762a1bSJed Brown suffix: bdm1_2d_2
1279c4762a1bSJed Brown requires: triangle
1280c4762a1bSJed Brown args: -petscspace_degree 1 -petscdualspace_type bdm \
1281c4762a1bSJed Brown -num_comp 2 -qorder 1 -porder 2
1282c4762a1bSJed Brown
1283c4762a1bSJed Brown # 2D BDM_1 on a quadrilateral
1284c4762a1bSJed Brown test:
1285c4762a1bSJed Brown suffix: bdm1q_2d_0
1286c4762a1bSJed Brown requires: triangle
1287c4762a1bSJed Brown args: -petscspace_degree 1 -petscdualspace_type bdm \
12883f27d899SToby Isaac -petscdualspace_lagrange_tensor 1 \
128930602db0SMatthew G. Knepley -dm_plex_simplex 0 -num_comp 2 -qorder 1 -convergence
1290c4762a1bSJed Brown test:
1291c4762a1bSJed Brown suffix: bdm1q_2d_1
1292c4762a1bSJed Brown requires: triangle
1293c4762a1bSJed Brown args: -petscspace_degree 1 -petscdualspace_type bdm \
12943f27d899SToby Isaac -petscdualspace_lagrange_tensor 1 \
129530602db0SMatthew G. Knepley -dm_plex_simplex 0 -num_comp 2 -qorder 1 -porder 1
1296c4762a1bSJed Brown test:
1297c4762a1bSJed Brown suffix: bdm1q_2d_2
1298c4762a1bSJed Brown requires: triangle
1299c4762a1bSJed Brown args: -petscspace_degree 1 -petscdualspace_type bdm \
13003f27d899SToby Isaac -petscdualspace_lagrange_tensor 1 \
130130602db0SMatthew G. Knepley -dm_plex_simplex 0 -num_comp 2 -qorder 1 -porder 2
1302c4762a1bSJed Brown
1303c4762a1bSJed Brown # Test high order quadrature
1304c4762a1bSJed Brown test:
1305c4762a1bSJed Brown suffix: p1_quad_2
1306c4762a1bSJed Brown requires: triangle
1307c4762a1bSJed Brown args: -petscspace_degree 1 -qorder 2 -porder 1
1308c4762a1bSJed Brown test:
1309c4762a1bSJed Brown suffix: p1_quad_5
1310c4762a1bSJed Brown requires: triangle
1311c4762a1bSJed Brown args: -petscspace_degree 1 -qorder 5 -porder 1
1312c4762a1bSJed Brown test:
1313c4762a1bSJed Brown suffix: p2_quad_3
1314c4762a1bSJed Brown requires: triangle
1315c4762a1bSJed Brown args: -petscspace_degree 2 -qorder 3 -porder 2
1316c4762a1bSJed Brown test:
1317c4762a1bSJed Brown suffix: p2_quad_5
1318c4762a1bSJed Brown requires: triangle
1319c4762a1bSJed Brown args: -petscspace_degree 2 -qorder 5 -porder 2
1320c4762a1bSJed Brown test:
1321c4762a1bSJed Brown suffix: q1_quad_2
132230602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 2 -porder 1
1323c4762a1bSJed Brown test:
1324c4762a1bSJed Brown suffix: q1_quad_5
132530602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 5 -porder 1
1326c4762a1bSJed Brown test:
1327c4762a1bSJed Brown suffix: q2_quad_3
132830602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 3 -porder 1
1329c4762a1bSJed Brown test:
1330c4762a1bSJed Brown suffix: q2_quad_5
133130602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 5 -porder 1
1332c4762a1bSJed Brown
1333c4762a1bSJed Brown # Nonconforming tests
1334c4762a1bSJed Brown test:
1335c4762a1bSJed Brown suffix: constraints
133630602db0SMatthew G. Knepley args: -dm_coord_space 0 -dm_plex_simplex 0 -petscspace_type tensor -petscspace_degree 1 -qorder 0 -constraints
1337c4762a1bSJed Brown test:
1338c4762a1bSJed Brown suffix: nonconforming_tensor_2
1339c4762a1bSJed Brown nsize: 4
134030602db0SMatthew G. Knepley args: -dist_dm_distribute -test_fe_jacobian -test_injector -petscpartitioner_type simple -tree -dm_plex_simplex 0 -dm_plex_max_projection_height 1 -petscspace_type tensor -petscspace_degree 2 -qorder 2 -dm_view ascii::ASCII_INFO_DETAIL
1341c4762a1bSJed Brown test:
1342c4762a1bSJed Brown suffix: nonconforming_tensor_3
1343c4762a1bSJed Brown nsize: 4
134430602db0SMatthew G. Knepley args: -dist_dm_distribute -test_fe_jacobian -petscpartitioner_type simple -tree -dm_plex_simplex 0 -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -dm_plex_max_projection_height 2 -petscspace_type tensor -petscspace_degree 1 -qorder 1 -dm_view ascii::ASCII_INFO_DETAIL
1345c4762a1bSJed Brown test:
1346c4762a1bSJed Brown suffix: nonconforming_tensor_2_fv
1347c4762a1bSJed Brown nsize: 4
134830602db0SMatthew G. Knepley args: -dist_dm_distribute -test_fv_grad -test_injector -petsclimiter_type none -petscpartitioner_type simple -tree -dm_plex_simplex 0 -num_comp 2
1349c4762a1bSJed Brown test:
1350c4762a1bSJed Brown suffix: nonconforming_tensor_3_fv
1351c4762a1bSJed Brown nsize: 4
135230602db0SMatthew G. Knepley args: -dist_dm_distribute -test_fv_grad -test_injector -petsclimiter_type none -petscpartitioner_type simple -tree -dm_plex_simplex 0 -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -num_comp 3
1353c4762a1bSJed Brown test:
1354c4762a1bSJed Brown suffix: nonconforming_tensor_2_hi
1355c4762a1bSJed Brown requires: !single
1356c4762a1bSJed Brown nsize: 4
135730602db0SMatthew G. Knepley args: -dist_dm_distribute -test_fe_jacobian -petscpartitioner_type simple -tree -dm_plex_simplex 0 -dm_plex_max_projection_height 1 -petscspace_type tensor -petscspace_degree 4 -qorder 4
1358c4762a1bSJed Brown test:
1359c4762a1bSJed Brown suffix: nonconforming_tensor_3_hi
1360c4762a1bSJed Brown requires: !single skip
1361c4762a1bSJed Brown nsize: 4
136230602db0SMatthew G. Knepley args: -dist_dm_distribute -test_fe_jacobian -petscpartitioner_type simple -tree -dm_plex_simplex 0 -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -dm_plex_max_projection_height 2 -petscspace_type tensor -petscspace_degree 4 -qorder 4
1363c4762a1bSJed Brown test:
1364c4762a1bSJed Brown suffix: nonconforming_simplex_2
1365c4762a1bSJed Brown requires: triangle
1366c4762a1bSJed Brown nsize: 4
136730602db0SMatthew G. Knepley args: -dist_dm_distribute -test_fe_jacobian -test_injector -petscpartitioner_type simple -tree -dm_plex_max_projection_height 1 -petscspace_degree 2 -qorder 2 -dm_view ascii::ASCII_INFO_DETAIL
1368c4762a1bSJed Brown test:
1369c4762a1bSJed Brown suffix: nonconforming_simplex_2_hi
1370c4762a1bSJed Brown requires: triangle !single
1371c4762a1bSJed Brown nsize: 4
137230602db0SMatthew G. Knepley args: -dist_dm_distribute -test_fe_jacobian -petscpartitioner_type simple -tree -dm_plex_max_projection_height 1 -petscspace_degree 4 -qorder 4
1373c4762a1bSJed Brown test:
1374c4762a1bSJed Brown suffix: nonconforming_simplex_2_fv
1375c4762a1bSJed Brown requires: triangle
1376c4762a1bSJed Brown nsize: 4
137730602db0SMatthew G. Knepley args: -dist_dm_distribute -test_fv_grad -test_injector -petsclimiter_type none -petscpartitioner_type simple -tree -num_comp 2
1378c4762a1bSJed Brown test:
1379c4762a1bSJed Brown suffix: nonconforming_simplex_3
1380c4762a1bSJed Brown requires: ctetgen
1381c4762a1bSJed Brown nsize: 4
138230602db0SMatthew G. Knepley args: -dist_dm_distribute -test_fe_jacobian -test_injector -petscpartitioner_type simple -tree -dm_plex_dim 3 -dm_plex_max_projection_height 2 -petscspace_degree 2 -qorder 2 -dm_view ascii::ASCII_INFO_DETAIL
1383c4762a1bSJed Brown test:
1384c4762a1bSJed Brown suffix: nonconforming_simplex_3_hi
1385c4762a1bSJed Brown requires: ctetgen skip
1386c4762a1bSJed Brown nsize: 4
138730602db0SMatthew G. Knepley args: -dist_dm_distribute -test_fe_jacobian -petscpartitioner_type simple -tree -dm_plex_dim 3 -dm_plex_max_projection_height 2 -petscspace_degree 4 -qorder 4
1388c4762a1bSJed Brown test:
1389c4762a1bSJed Brown suffix: nonconforming_simplex_3_fv
1390c4762a1bSJed Brown requires: ctetgen
1391c4762a1bSJed Brown nsize: 4
139230602db0SMatthew G. Knepley args: -dist_dm_distribute -test_fv_grad -test_injector -petsclimiter_type none -petscpartitioner_type simple -tree -dm_plex_dim 3 -num_comp 3
1393c4762a1bSJed Brown
1394d21efd2eSMatthew G. Knepley # 3D WXY on a triangular prism
1395d21efd2eSMatthew G. Knepley test:
1396d21efd2eSMatthew G. Knepley suffix: wxy_0
1397d21efd2eSMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism -qorder 2 -porder 0 \
1398e239af90SMatthew G. Knepley -petscspace_type sum \
1399e239af90SMatthew G. Knepley -petscspace_variables 3 \
1400e239af90SMatthew G. Knepley -petscspace_components 3 \
1401e239af90SMatthew G. Knepley -petscspace_sum_spaces 2 \
1402e239af90SMatthew G. Knepley -petscspace_sum_concatenate false \
1403e239af90SMatthew G. Knepley -sumcomp_0_petscspace_variables 3 \
1404e239af90SMatthew G. Knepley -sumcomp_0_petscspace_components 3 \
1405e239af90SMatthew G. Knepley -sumcomp_0_petscspace_degree 1 \
1406e239af90SMatthew G. Knepley -sumcomp_1_petscspace_variables 3 \
1407e239af90SMatthew G. Knepley -sumcomp_1_petscspace_components 3 \
1408e239af90SMatthew G. Knepley -sumcomp_1_petscspace_type wxy \
1409e239af90SMatthew G. Knepley -petscdualspace_refcell triangular_prism \
1410e239af90SMatthew G. Knepley -petscdualspace_form_degree 0 \
1411e239af90SMatthew G. Knepley -petscdualspace_order 1 \
1412e239af90SMatthew G. Knepley -petscdualspace_components 3
1413d21efd2eSMatthew G. Knepley
14146c2660fcSDarsh Nathawani # 2D RT_0 on a triangle
14156c2660fcSDarsh Nathawani test:
14166c2660fcSDarsh Nathawani suffix: rt0_2d_tri
14176c2660fcSDarsh Nathawani requires: triangle
14186c2660fcSDarsh Nathawani args: -qorder 1 -porder 0 -RT \
14196c2660fcSDarsh Nathawani -petscspace_type ptrimmed \
14206c2660fcSDarsh Nathawani -petscspace_components 2 \
14216c2660fcSDarsh Nathawani -petscspace_ptrimmed_form_degree -1 \
14226c2660fcSDarsh Nathawani -petscdualspace_order 1 \
14236c2660fcSDarsh Nathawani -petscdualspace_form_degree -1 \
14246c2660fcSDarsh Nathawani -petscdualspace_lagrange_trimmed true
14256c2660fcSDarsh Nathawani
14266c2660fcSDarsh Nathawani # 2D RT_0 on a quadrilateral
14276c2660fcSDarsh Nathawani test:
14286c2660fcSDarsh Nathawani suffix: rt0_2d_quad
14296c2660fcSDarsh Nathawani requires: triangle
14306c2660fcSDarsh Nathawani args: -dm_plex_simplex 0 -qorder 1 -porder 0 -RT \
14316c2660fcSDarsh Nathawani -petscspace_degree 1 \
14326c2660fcSDarsh Nathawani -petscspace_type sum \
14336c2660fcSDarsh Nathawani -petscspace_variables 2 \
14346c2660fcSDarsh Nathawani -petscspace_components 2 \
14356c2660fcSDarsh Nathawani -petscspace_sum_spaces 2 \
14366c2660fcSDarsh Nathawani -petscspace_sum_concatenate true \
14376c2660fcSDarsh Nathawani -sumcomp_0_petscspace_variables 2 \
14386c2660fcSDarsh Nathawani -sumcomp_0_petscspace_type tensor \
14396c2660fcSDarsh Nathawani -sumcomp_0_petscspace_tensor_spaces 2 \
14406c2660fcSDarsh Nathawani -sumcomp_0_petscspace_tensor_uniform false \
14416c2660fcSDarsh Nathawani -sumcomp_0_tensorcomp_0_petscspace_degree 1 \
14426c2660fcSDarsh Nathawani -sumcomp_0_tensorcomp_1_petscspace_degree 0 \
14436c2660fcSDarsh Nathawani -sumcomp_1_petscspace_variables 2 \
14446c2660fcSDarsh Nathawani -sumcomp_1_petscspace_type tensor \
14456c2660fcSDarsh Nathawani -sumcomp_1_petscspace_tensor_spaces 2 \
14466c2660fcSDarsh Nathawani -sumcomp_1_petscspace_tensor_uniform false \
14476c2660fcSDarsh Nathawani -sumcomp_1_tensorcomp_0_petscspace_degree 0 \
14486c2660fcSDarsh Nathawani -sumcomp_1_tensorcomp_1_petscspace_degree 1 \
14496c2660fcSDarsh Nathawani -petscdualspace_form_degree -1 \
14506c2660fcSDarsh Nathawani -petscdualspace_order 1 \
14516c2660fcSDarsh Nathawani -petscdualspace_lagrange_trimmed true
14526c2660fcSDarsh Nathawani
1453c4762a1bSJed Brown TEST*/
1454c4762a1bSJed Brown
1455c4762a1bSJed Brown /*
1456c4762a1bSJed Brown # 2D Q_2 on a quadrilaterial Plex
1457c4762a1bSJed Brown test:
1458c4762a1bSJed Brown suffix: q2_2d_plex_0
145930602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -convergence
1460c4762a1bSJed Brown test:
1461c4762a1bSJed Brown suffix: q2_2d_plex_1
146230602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 1
1463c4762a1bSJed Brown test:
1464c4762a1bSJed Brown suffix: q2_2d_plex_2
146530602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 2
1466c4762a1bSJed Brown test:
1467c4762a1bSJed Brown suffix: q2_2d_plex_3
146830602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 1 -shear_coords
1469c4762a1bSJed Brown test:
1470c4762a1bSJed Brown suffix: q2_2d_plex_4
147130602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 2 -shear_coords
1472c4762a1bSJed Brown test:
1473c4762a1bSJed Brown suffix: q2_2d_plex_5
147430602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_poly_tensor 1 -qorder 2 -porder 0 -non_affine_coords
1475c4762a1bSJed Brown test:
1476c4762a1bSJed Brown suffix: q2_2d_plex_6
147730602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_poly_tensor 1 -qorder 2 -porder 1 -non_affine_coords
1478c4762a1bSJed Brown test:
1479c4762a1bSJed Brown suffix: q2_2d_plex_7
148030602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_poly_tensor 1 -qorder 2 -porder 2 -non_affine_coords
1481c4762a1bSJed Brown
1482c4762a1bSJed Brown test:
1483c4762a1bSJed Brown suffix: p1d_2d_6
1484e788613bSJoe Wallwork requires: mmg
1485c4762a1bSJed Brown args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -convergence -conv_refine 0
1486c4762a1bSJed Brown test:
1487c4762a1bSJed Brown suffix: p1d_2d_7
1488e788613bSJoe Wallwork requires: mmg
1489c4762a1bSJed Brown args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 1 -conv_refine 0
1490c4762a1bSJed Brown test:
1491c4762a1bSJed Brown suffix: p1d_2d_8
1492e788613bSJoe Wallwork requires: mmg
1493c4762a1bSJed Brown args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 2 -conv_refine 0
1494c4762a1bSJed Brown */
1495