1a2424aedSMatthew G. Knepley static char help[] = "Tests for Gauss' Law\n\n"; 2a2424aedSMatthew G. Knepley 3a2424aedSMatthew G. Knepley /* We want to check the weak version of Gauss' Law, namely that 4a2424aedSMatthew G. Knepley 5a2424aedSMatthew G. Knepley \int_\Omega v div q - \int_\Gamma v (q \cdot n) = 0 6a2424aedSMatthew G. Knepley 7a2424aedSMatthew G. Knepley */ 8a2424aedSMatthew G. Knepley 9a2424aedSMatthew G. Knepley #include <petscdmplex.h> 10a2424aedSMatthew G. Knepley #include <petscsnes.h> 11a2424aedSMatthew G. Knepley #include <petscds.h> 12a2424aedSMatthew G. Knepley 13a2424aedSMatthew G. Knepley typedef struct { 14a2424aedSMatthew G. Knepley PetscInt degree; // The degree of the discretization 15a2424aedSMatthew G. Knepley PetscBool divFree; // True if the solution is divergence-free 16a2424aedSMatthew G. Knepley } AppCtx; 17a2424aedSMatthew G. Knepley 18a2424aedSMatthew G. Knepley static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 19a2424aedSMatthew G. Knepley { 20a2424aedSMatthew G. Knepley u[0] = 0.0; 21a2424aedSMatthew G. Knepley return PETSC_SUCCESS; 22a2424aedSMatthew G. Knepley } 23a2424aedSMatthew G. Knepley 24a2424aedSMatthew G. Knepley // div <x^d y^{d-1}, -x^{d-1} y^d> = 0 25a2424aedSMatthew G. Knepley static void solenoidal_2d(PetscInt n, const PetscReal x[], PetscScalar u[]) 26a2424aedSMatthew G. Knepley { 27a2424aedSMatthew G. Knepley u[0] = PetscPowRealInt(x[0], n) * PetscPowRealInt(x[1], n - 1); 28a2424aedSMatthew G. Knepley u[1] = -PetscPowRealInt(x[0], n - 1) * PetscPowRealInt(x[1], n); 29a2424aedSMatthew G. Knepley } 30a2424aedSMatthew G. Knepley // div <x^d y^{d-1} z^{d-1}, -2 x^{d-1} y^d z^{d-1}, x^{d-1} y^{d-1} z^d> = 0 31a2424aedSMatthew G. Knepley static void solenoidal_3d(PetscInt n, const PetscReal x[], PetscScalar u[]) 32a2424aedSMatthew G. Knepley { 33a2424aedSMatthew G. Knepley u[0] = PetscPowRealInt(x[0], n) * PetscPowRealInt(x[1], n - 1) * PetscPowRealInt(x[2], n - 1); 34a2424aedSMatthew G. Knepley u[1] = -2. * PetscPowRealInt(x[0], n - 1) * PetscPowRealInt(x[1], n) * PetscPowRealInt(x[2], n - 1); 35a2424aedSMatthew G. Knepley u[2] = PetscPowRealInt(x[0], n - 1) * PetscPowRealInt(x[1], n - 1) * PetscPowRealInt(x[2], n); 36a2424aedSMatthew G. Knepley } 37a2424aedSMatthew G. Knepley 38a2424aedSMatthew G. Knepley static PetscErrorCode solenoidal_totaldeg_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 39a2424aedSMatthew G. Knepley { 40a2424aedSMatthew G. Knepley const PetscInt deg = *(PetscInt *)ctx; 41a2424aedSMatthew G. Knepley const PetscInt n = deg / 2 + deg % 2; 42a2424aedSMatthew G. Knepley 43a2424aedSMatthew G. Knepley solenoidal_2d(n, x, u); 44a2424aedSMatthew G. Knepley return PETSC_SUCCESS; 45a2424aedSMatthew G. Knepley } 46a2424aedSMatthew G. Knepley 47a2424aedSMatthew G. Knepley static PetscErrorCode solenoidal_totaldeg_3d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 48a2424aedSMatthew G. Knepley { 49a2424aedSMatthew G. Knepley const PetscInt deg = *(PetscInt *)ctx; 50a2424aedSMatthew G. Knepley const PetscInt n = deg / 3 + (deg % 3 ? 1 : 0); 51a2424aedSMatthew G. Knepley 52a2424aedSMatthew G. Knepley solenoidal_3d(n, x, u); 53a2424aedSMatthew G. Knepley return PETSC_SUCCESS; 54a2424aedSMatthew G. Knepley } 55a2424aedSMatthew G. Knepley 56a2424aedSMatthew G. Knepley // This is in P_n^{-} 57a2424aedSMatthew G. Knepley static PetscErrorCode source_totaldeg(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 58a2424aedSMatthew G. Knepley { 59a2424aedSMatthew G. Knepley const PetscInt n = *(PetscInt *)ctx; 60a2424aedSMatthew G. Knepley 61a2424aedSMatthew G. Knepley for (PetscInt d = 0; d < dim; ++d) u[d] = PetscPowRealInt(x[d], n + 1); 62a2424aedSMatthew G. Knepley return PETSC_SUCCESS; 63a2424aedSMatthew G. Knepley } 64a2424aedSMatthew G. Knepley 65a2424aedSMatthew G. Knepley static void identity(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 66a2424aedSMatthew G. Knepley { 67a2424aedSMatthew G. Knepley const PetscInt deg = (PetscInt)PetscRealPart(constants[0]); 68a2424aedSMatthew G. Knepley PetscScalar p[3]; 69a2424aedSMatthew G. Knepley 70a2424aedSMatthew G. Knepley if (dim == 2) PetscCallVoid(solenoidal_totaldeg_2d(dim, t, x, uOff[1] - uOff[0], p, (void *)°)); 71a2424aedSMatthew G. Knepley else PetscCallVoid(solenoidal_totaldeg_3d(dim, t, x, uOff[1] - uOff[0], p, (void *)°)); 72a2424aedSMatthew G. Knepley for (PetscInt c = 0; c < dim; ++c) f0[c] = -u[c] + p[c]; 73a2424aedSMatthew G. Knepley } 74a2424aedSMatthew G. Knepley 75a2424aedSMatthew G. Knepley static void zero_bd(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 76a2424aedSMatthew G. Knepley { 77a2424aedSMatthew G. Knepley for (PetscInt d = 0; d < dim; ++d) f0[0] = 0.; 78a2424aedSMatthew G. Knepley } 79a2424aedSMatthew G. Knepley 80a2424aedSMatthew G. Knepley static void flux(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 81a2424aedSMatthew G. Knepley { 82a2424aedSMatthew G. Knepley for (PetscInt d = 0; d < dim; ++d) f0[0] -= u[d] * n[d]; 83a2424aedSMatthew G. Knepley } 84a2424aedSMatthew G. Knepley 85a2424aedSMatthew G. Knepley static void divergence(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 86a2424aedSMatthew G. Knepley { 87a2424aedSMatthew G. Knepley for (PetscInt d = 0; d < dim; ++d) f0[0] += u_x[d * dim + d]; 88a2424aedSMatthew G. Knepley } 89a2424aedSMatthew G. Knepley 90a2424aedSMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm) 91a2424aedSMatthew G. Knepley { 92a2424aedSMatthew G. Knepley PetscFunctionBegin; 93a2424aedSMatthew G. Knepley PetscCall(DMCreate(comm, dm)); 94a2424aedSMatthew G. Knepley PetscCall(DMSetType(*dm, DMPLEX)); 95a2424aedSMatthew G. Knepley PetscCall(DMSetFromOptions(*dm)); 96a2424aedSMatthew G. Knepley PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 97a2424aedSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 98a2424aedSMatthew G. Knepley } 99a2424aedSMatthew G. Knepley 100a2424aedSMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 101a2424aedSMatthew G. Knepley { 102a2424aedSMatthew G. Knepley options->degree = -1; 103a2424aedSMatthew G. Knepley 104a2424aedSMatthew G. Knepley PetscFunctionBeginUser; 105a2424aedSMatthew G. Knepley PetscOptionsBegin(comm, "", "Gauss' Law Test Options", "DMPLEX"); 106a2424aedSMatthew G. Knepley PetscOptionsEnd(); 107a2424aedSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 108a2424aedSMatthew G. Knepley } 109a2424aedSMatthew G. Knepley 110a2424aedSMatthew G. Knepley static PetscErrorCode CreateDiscretization(DM dm, AppCtx *user) 111a2424aedSMatthew G. Knepley { 112a2424aedSMatthew G. Knepley PetscFE feq, fep; 113a2424aedSMatthew G. Knepley PetscSpace sp; 114a2424aedSMatthew G. Knepley PetscQuadrature quad, fquad; 115a2424aedSMatthew G. Knepley PetscDS ds; 116a2424aedSMatthew G. Knepley DMLabel label; 117a2424aedSMatthew G. Knepley DMPolytopeType ct; 118a2424aedSMatthew G. Knepley PetscInt dim, cStart, minDeg, maxDeg; 119a2424aedSMatthew G. Knepley PetscBool isTrimmed, isSum; 120a2424aedSMatthew G. Knepley 121a2424aedSMatthew G. Knepley PetscFunctionBeginUser; 122a2424aedSMatthew G. Knepley PetscCall(DMGetDimension(dm, &dim)); 123a2424aedSMatthew G. Knepley PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 124a2424aedSMatthew G. Knepley PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 125a2424aedSMatthew G. Knepley PetscCall(PetscFECreateByCell(PETSC_COMM_SELF, dim, dim, ct, "field_", -1, &feq)); 126a2424aedSMatthew G. Knepley PetscCall(DMSetField(dm, 0, NULL, (PetscObject)feq)); 127a2424aedSMatthew G. Knepley PetscCall(PetscFEGetQuadrature(feq, &quad)); 128a2424aedSMatthew G. Knepley PetscCall(PetscFEGetFaceQuadrature(feq, &fquad)); 129a2424aedSMatthew G. Knepley PetscCall(PetscFEGetBasisSpace(feq, &sp)); 130a2424aedSMatthew G. Knepley PetscCall(PetscSpaceGetDegree(sp, &minDeg, &maxDeg)); 131a2424aedSMatthew G. Knepley PetscCall(PetscObjectTypeCompare((PetscObject)sp, PETSCSPACEPTRIMMED, &isTrimmed)); 132a2424aedSMatthew G. Knepley PetscCall(PetscObjectTypeCompare((PetscObject)sp, PETSCSPACESUM, &isSum)); 133a2424aedSMatthew G. Knepley if (isSum) { 134a2424aedSMatthew G. Knepley PetscSpace subsp, xsp, ysp; 135a2424aedSMatthew G. Knepley PetscInt xdeg, ydeg; 136a2424aedSMatthew G. Knepley PetscBool isTensor; 137a2424aedSMatthew G. Knepley 138a2424aedSMatthew G. Knepley PetscCall(PetscSpaceSumGetSubspace(sp, 0, &subsp)); 139a2424aedSMatthew G. Knepley PetscCall(PetscObjectTypeCompare((PetscObject)subsp, PETSCSPACETENSOR, &isTensor)); 140a2424aedSMatthew G. Knepley if (isTensor) { 141a2424aedSMatthew G. Knepley PetscCall(PetscSpaceTensorGetSubspace(subsp, 0, &xsp)); 142a2424aedSMatthew G. Knepley PetscCall(PetscSpaceTensorGetSubspace(subsp, 1, &ysp)); 143a2424aedSMatthew G. Knepley PetscCall(PetscSpaceGetDegree(xsp, &xdeg, NULL)); 144a2424aedSMatthew G. Knepley PetscCall(PetscSpaceGetDegree(ysp, &ydeg, NULL)); 145a2424aedSMatthew G. Knepley isTrimmed = xdeg != ydeg ? PETSC_TRUE : PETSC_FALSE; 146a2424aedSMatthew G. Knepley } 147a2424aedSMatthew G. Knepley } 148a2424aedSMatthew G. Knepley user->degree = minDeg; 149a2424aedSMatthew G. Knepley if (isTrimmed) user->divFree = PETSC_FALSE; 150a2424aedSMatthew G. Knepley else user->divFree = PETSC_TRUE; 151a2424aedSMatthew G. Knepley PetscCheck(!user->divFree || user->degree, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Degree 0 solution not available"); 152a2424aedSMatthew G. Knepley PetscCall(PetscFEDestroy(&feq)); 153a2424aedSMatthew G. Knepley PetscCall(PetscFECreateByCell(PETSC_COMM_SELF, dim, 1, ct, "pot_", -1, &fep)); 154a2424aedSMatthew G. Knepley PetscCall(DMSetField(dm, 1, NULL, (PetscObject)fep)); 155a2424aedSMatthew G. Knepley PetscCall(PetscFESetQuadrature(fep, quad)); 156a2424aedSMatthew G. Knepley PetscCall(PetscFESetFaceQuadrature(fep, fquad)); 157a2424aedSMatthew G. Knepley PetscCall(PetscFEDestroy(&fep)); 158a2424aedSMatthew G. Knepley PetscCall(DMCreateDS(dm)); 159a2424aedSMatthew G. Knepley 160a2424aedSMatthew G. Knepley PetscCall(DMGetDS(dm, &ds)); 161a2424aedSMatthew G. Knepley PetscCall(PetscDSSetResidual(ds, 0, identity, NULL)); 162a2424aedSMatthew G. Knepley PetscCall(PetscDSSetResidual(ds, 1, divergence, NULL)); 163a2424aedSMatthew G. Knepley if (user->divFree) { 164a2424aedSMatthew G. Knepley if (dim == 2) PetscCall(PetscDSSetExactSolution(ds, 0, solenoidal_totaldeg_2d, &user->degree)); 165a2424aedSMatthew G. Knepley else PetscCall(PetscDSSetExactSolution(ds, 0, solenoidal_totaldeg_3d, &user->degree)); 166a2424aedSMatthew G. Knepley } else { 167a2424aedSMatthew G. Knepley PetscCall(PetscDSSetExactSolution(ds, 0, source_totaldeg, &user->degree)); 168a2424aedSMatthew G. Knepley } 169a2424aedSMatthew G. Knepley PetscCall(PetscDSSetExactSolution(ds, 1, zero, &user->degree)); 170a2424aedSMatthew G. Knepley PetscCall(DMGetLabel(dm, "marker", &label)); 171a2424aedSMatthew G. Knepley 172a2424aedSMatthew G. Knepley // TODO Can we also test the boundary residual integration? 173a2424aedSMatthew G. Knepley //PetscWeakForm wf; 174a2424aedSMatthew G. Knepley //PetscInt bd, id = 1; 175a2424aedSMatthew G. Knepley //PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "boundary", label, 1, &id, 1, 0, NULL, NULL, NULL, user, &bd)); 176a2424aedSMatthew G. Knepley //PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 177a2424aedSMatthew G. Knepley //PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, 1, 0, 0, flux, 0, NULL)); 178a2424aedSMatthew G. Knepley 179a2424aedSMatthew G. Knepley { 180a2424aedSMatthew G. Knepley PetscScalar constants[1]; 181a2424aedSMatthew G. Knepley 182a2424aedSMatthew G. Knepley constants[0] = user->degree; 183a2424aedSMatthew G. Knepley PetscCall(PetscDSSetConstants(ds, 1, constants)); 184a2424aedSMatthew G. Knepley } 185a2424aedSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 186a2424aedSMatthew G. Knepley } 187a2424aedSMatthew G. Knepley 188a2424aedSMatthew G. Knepley int main(int argc, char **argv) 189a2424aedSMatthew G. Knepley { 190a2424aedSMatthew G. Knepley DM dm; 191a2424aedSMatthew G. Knepley SNES snes; 192a2424aedSMatthew G. Knepley Vec u; 193a2424aedSMatthew G. Knepley PetscReal error[2], residual; 194a2424aedSMatthew G. Knepley PetscScalar source[2], outflow[2]; 195a2424aedSMatthew G. Knepley AppCtx user; 196a2424aedSMatthew G. Knepley 197a2424aedSMatthew G. Knepley PetscFunctionBeginUser; 198a2424aedSMatthew G. Knepley PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 199a2424aedSMatthew G. Knepley PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 200a2424aedSMatthew G. Knepley PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm)); 201a2424aedSMatthew G. Knepley PetscCall(CreateDiscretization(dm, &user)); 202a2424aedSMatthew G. Knepley PetscCall(DMGetGlobalVector(dm, &u)); 203a2424aedSMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)u, "solution")); 204a2424aedSMatthew G. Knepley PetscCall(DMComputeExactSolution(dm, 0., u, NULL)); 205a2424aedSMatthew G. Knepley 206a2424aedSMatthew G. Knepley PetscCall(SNESCreate(PetscObjectComm((PetscObject)dm), &snes)); 207a2424aedSMatthew G. Knepley PetscCall(SNESSetDM(snes, dm)); 208a2424aedSMatthew G. Knepley PetscCall(DMPlexSetSNESLocalFEM(dm, PETSC_FALSE, &user)); 209a2424aedSMatthew G. Knepley PetscCall(SNESSetFromOptions(snes)); 210a2424aedSMatthew G. Knepley PetscCall(DMSNESCheckDiscretization(snes, dm, 0., u, PETSC_DETERMINE, error)); 211a2424aedSMatthew G. Knepley PetscCheck(PetscAbsReal(error[0]) < PETSC_SMALL, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Exact solution does not fit into FEM space: %g should be zero", (double)error[0]); 212a2424aedSMatthew G. Knepley if (user.divFree) { 213a2424aedSMatthew G. Knepley PetscCall(DMSNESCheckResidual(snes, dm, u, PETSC_DETERMINE, &residual)); 214a2424aedSMatthew G. Knepley PetscCheck(PetscAbsReal(residual) < PETSC_SMALL, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Exact solution is not divergence-free: %g should be zero", (double)residual); 215a2424aedSMatthew G. Knepley } else { 216a2424aedSMatthew G. Knepley PetscDS ds; 217a2424aedSMatthew G. Knepley 218a2424aedSMatthew G. Knepley PetscCall(DMGetDS(dm, &ds)); 219a2424aedSMatthew G. Knepley PetscCall(PetscDSSetObjective(ds, 1, divergence)); 220a2424aedSMatthew G. Knepley PetscCall(DMPlexComputeIntegralFEM(dm, u, source, &user)); 221a2424aedSMatthew G. Knepley } 222a2424aedSMatthew G. Knepley PetscCall(SNESDestroy(&snes)); 223a2424aedSMatthew G. Knepley 224*2192575eSBarry Smith PetscBdPointFn *funcs[] = {zero_bd, flux}; 225a2424aedSMatthew G. Knepley DMLabel label; 226a2424aedSMatthew G. Knepley PetscInt id = 1; 227a2424aedSMatthew G. Knepley 228a2424aedSMatthew G. Knepley PetscCall(DMGetLabel(dm, "marker", &label)); 229a2424aedSMatthew G. Knepley PetscCall(DMPlexComputeBdIntegral(dm, u, label, 1, &id, funcs, outflow, &user)); 230a2424aedSMatthew G. Knepley if (user.divFree) PetscCheck(PetscAbsScalar(outflow[1]) < PETSC_SMALL, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Outflow %g should be zero for a divergence-free field", (double)PetscRealPart(outflow[1])); 231a2424aedSMatthew G. Knepley else PetscCheck(PetscAbsScalar(source[1] + outflow[1]) < PETSC_SMALL, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Outflow %g should oppose source %g", (double)PetscRealPart(outflow[1]), (double)PetscRealPart(source[1])); 232a2424aedSMatthew G. Knepley 233a2424aedSMatthew G. Knepley PetscCall(DMRestoreGlobalVector(dm, &u)); 234a2424aedSMatthew G. Knepley PetscCall(DMDestroy(&dm)); 235a2424aedSMatthew G. Knepley PetscCall(PetscFinalize()); 236a2424aedSMatthew G. Knepley return 0; 237a2424aedSMatthew G. Knepley } 238a2424aedSMatthew G. Knepley 239a2424aedSMatthew G. Knepley /*TEST 240a2424aedSMatthew G. Knepley 241a2424aedSMatthew G. Knepley testset: 242a2424aedSMatthew G. Knepley suffix: p 243a2424aedSMatthew G. Knepley requires: triangle ctetgen 244a2424aedSMatthew G. Knepley args: -dm_plex_dim {{2 3}} -dm_plex_box_faces 2,2,2 245a2424aedSMatthew G. Knepley 246a2424aedSMatthew G. Knepley test: 247a2424aedSMatthew G. Knepley suffix: 1 248a2424aedSMatthew G. Knepley args: -field_petscspace_degree 1 -pot_petscspace_degree 1 249a2424aedSMatthew G. Knepley test: 250a2424aedSMatthew G. Knepley suffix: 2 251a2424aedSMatthew G. Knepley args: -field_petscspace_degree 2 -pot_petscspace_degree 2 252a2424aedSMatthew G. Knepley test: 253a2424aedSMatthew G. Knepley suffix: 3 254a2424aedSMatthew G. Knepley args: -field_petscspace_degree 3 -pot_petscspace_degree 3 255a2424aedSMatthew G. Knepley test: 256a2424aedSMatthew G. Knepley suffix: 4 257a2424aedSMatthew G. Knepley args: -field_petscspace_degree 4 -pot_petscspace_degree 4 258a2424aedSMatthew G. Knepley 259a2424aedSMatthew G. Knepley testset: 260a2424aedSMatthew G. Knepley suffix: q 261a2424aedSMatthew G. Knepley args: -dm_plex_dim {{2 3}} -dm_plex_simplex 0 -dm_plex_box_faces 2,2 262a2424aedSMatthew G. Knepley 263a2424aedSMatthew G. Knepley test: 264a2424aedSMatthew G. Knepley suffix: 1 265a2424aedSMatthew G. Knepley args: -field_petscspace_degree 1 -pot_petscspace_degree 1 266a2424aedSMatthew G. Knepley test: 267a2424aedSMatthew G. Knepley suffix: 2 268a2424aedSMatthew G. Knepley args: -field_petscspace_degree 2 -pot_petscspace_degree 2 269a2424aedSMatthew G. Knepley test: 270a2424aedSMatthew G. Knepley suffix: 3 271a2424aedSMatthew G. Knepley args: -field_petscspace_degree 3 -pot_petscspace_degree 3 272a2424aedSMatthew G. Knepley test: 273a2424aedSMatthew G. Knepley suffix: 4 274a2424aedSMatthew G. Knepley args: -field_petscspace_degree 4 -pot_petscspace_degree 4 275a2424aedSMatthew G. Knepley 276a2424aedSMatthew G. Knepley testset: 277a2424aedSMatthew G. Knepley suffix: bdm 278a2424aedSMatthew G. Knepley requires: triangle ctetgen 279a2424aedSMatthew G. Knepley args: -dm_plex_dim 2 -dm_plex_box_faces 2,2 280a2424aedSMatthew G. Knepley 281a2424aedSMatthew G. Knepley test: 282a2424aedSMatthew G. Knepley suffix: 1 283a2424aedSMatthew G. Knepley args: -pot_petscspace_degree 0 -pot_petscdualspace_lagrange_continuity 0 \ 284a2424aedSMatthew G. Knepley -field_petscspace_degree 1 -field_petscdualspace_type bdm \ 285a2424aedSMatthew G. Knepley -field_petscfe_default_quadrature_order 2 286a2424aedSMatthew G. Knepley 287a2424aedSMatthew G. Knepley testset: 288a2424aedSMatthew G. Knepley suffix: rt 289a2424aedSMatthew G. Knepley requires: triangle ctetgen 290a2424aedSMatthew G. Knepley args: -dm_plex_dim 2 -dm_plex_box_faces 2,2 291a2424aedSMatthew G. Knepley 292a2424aedSMatthew G. Knepley test: 293a2424aedSMatthew G. Knepley suffix: 1 294a2424aedSMatthew G. Knepley args: -pot_petscspace_degree 0 -pot_petscdualspace_lagrange_continuity 0 \ 295a2424aedSMatthew G. Knepley -field_petscspace_type ptrimmed \ 296a2424aedSMatthew G. Knepley -field_petscspace_components 2 \ 297a2424aedSMatthew G. Knepley -field_petscspace_ptrimmed_form_degree -1 \ 298a2424aedSMatthew G. Knepley -field_petscdualspace_order 1 \ 299a2424aedSMatthew G. Knepley -field_petscdualspace_form_degree -1 \ 300a2424aedSMatthew G. Knepley -field_petscdualspace_lagrange_trimmed true \ 301a2424aedSMatthew G. Knepley -field_petscfe_default_quadrature_order 2 302a2424aedSMatthew G. Knepley 303a2424aedSMatthew G. Knepley testset: 304a2424aedSMatthew G. Knepley suffix: rtq 305a2424aedSMatthew G. Knepley requires: triangle ctetgen 306a2424aedSMatthew G. Knepley args: -dm_plex_dim 2 -dm_plex_simplex 0 -dm_plex_box_faces 2,2 307a2424aedSMatthew G. Knepley 308a2424aedSMatthew G. Knepley test: 309a2424aedSMatthew G. Knepley suffix: 1 310a2424aedSMatthew G. Knepley args: -pot_petscspace_degree 0 -pot_petscdualspace_lagrange_continuity 0 \ 311a2424aedSMatthew G. Knepley -field_petscspace_degree 1 \ 312a2424aedSMatthew G. Knepley -field_petscspace_type sum \ 313a2424aedSMatthew G. Knepley -field_petscspace_variables 2 \ 314a2424aedSMatthew G. Knepley -field_petscspace_components 2 \ 315a2424aedSMatthew G. Knepley -field_petscspace_sum_spaces 2 \ 316a2424aedSMatthew G. Knepley -field_petscspace_sum_concatenate true \ 317a2424aedSMatthew G. Knepley -field_sumcomp_0_petscspace_variables 2 \ 318a2424aedSMatthew G. Knepley -field_sumcomp_0_petscspace_type tensor \ 319a2424aedSMatthew G. Knepley -field_sumcomp_0_petscspace_tensor_spaces 2 \ 320a2424aedSMatthew G. Knepley -field_sumcomp_0_petscspace_tensor_uniform false \ 321a2424aedSMatthew G. Knepley -field_sumcomp_0_tensorcomp_0_petscspace_degree 1 \ 322a2424aedSMatthew G. Knepley -field_sumcomp_0_tensorcomp_1_petscspace_degree 0 \ 323a2424aedSMatthew G. Knepley -field_sumcomp_1_petscspace_variables 2 \ 324a2424aedSMatthew G. Knepley -field_sumcomp_1_petscspace_type tensor \ 325a2424aedSMatthew G. Knepley -field_sumcomp_1_petscspace_tensor_spaces 2 \ 326a2424aedSMatthew G. Knepley -field_sumcomp_1_petscspace_tensor_uniform false \ 327a2424aedSMatthew G. Knepley -field_sumcomp_1_tensorcomp_0_petscspace_degree 0 \ 328a2424aedSMatthew G. Knepley -field_sumcomp_1_tensorcomp_1_petscspace_degree 1 \ 329a2424aedSMatthew G. Knepley -field_petscdualspace_order 1 \ 330a2424aedSMatthew G. Knepley -field_petscdualspace_form_degree -1 \ 331a2424aedSMatthew G. Knepley -field_petscdualspace_lagrange_trimmed true \ 332a2424aedSMatthew G. Knepley -field_petscfe_default_quadrature_order 2 333a2424aedSMatthew G. Knepley 334a2424aedSMatthew G. Knepley TEST*/ 335