1*d21efd2eSMatthew G. Knepley static const char help[] = "Tests for determining whether a new finite element works"; 2*d21efd2eSMatthew G. Knepley 3*d21efd2eSMatthew G. Knepley /* 4*d21efd2eSMatthew G. Knepley Use -interpolation_view and -l2_projection_view to look at the interpolants. 5*d21efd2eSMatthew G. Knepley */ 6*d21efd2eSMatthew G. Knepley 7*d21efd2eSMatthew G. Knepley #include <petscdmplex.h> 8*d21efd2eSMatthew G. Knepley #include <petscfe.h> 9*d21efd2eSMatthew G. Knepley #include <petscds.h> 10*d21efd2eSMatthew G. Knepley #include <petscsnes.h> 11*d21efd2eSMatthew G. Knepley 12*d21efd2eSMatthew G. Knepley static void constant(PetscInt dim, PetscInt Nf, PetscInt NfAux, 13*d21efd2eSMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 14*d21efd2eSMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 15*d21efd2eSMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 16*d21efd2eSMatthew G. Knepley { 17*d21efd2eSMatthew G. Knepley const PetscInt Nc = uOff[1] - uOff[0]; 18*d21efd2eSMatthew G. Knepley for (PetscInt c = 0; c < Nc; ++c) f0[c] += 5.; 19*d21efd2eSMatthew G. Knepley } 20*d21efd2eSMatthew G. Knepley 21*d21efd2eSMatthew G. Knepley static void linear(PetscInt dim, PetscInt Nf, PetscInt NfAux, 22*d21efd2eSMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 23*d21efd2eSMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 24*d21efd2eSMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 25*d21efd2eSMatthew G. Knepley { 26*d21efd2eSMatthew G. Knepley const PetscInt Nc = uOff[1] - uOff[0]; 27*d21efd2eSMatthew G. Knepley for (PetscInt c = 0; c < Nc; ++c) f0[c] += 5.*x[c]; 28*d21efd2eSMatthew G. Knepley } 29*d21efd2eSMatthew G. Knepley 30*d21efd2eSMatthew G. Knepley static void quadratic(PetscInt dim, PetscInt Nf, PetscInt NfAux, 31*d21efd2eSMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 32*d21efd2eSMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 33*d21efd2eSMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 34*d21efd2eSMatthew G. Knepley { 35*d21efd2eSMatthew G. Knepley const PetscInt Nc = uOff[1] - uOff[0]; 36*d21efd2eSMatthew G. Knepley for (PetscInt c = 0; c < Nc; ++c) f0[c] += 5.*x[c]*x[c]; 37*d21efd2eSMatthew G. Knepley } 38*d21efd2eSMatthew G. Knepley 39*d21efd2eSMatthew G. Knepley static void trig(PetscInt dim, PetscInt Nf, PetscInt NfAux, 40*d21efd2eSMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 41*d21efd2eSMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 42*d21efd2eSMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 43*d21efd2eSMatthew G. Knepley { 44*d21efd2eSMatthew G. Knepley const PetscInt Nc = uOff[1] - uOff[0]; 45*d21efd2eSMatthew G. Knepley for (PetscInt c = 0; c < Nc; ++c) f0[c] += PetscCosReal(2.*PETSC_PI*x[c]); 46*d21efd2eSMatthew G. Knepley } 47*d21efd2eSMatthew G. Knepley 48*d21efd2eSMatthew G. Knepley static const char *names[] = {"constant", "linear", "quadratic", "trig"}; 49*d21efd2eSMatthew G. Knepley static PetscPointFunc functions[] = { constant, linear, quadratic, trig}; 50*d21efd2eSMatthew G. Knepley 51*d21efd2eSMatthew G. Knepley typedef struct { 52*d21efd2eSMatthew G. Knepley PetscPointFunc exactSol; 53*d21efd2eSMatthew G. Knepley } AppCtx; 54*d21efd2eSMatthew G. Knepley 55*d21efd2eSMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 56*d21efd2eSMatthew G. Knepley { 57*d21efd2eSMatthew G. Knepley char name[PETSC_MAX_PATH_LEN] = "constant"; 58*d21efd2eSMatthew G. Knepley PetscInt Nfunc = sizeof(names)/sizeof(char *), i; 59*d21efd2eSMatthew G. Knepley PetscErrorCode ierr; 60*d21efd2eSMatthew G. Knepley 61*d21efd2eSMatthew G. Knepley PetscFunctionBeginUser; 62*d21efd2eSMatthew G. Knepley options->exactSol = NULL; 63*d21efd2eSMatthew G. Knepley 64*d21efd2eSMatthew G. Knepley ierr = PetscOptionsBegin(comm, "", "FE Test Options", "PETSCFE");CHKERRQ(ierr); 65*d21efd2eSMatthew G. Knepley ierr = PetscOptionsString("-func", "Function to project into space", "", name, name, PETSC_MAX_PATH_LEN, NULL);CHKERRQ(ierr); 66*d21efd2eSMatthew G. Knepley ierr = PetscOptionsEnd();CHKERRQ(ierr); 67*d21efd2eSMatthew G. Knepley 68*d21efd2eSMatthew G. Knepley for (i = 0; i < Nfunc; ++i) { 69*d21efd2eSMatthew G. Knepley PetscBool flg; 70*d21efd2eSMatthew G. Knepley 71*d21efd2eSMatthew G. Knepley ierr = PetscStrcmp(name, names[i], &flg);CHKERRQ(ierr); 72*d21efd2eSMatthew G. Knepley if (flg) {options->exactSol = functions[i]; break;} 73*d21efd2eSMatthew G. Knepley } 74*d21efd2eSMatthew G. Knepley if (!options->exactSol) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Invalid test function %s", name); 75*d21efd2eSMatthew G. Knepley PetscFunctionReturn(0); 76*d21efd2eSMatthew G. Knepley } 77*d21efd2eSMatthew G. Knepley 78*d21efd2eSMatthew G. Knepley /* The exact solution is the negative of the f0 contribution */ 79*d21efd2eSMatthew G. Knepley static PetscErrorCode exactSolution(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 80*d21efd2eSMatthew G. Knepley { 81*d21efd2eSMatthew G. Knepley AppCtx *user = (AppCtx *) ctx; 82*d21efd2eSMatthew G. Knepley PetscInt uOff[2] = {0., Nc}; 83*d21efd2eSMatthew G. Knepley 84*d21efd2eSMatthew G. Knepley user->exactSol(dim, 1, 0, uOff, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, time, x, 0, NULL, u); 85*d21efd2eSMatthew G. Knepley for (PetscInt c = 0; c < Nc; ++c) u[c] *= -1.; 86*d21efd2eSMatthew G. Knepley return 0; 87*d21efd2eSMatthew G. Knepley } 88*d21efd2eSMatthew G. Knepley 89*d21efd2eSMatthew G. Knepley static void f0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 90*d21efd2eSMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 91*d21efd2eSMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 92*d21efd2eSMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 93*d21efd2eSMatthew G. Knepley { 94*d21efd2eSMatthew G. Knepley const PetscInt Nc = uOff[1] - uOff[0]; 95*d21efd2eSMatthew G. Knepley for (PetscInt c = 0; c < Nc; ++c) f0[c] += u[c]; 96*d21efd2eSMatthew G. Knepley } 97*d21efd2eSMatthew G. Knepley 98*d21efd2eSMatthew G. Knepley static void g0(PetscInt dim, PetscInt Nf, PetscInt NfAux, 99*d21efd2eSMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 100*d21efd2eSMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 101*d21efd2eSMatthew G. Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 102*d21efd2eSMatthew G. Knepley { 103*d21efd2eSMatthew G. Knepley const PetscInt Nc = uOff[1] - uOff[0]; 104*d21efd2eSMatthew G. Knepley for (PetscInt c = 0; c < Nc; ++c) g0[c*Nc+c] = 1.0; 105*d21efd2eSMatthew G. Knepley } 106*d21efd2eSMatthew G. Knepley 107*d21efd2eSMatthew G. Knepley static PetscErrorCode SetupProblem(DM dm, AppCtx *user) 108*d21efd2eSMatthew G. Knepley { 109*d21efd2eSMatthew G. Knepley PetscDS ds; 110*d21efd2eSMatthew G. Knepley PetscWeakForm wf; 111*d21efd2eSMatthew G. Knepley PetscErrorCode ierr; 112*d21efd2eSMatthew G. Knepley 113*d21efd2eSMatthew G. Knepley PetscFunctionBeginUser; 114*d21efd2eSMatthew G. Knepley ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 115*d21efd2eSMatthew G. Knepley ierr = PetscDSGetWeakForm(ds, &wf);CHKERRQ(ierr); 116*d21efd2eSMatthew G. Knepley ierr = PetscWeakFormSetIndexResidual(wf, NULL, 0, 0, 0, 0, f0, 0, NULL);CHKERRQ(ierr); 117*d21efd2eSMatthew G. Knepley ierr = PetscWeakFormSetIndexResidual(wf, NULL, 0, 0, 0, 1, user->exactSol, 0, NULL);CHKERRQ(ierr); 118*d21efd2eSMatthew G. Knepley ierr = PetscWeakFormSetIndexJacobian(wf, NULL, 0, 0, 0, 0, 0, g0, 0, NULL, 0, NULL, 0, NULL);CHKERRQ(ierr); 119*d21efd2eSMatthew G. Knepley ierr = PetscDSSetExactSolution(ds, 0, exactSolution, user);CHKERRQ(ierr); 120*d21efd2eSMatthew G. Knepley PetscFunctionReturn(0); 121*d21efd2eSMatthew G. Knepley } 122*d21efd2eSMatthew G. Knepley 123*d21efd2eSMatthew G. Knepley static PetscErrorCode SetupDiscretization(DM dm, const char name[], AppCtx *user) 124*d21efd2eSMatthew G. Knepley { 125*d21efd2eSMatthew G. Knepley DM cdm = dm; 126*d21efd2eSMatthew G. Knepley PetscFE fe; 127*d21efd2eSMatthew G. Knepley char prefix[PETSC_MAX_PATH_LEN]; 128*d21efd2eSMatthew G. Knepley PetscErrorCode ierr; 129*d21efd2eSMatthew G. Knepley 130*d21efd2eSMatthew G. Knepley PetscFunctionBeginUser; 131*d21efd2eSMatthew G. Knepley ierr = PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name);CHKERRQ(ierr); 132*d21efd2eSMatthew G. Knepley ierr = DMCreateFEDefault(dm, 1, name ? prefix : NULL, -1, &fe);CHKERRQ(ierr); 133*d21efd2eSMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) fe, name ? name : "Solution");CHKERRQ(ierr); 134*d21efd2eSMatthew G. Knepley /* Set discretization and boundary conditions for each mesh */ 135*d21efd2eSMatthew G. Knepley ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr); 136*d21efd2eSMatthew G. Knepley ierr = DMCreateDS(dm);CHKERRQ(ierr); 137*d21efd2eSMatthew G. Knepley ierr = SetupProblem(dm, user);CHKERRQ(ierr); 138*d21efd2eSMatthew G. Knepley while (cdm) { 139*d21efd2eSMatthew G. Knepley ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr); 140*d21efd2eSMatthew G. Knepley ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr); 141*d21efd2eSMatthew G. Knepley } 142*d21efd2eSMatthew G. Knepley ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 143*d21efd2eSMatthew G. Knepley PetscFunctionReturn(0); 144*d21efd2eSMatthew G. Knepley } 145*d21efd2eSMatthew G. Knepley 146*d21efd2eSMatthew G. Knepley /* This test tells us whether the given function is contained in the approximation space */ 147*d21efd2eSMatthew G. Knepley static PetscErrorCode CheckInterpolation(DM dm, AppCtx *user) 148*d21efd2eSMatthew G. Knepley { 149*d21efd2eSMatthew G. Knepley PetscSimplePointFunc exactSol[1]; 150*d21efd2eSMatthew G. Knepley void *exactCtx[1]; 151*d21efd2eSMatthew G. Knepley PetscDS ds; 152*d21efd2eSMatthew G. Knepley Vec u; 153*d21efd2eSMatthew G. Knepley PetscReal error, tol = PETSC_SMALL; 154*d21efd2eSMatthew G. Knepley MPI_Comm comm; 155*d21efd2eSMatthew G. Knepley PetscErrorCode ierr; 156*d21efd2eSMatthew G. Knepley 157*d21efd2eSMatthew G. Knepley PetscFunctionBeginUser; 158*d21efd2eSMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 159*d21efd2eSMatthew G. Knepley ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 160*d21efd2eSMatthew G. Knepley ierr = DMGetGlobalVector(dm, &u);CHKERRQ(ierr); 161*d21efd2eSMatthew G. Knepley ierr = PetscDSGetExactSolution(ds, 0, &exactSol[0], &exactCtx[0]);CHKERRQ(ierr); 162*d21efd2eSMatthew G. Knepley ierr = DMProjectFunction(dm, 0.0, exactSol, exactCtx, INSERT_ALL_VALUES, u);CHKERRQ(ierr); 163*d21efd2eSMatthew G. Knepley ierr = VecViewFromOptions(u, NULL, "-interpolation_view");CHKERRQ(ierr); 164*d21efd2eSMatthew G. Knepley ierr = DMComputeL2Diff(dm, 0.0, exactSol, exactCtx, u, &error);CHKERRQ(ierr); 165*d21efd2eSMatthew G. Knepley ierr = DMRestoreGlobalVector(dm, &u);CHKERRQ(ierr); 166*d21efd2eSMatthew G. Knepley if (error > tol) {ierr = PetscPrintf(comm, "Interpolation tests FAIL at tolerance %g error %g\n", (double) tol, (double) error);CHKERRQ(ierr);} 167*d21efd2eSMatthew G. Knepley else {ierr = PetscPrintf(comm, "Interpolation tests pass at tolerance %g\n", (double) tol);CHKERRQ(ierr);} 168*d21efd2eSMatthew G. Knepley PetscFunctionReturn(0); 169*d21efd2eSMatthew G. Knepley } 170*d21efd2eSMatthew G. Knepley 171*d21efd2eSMatthew G. Knepley /* This test tells us whether the element is unisolvent (the mass matrix has full rank), and what rate of convergence we achieve */ 172*d21efd2eSMatthew G. Knepley static PetscErrorCode CheckL2Projection(DM dm, AppCtx *user) 173*d21efd2eSMatthew G. Knepley { 174*d21efd2eSMatthew G. Knepley PetscSimplePointFunc exactSol[1]; 175*d21efd2eSMatthew G. Knepley void *exactCtx[1]; 176*d21efd2eSMatthew G. Knepley SNES snes; 177*d21efd2eSMatthew G. Knepley PetscDS ds; 178*d21efd2eSMatthew G. Knepley Vec u; 179*d21efd2eSMatthew G. Knepley PetscReal error, tol = PETSC_SMALL; 180*d21efd2eSMatthew G. Knepley MPI_Comm comm; 181*d21efd2eSMatthew G. Knepley PetscErrorCode ierr; 182*d21efd2eSMatthew G. Knepley 183*d21efd2eSMatthew G. Knepley PetscFunctionBeginUser; 184*d21efd2eSMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 185*d21efd2eSMatthew G. Knepley ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 186*d21efd2eSMatthew G. Knepley ierr = DMGetGlobalVector(dm, &u);CHKERRQ(ierr); 187*d21efd2eSMatthew G. Knepley ierr = PetscDSGetExactSolution(ds, 0, &exactSol[0], &exactCtx[0]);CHKERRQ(ierr); 188*d21efd2eSMatthew G. Knepley ierr = SNESCreate(comm, &snes);CHKERRQ(ierr); 189*d21efd2eSMatthew G. Knepley ierr = SNESSetDM(snes, dm);CHKERRQ(ierr); 190*d21efd2eSMatthew G. Knepley ierr = VecSet(u, 0.0);CHKERRQ(ierr); 191*d21efd2eSMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) u, "solution");CHKERRQ(ierr); 192*d21efd2eSMatthew G. Knepley ierr = DMPlexSetSNESLocalFEM(dm, user, user, user);CHKERRQ(ierr); 193*d21efd2eSMatthew G. Knepley ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 194*d21efd2eSMatthew G. Knepley ierr = DMSNESCheckFromOptions(snes, u);CHKERRQ(ierr); 195*d21efd2eSMatthew G. Knepley ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr); 196*d21efd2eSMatthew G. Knepley ierr = SNESDestroy(&snes);CHKERRQ(ierr); 197*d21efd2eSMatthew G. Knepley ierr = VecViewFromOptions(u, NULL, "-l2_projection_view");CHKERRQ(ierr); 198*d21efd2eSMatthew G. Knepley ierr = DMComputeL2Diff(dm, 0.0, exactSol, exactCtx, u, &error);CHKERRQ(ierr); 199*d21efd2eSMatthew G. Knepley ierr = DMRestoreGlobalVector(dm, &u);CHKERRQ(ierr); 200*d21efd2eSMatthew G. Knepley if (error > tol) {ierr = PetscPrintf(comm, "L2 projection tests FAIL at tolerance %g error %g\n", (double) tol, (double) error);CHKERRQ(ierr);} 201*d21efd2eSMatthew G. Knepley else {ierr = PetscPrintf(comm, "L2 projection tests pass at tolerance %g\n", (double) tol);CHKERRQ(ierr);} 202*d21efd2eSMatthew G. Knepley PetscFunctionReturn(0); 203*d21efd2eSMatthew G. Knepley } 204*d21efd2eSMatthew G. Knepley 205*d21efd2eSMatthew G. Knepley int main(int argc, char **argv) 206*d21efd2eSMatthew G. Knepley { 207*d21efd2eSMatthew G. Knepley DM dm; 208*d21efd2eSMatthew G. Knepley AppCtx user; 209*d21efd2eSMatthew G. Knepley PetscMPIInt size; 210*d21efd2eSMatthew G. Knepley PetscErrorCode ierr; 211*d21efd2eSMatthew G. Knepley 212*d21efd2eSMatthew G. Knepley ierr = PetscInitialize(&argc, &argv, NULL, help); if (ierr) return ierr; 213*d21efd2eSMatthew G. Knepley ierr = MPI_Comm_size(PETSC_COMM_WORLD, &size);CHKERRMPI(ierr); 214*d21efd2eSMatthew G. Knepley if (size > 1) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "This is a uniprocessor example only."); 215*d21efd2eSMatthew G. Knepley ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 216*d21efd2eSMatthew G. Knepley ierr = DMCreate(PETSC_COMM_WORLD, &dm);CHKERRQ(ierr); 217*d21efd2eSMatthew G. Knepley ierr = DMSetType(dm, DMPLEX);CHKERRQ(ierr); 218*d21efd2eSMatthew G. Knepley ierr = DMSetFromOptions(dm);CHKERRQ(ierr); 219*d21efd2eSMatthew G. Knepley ierr = SetupDiscretization(dm, NULL, &user);CHKERRQ(ierr); 220*d21efd2eSMatthew G. Knepley 221*d21efd2eSMatthew G. Knepley ierr = CheckInterpolation(dm, &user);CHKERRQ(ierr); 222*d21efd2eSMatthew G. Knepley ierr = CheckL2Projection(dm, &user);CHKERRQ(ierr); 223*d21efd2eSMatthew G. Knepley 224*d21efd2eSMatthew G. Knepley ierr = DMDestroy(&dm);CHKERRQ(ierr); 225*d21efd2eSMatthew G. Knepley ierr = PetscFinalize(); 226*d21efd2eSMatthew G. Knepley return ierr; 227*d21efd2eSMatthew G. Knepley } 228*d21efd2eSMatthew G. Knepley 229*d21efd2eSMatthew G. Knepley /*TEST 230*d21efd2eSMatthew G. Knepley 231*d21efd2eSMatthew G. Knepley testset: 232*d21efd2eSMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -petscspace_degree 1\ 233*d21efd2eSMatthew G. Knepley -snes_error_if_not_converged -ksp_error_if_not_converged -pc_type lu 234*d21efd2eSMatthew G. Knepley 235*d21efd2eSMatthew G. Knepley test: 236*d21efd2eSMatthew G. Knepley suffix: p1_0 237*d21efd2eSMatthew G. Knepley args: -func {{constant linear}} 238*d21efd2eSMatthew G. Knepley 239*d21efd2eSMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 4 gives convergence rate 2.0 240*d21efd2eSMatthew G. Knepley test: 241*d21efd2eSMatthew G. Knepley suffix: p1_1 242*d21efd2eSMatthew G. Knepley args: -func {{quadratic trig}} \ 243*d21efd2eSMatthew G. Knepley -snes_convergence_estimate -convest_num_refine 2 244*d21efd2eSMatthew G. Knepley 245*d21efd2eSMatthew G. Knepley testset: 246*d21efd2eSMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism \ 247*d21efd2eSMatthew G. Knepley -petscspace_type sum \ 248*d21efd2eSMatthew G. Knepley -petscspace_variables 3 \ 249*d21efd2eSMatthew G. Knepley -petscspace_components 3 \ 250*d21efd2eSMatthew G. Knepley -petscspace_sum_spaces 2 \ 251*d21efd2eSMatthew G. Knepley -petscspace_sum_concatenate false \ 252*d21efd2eSMatthew G. Knepley -sumcomp_0_petscspace_variables 3 \ 253*d21efd2eSMatthew G. Knepley -sumcomp_0_petscspace_components 3 \ 254*d21efd2eSMatthew G. Knepley -sumcomp_0_petscspace_degree 1 \ 255*d21efd2eSMatthew G. Knepley -sumcomp_1_petscspace_variables 3 \ 256*d21efd2eSMatthew G. Knepley -sumcomp_1_petscspace_components 3 \ 257*d21efd2eSMatthew G. Knepley -sumcomp_1_petscspace_type wxy \ 258*d21efd2eSMatthew G. Knepley -petscdualspace_form_degree 0 \ 259*d21efd2eSMatthew G. Knepley -petscdualspace_order 1 \ 260*d21efd2eSMatthew G. Knepley -petscdualspace_components 3 \ 261*d21efd2eSMatthew G. Knepley -snes_error_if_not_converged -ksp_error_if_not_converged -pc_type lu 262*d21efd2eSMatthew G. Knepley 263*d21efd2eSMatthew G. Knepley test: 264*d21efd2eSMatthew G. Knepley suffix: wxy_0 265*d21efd2eSMatthew G. Knepley args: -func constant 266*d21efd2eSMatthew G. Knepley 267*d21efd2eSMatthew G. Knepley test: 268*d21efd2eSMatthew G. Knepley suffix: wxy_1 269*d21efd2eSMatthew G. Knepley args: -func linear 270*d21efd2eSMatthew G. Knepley 271*d21efd2eSMatthew G. Knepley TEST*/ 272