xref: /petsc/src/dm/dt/fe/tests/ex3.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1d21efd2eSMatthew G. Knepley static const char help[] = "Tests for determining whether a new finite element works";
2d21efd2eSMatthew G. Knepley 
3d21efd2eSMatthew G. Knepley /*
4d21efd2eSMatthew G. Knepley   Use -interpolation_view and -l2_projection_view to look at the interpolants.
5d21efd2eSMatthew G. Knepley */
6d21efd2eSMatthew G. Knepley 
7d21efd2eSMatthew G. Knepley #include <petscdmplex.h>
8d21efd2eSMatthew G. Knepley #include <petscfe.h>
9d21efd2eSMatthew G. Knepley #include <petscds.h>
10d21efd2eSMatthew G. Knepley #include <petscsnes.h>
11d21efd2eSMatthew G. Knepley 
12*9371c9d4SSatish Balay static void constant(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[]) {
13d21efd2eSMatthew G. Knepley   const PetscInt Nc = uOff[1] - uOff[0];
14d21efd2eSMatthew G. Knepley   for (PetscInt c = 0; c < Nc; ++c) f0[c] += 5.;
15d21efd2eSMatthew G. Knepley }
16d21efd2eSMatthew G. Knepley 
17*9371c9d4SSatish Balay static void linear(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[]) {
18d21efd2eSMatthew G. Knepley   const PetscInt Nc = uOff[1] - uOff[0];
19d21efd2eSMatthew G. Knepley   for (PetscInt c = 0; c < Nc; ++c) f0[c] += 5. * x[c];
20d21efd2eSMatthew G. Knepley }
21d21efd2eSMatthew G. Knepley 
22*9371c9d4SSatish Balay static void quadratic(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[]) {
23d21efd2eSMatthew G. Knepley   const PetscInt Nc = uOff[1] - uOff[0];
24d21efd2eSMatthew G. Knepley   for (PetscInt c = 0; c < Nc; ++c) f0[c] += 5. * x[c] * x[c];
25d21efd2eSMatthew G. Knepley }
26d21efd2eSMatthew G. Knepley 
27*9371c9d4SSatish Balay static void trig(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[]) {
28d21efd2eSMatthew G. Knepley   const PetscInt Nc = uOff[1] - uOff[0];
29d21efd2eSMatthew G. Knepley   for (PetscInt c = 0; c < Nc; ++c) f0[c] += PetscCosReal(2. * PETSC_PI * x[c]);
30d21efd2eSMatthew G. Knepley }
31d21efd2eSMatthew G. Knepley 
32e239af90SMatthew G. Knepley /*
33e239af90SMatthew G. Knepley  The prime basis for the Wheeler-Yotov-Xue prism.
34e239af90SMatthew G. Knepley  */
35*9371c9d4SSatish Balay static void prime(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[]) {
36e239af90SMatthew G. Knepley   PetscReal x = X[0], y = X[1], z = X[2], b = 1 + x + y + z;
37e239af90SMatthew G. Knepley   f0[0] += b + 2.0 * x * z + 2.0 * y * z + x * y + x * x;
38e239af90SMatthew G. Knepley   f0[1] += b + 2.0 * x * z + 2.0 * y * z + x * y + y * y;
39e239af90SMatthew G. Knepley   f0[2] += b - 3.0 * x * z - 3.0 * y * z - 2.0 * z * z;
40e239af90SMatthew G. Knepley }
41e239af90SMatthew G. Knepley 
42e239af90SMatthew G. Knepley static const char    *names[]     = {"constant", "linear", "quadratic", "trig", "prime"};
43e239af90SMatthew G. Knepley static PetscPointFunc functions[] = {constant, linear, quadratic, trig, prime};
44d21efd2eSMatthew G. Knepley 
45d21efd2eSMatthew G. Knepley typedef struct {
46d21efd2eSMatthew G. Knepley   PetscPointFunc exactSol;
47e239af90SMatthew G. Knepley   PetscReal      shear, flatten;
48d21efd2eSMatthew G. Knepley } AppCtx;
49d21efd2eSMatthew G. Knepley 
50*9371c9d4SSatish Balay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) {
51d21efd2eSMatthew G. Knepley   char     name[PETSC_MAX_PATH_LEN] = "constant";
52dd39110bSPierre Jolivet   PetscInt Nfunc                    = PETSC_STATIC_ARRAY_LENGTH(names), i;
53d21efd2eSMatthew G. Knepley 
54d21efd2eSMatthew G. Knepley   PetscFunctionBeginUser;
55d21efd2eSMatthew G. Knepley   options->exactSol = NULL;
56e239af90SMatthew G. Knepley   options->shear    = 0.;
57e239af90SMatthew G. Knepley   options->flatten  = 1.;
58d21efd2eSMatthew G. Knepley 
59d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "FE Test Options", "PETSCFE");
609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-func", "Function to project into space", "", name, name, PETSC_MAX_PATH_LEN, NULL));
619566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-shear", "Factor by which to shear along the x-direction", "", options->shear, &(options->shear), NULL));
629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-flatten", "Factor by which to flatten", "", options->flatten, &(options->flatten), NULL));
63d0609cedSBarry Smith   PetscOptionsEnd();
64d21efd2eSMatthew G. Knepley 
65d21efd2eSMatthew G. Knepley   for (i = 0; i < Nfunc; ++i) {
66d21efd2eSMatthew G. Knepley     PetscBool flg;
67d21efd2eSMatthew G. Knepley 
689566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(name, names[i], &flg));
69*9371c9d4SSatish Balay     if (flg) {
70*9371c9d4SSatish Balay       options->exactSol = functions[i];
71*9371c9d4SSatish Balay       break;
72*9371c9d4SSatish Balay     }
73d21efd2eSMatthew G. Knepley   }
74e239af90SMatthew G. Knepley   PetscCheck(options->exactSol, comm, PETSC_ERR_ARG_WRONG, "Invalid test function %s", name);
75d21efd2eSMatthew G. Knepley   PetscFunctionReturn(0);
76d21efd2eSMatthew G. Knepley }
77d21efd2eSMatthew G. Knepley 
78d21efd2eSMatthew G. Knepley /* The exact solution is the negative of the f0 contribution */
79*9371c9d4SSatish Balay static PetscErrorCode exactSolution(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
80d21efd2eSMatthew G. Knepley   AppCtx  *user    = (AppCtx *)ctx;
81e239af90SMatthew G. Knepley   PetscInt uOff[2] = {0, Nc};
82d21efd2eSMatthew G. Knepley 
83d21efd2eSMatthew G. Knepley   user->exactSol(dim, 1, 0, uOff, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, time, x, 0, NULL, u);
84d21efd2eSMatthew G. Knepley   for (PetscInt c = 0; c < Nc; ++c) u[c] *= -1.;
85d21efd2eSMatthew G. Knepley   return 0;
86d21efd2eSMatthew G. Knepley }
87d21efd2eSMatthew G. Knepley 
88*9371c9d4SSatish Balay static void f0(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[]) {
89d21efd2eSMatthew G. Knepley   const PetscInt Nc = uOff[1] - uOff[0];
90d21efd2eSMatthew G. Knepley   for (PetscInt c = 0; c < Nc; ++c) f0[c] += u[c];
91d21efd2eSMatthew G. Knepley }
92d21efd2eSMatthew G. Knepley 
93*9371c9d4SSatish Balay static void g0(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[]) {
94d21efd2eSMatthew G. Knepley   const PetscInt Nc = uOff[1] - uOff[0];
95d21efd2eSMatthew G. Knepley   for (PetscInt c = 0; c < Nc; ++c) g0[c * Nc + c] = 1.0;
96d21efd2eSMatthew G. Knepley }
97d21efd2eSMatthew G. Knepley 
98*9371c9d4SSatish Balay static PetscErrorCode SetupProblem(DM dm, AppCtx *user) {
99d21efd2eSMatthew G. Knepley   PetscDS       ds;
100d21efd2eSMatthew G. Knepley   PetscWeakForm wf;
101d21efd2eSMatthew G. Knepley 
102d21efd2eSMatthew G. Knepley   PetscFunctionBeginUser;
1039566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
1049566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakForm(ds, &wf));
1059566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexResidual(wf, NULL, 0, 0, 0, 0, f0, 0, NULL));
1069566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexResidual(wf, NULL, 0, 0, 0, 1, user->exactSol, 0, NULL));
1079566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexJacobian(wf, NULL, 0, 0, 0, 0, 0, g0, 0, NULL, 0, NULL, 0, NULL));
1089566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolution(ds, 0, exactSolution, user));
109d21efd2eSMatthew G. Knepley   PetscFunctionReturn(0);
110d21efd2eSMatthew G. Knepley }
111d21efd2eSMatthew G. Knepley 
112*9371c9d4SSatish Balay static PetscErrorCode SetupDiscretization(DM dm, const char name[], AppCtx *user) {
113d21efd2eSMatthew G. Knepley   DM      cdm = dm;
114d21efd2eSMatthew G. Knepley   PetscFE fe;
115d21efd2eSMatthew G. Knepley   char    prefix[PETSC_MAX_PATH_LEN];
116d21efd2eSMatthew G. Knepley 
117d21efd2eSMatthew G. Knepley   PetscFunctionBeginUser;
11863a3b9bcSJacob Faibussowitsch   if (name) PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%s_", name));
1199566063dSJacob Faibussowitsch   PetscCall(DMCreateFEDefault(dm, 1, name ? prefix : NULL, -1, &fe));
1209566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)fe, name ? name : "Solution"));
121d21efd2eSMatthew G. Knepley   /* Set discretization and boundary conditions for each mesh */
1229566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
1239566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
1249566063dSJacob Faibussowitsch   PetscCall(SetupProblem(dm, user));
125d21efd2eSMatthew G. Knepley   while (cdm) {
1269566063dSJacob Faibussowitsch     PetscCall(DMCopyDisc(dm, cdm));
1279566063dSJacob Faibussowitsch     PetscCall(DMGetCoarseDM(cdm, &cdm));
128d21efd2eSMatthew G. Knepley   }
1299566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
130d21efd2eSMatthew G. Knepley   PetscFunctionReturn(0);
131d21efd2eSMatthew G. Knepley }
132d21efd2eSMatthew G. Knepley 
133d21efd2eSMatthew G. Knepley /* This test tells us whether the given function is contained in the approximation space */
134*9371c9d4SSatish Balay static PetscErrorCode CheckInterpolation(DM dm, AppCtx *user) {
135d21efd2eSMatthew G. Knepley   PetscSimplePointFunc exactSol[1];
136d21efd2eSMatthew G. Knepley   void                *exactCtx[1];
137d21efd2eSMatthew G. Knepley   PetscDS              ds;
138d21efd2eSMatthew G. Knepley   Vec                  u;
139d21efd2eSMatthew G. Knepley   PetscReal            error, tol = PETSC_SMALL;
140d21efd2eSMatthew G. Knepley   MPI_Comm             comm;
141d21efd2eSMatthew G. Knepley 
142d21efd2eSMatthew G. Knepley   PetscFunctionBeginUser;
1439566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
1449566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
1459566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dm, &u));
1469566063dSJacob Faibussowitsch   PetscCall(PetscDSGetExactSolution(ds, 0, &exactSol[0], &exactCtx[0]));
1479566063dSJacob Faibussowitsch   PetscCall(DMProjectFunction(dm, 0.0, exactSol, exactCtx, INSERT_ALL_VALUES, u));
1489566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(u, NULL, "-interpolation_view"));
1499566063dSJacob Faibussowitsch   PetscCall(DMComputeL2Diff(dm, 0.0, exactSol, exactCtx, u, &error));
1509566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dm, &u));
1519566063dSJacob Faibussowitsch   if (error > tol) PetscCall(PetscPrintf(comm, "Interpolation tests FAIL at tolerance %g error %g\n", (double)tol, (double)error));
1529566063dSJacob Faibussowitsch   else PetscCall(PetscPrintf(comm, "Interpolation tests pass at tolerance %g\n", (double)tol));
153d21efd2eSMatthew G. Knepley   PetscFunctionReturn(0);
154d21efd2eSMatthew G. Knepley }
155d21efd2eSMatthew G. Knepley 
156d21efd2eSMatthew G. Knepley /* This test tells us whether the element is unisolvent (the mass matrix has full rank), and what rate of convergence we achieve */
157*9371c9d4SSatish Balay static PetscErrorCode CheckL2Projection(DM dm, AppCtx *user) {
158d21efd2eSMatthew G. Knepley   PetscSimplePointFunc exactSol[1];
159d21efd2eSMatthew G. Knepley   void                *exactCtx[1];
160d21efd2eSMatthew G. Knepley   SNES                 snes;
161d21efd2eSMatthew G. Knepley   PetscDS              ds;
162d21efd2eSMatthew G. Knepley   Vec                  u;
163d21efd2eSMatthew G. Knepley   PetscReal            error, tol = PETSC_SMALL;
164d21efd2eSMatthew G. Knepley   MPI_Comm             comm;
165d21efd2eSMatthew G. Knepley 
166d21efd2eSMatthew G. Knepley   PetscFunctionBeginUser;
1679566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
1689566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
1699566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dm, &u));
1709566063dSJacob Faibussowitsch   PetscCall(PetscDSGetExactSolution(ds, 0, &exactSol[0], &exactCtx[0]));
1719566063dSJacob Faibussowitsch   PetscCall(SNESCreate(comm, &snes));
1729566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes, dm));
1739566063dSJacob Faibussowitsch   PetscCall(VecSet(u, 0.0));
1749566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)u, "solution"));
1759566063dSJacob Faibussowitsch   PetscCall(DMPlexSetSNESLocalFEM(dm, user, user, user));
1769566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
1779566063dSJacob Faibussowitsch   PetscCall(DMSNESCheckFromOptions(snes, u));
1789566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, NULL, u));
1799566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
1809566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(u, NULL, "-l2_projection_view"));
1819566063dSJacob Faibussowitsch   PetscCall(DMComputeL2Diff(dm, 0.0, exactSol, exactCtx, u, &error));
1829566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dm, &u));
1839566063dSJacob Faibussowitsch   if (error > tol) PetscCall(PetscPrintf(comm, "L2 projection tests FAIL at tolerance %g error %g\n", (double)tol, (double)error));
1849566063dSJacob Faibussowitsch   else PetscCall(PetscPrintf(comm, "L2 projection tests pass at tolerance %g\n", (double)tol));
185d21efd2eSMatthew G. Knepley   PetscFunctionReturn(0);
186d21efd2eSMatthew G. Knepley }
187d21efd2eSMatthew G. Knepley 
188e239af90SMatthew G. Knepley /* Distorts the mesh by shearing in the x-direction and flattening, factors provided in the options. */
189*9371c9d4SSatish Balay static PetscErrorCode DistortMesh(DM dm, AppCtx *user) {
190e239af90SMatthew G. Knepley   Vec          coordinates;
191e239af90SMatthew G. Knepley   PetscScalar *ca;
192e239af90SMatthew G. Knepley   PetscInt     dE, n, i;
193e239af90SMatthew G. Knepley 
194e239af90SMatthew G. Knepley   PetscFunctionBeginUser;
1959566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &dE));
1969566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinates(dm, &coordinates));
1979566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(coordinates, &n));
1989566063dSJacob Faibussowitsch   PetscCall(VecGetArray(coordinates, &ca));
199e239af90SMatthew G. Knepley   for (i = 0; i < (n / dE); ++i) {
200e239af90SMatthew G. Knepley     ca[i * dE + 0] += user->shear * ca[i * dE + 0];
201e239af90SMatthew G. Knepley     ca[i * dE + 1] *= user->flatten;
202e239af90SMatthew G. Knepley   }
2039566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(coordinates, &ca));
204e239af90SMatthew G. Knepley   PetscFunctionReturn(0);
205e239af90SMatthew G. Knepley }
206e239af90SMatthew G. Knepley 
207*9371c9d4SSatish Balay int main(int argc, char **argv) {
208d21efd2eSMatthew G. Knepley   DM          dm;
209d21efd2eSMatthew G. Knepley   AppCtx      user;
210d21efd2eSMatthew G. Knepley   PetscMPIInt size;
211d21efd2eSMatthew G. Knepley 
212327415f7SBarry Smith   PetscFunctionBeginUser;
2139566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
2149566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
215e239af90SMatthew G. Knepley   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_SUP, "This is a uniprocessor example only.");
2169566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
2179566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
2189566063dSJacob Faibussowitsch   PetscCall(DMSetType(dm, DMPLEX));
2199566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
2209566063dSJacob Faibussowitsch   PetscCall(DistortMesh(dm, &user));
2219566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
2229566063dSJacob Faibussowitsch   PetscCall(SetupDiscretization(dm, NULL, &user));
223d21efd2eSMatthew G. Knepley 
2249566063dSJacob Faibussowitsch   PetscCall(CheckInterpolation(dm, &user));
2259566063dSJacob Faibussowitsch   PetscCall(CheckL2Projection(dm, &user));
226d21efd2eSMatthew G. Knepley 
2279566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
2289566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
229b122ec5aSJacob Faibussowitsch   return 0;
230d21efd2eSMatthew G. Knepley }
231d21efd2eSMatthew G. Knepley 
232d21efd2eSMatthew G. Knepley /*TEST
233d21efd2eSMatthew G. Knepley 
234d21efd2eSMatthew G. Knepley   testset:
235d21efd2eSMatthew G. Knepley     args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -petscspace_degree 1\
236d21efd2eSMatthew G. Knepley           -snes_error_if_not_converged -ksp_error_if_not_converged -pc_type lu
237d21efd2eSMatthew G. Knepley 
238d21efd2eSMatthew G. Knepley     test:
239d21efd2eSMatthew G. Knepley       suffix: p1_0
240d21efd2eSMatthew G. Knepley       args: -func {{constant linear}}
241d21efd2eSMatthew G. Knepley 
242d21efd2eSMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 4 gives convergence rate 2.0
243d21efd2eSMatthew G. Knepley     test:
244d21efd2eSMatthew G. Knepley       suffix: p1_1
245d21efd2eSMatthew G. Knepley       args: -func {{quadratic trig}} \
246d21efd2eSMatthew G. Knepley             -snes_convergence_estimate -convest_num_refine 2
247d21efd2eSMatthew G. Knepley 
248d21efd2eSMatthew G. Knepley   testset:
249e239af90SMatthew G. Knepley     requires: !complex double
250d21efd2eSMatthew G. Knepley     args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism \
251d21efd2eSMatthew G. Knepley             -petscspace_type sum \
252d21efd2eSMatthew G. Knepley             -petscspace_variables 3 \
253d21efd2eSMatthew G. Knepley             -petscspace_components 3 \
254d21efd2eSMatthew G. Knepley             -petscspace_sum_spaces 2 \
255d21efd2eSMatthew G. Knepley             -petscspace_sum_concatenate false \
256d21efd2eSMatthew G. Knepley               -sumcomp_0_petscspace_variables 3 \
257d21efd2eSMatthew G. Knepley               -sumcomp_0_petscspace_components 3 \
258d21efd2eSMatthew G. Knepley               -sumcomp_0_petscspace_degree 1 \
259d21efd2eSMatthew G. Knepley               -sumcomp_1_petscspace_variables 3 \
260d21efd2eSMatthew G. Knepley               -sumcomp_1_petscspace_components 3 \
261d21efd2eSMatthew G. Knepley               -sumcomp_1_petscspace_type wxy \
262d21efd2eSMatthew G. Knepley             -petscdualspace_form_degree 0 \
263d21efd2eSMatthew G. Knepley             -petscdualspace_order 1 \
264d21efd2eSMatthew G. Knepley             -petscdualspace_components 3 \
265d21efd2eSMatthew G. Knepley           -snes_error_if_not_converged -ksp_error_if_not_converged -pc_type lu
266d21efd2eSMatthew G. Knepley 
267d21efd2eSMatthew G. Knepley     test:
268d21efd2eSMatthew G. Knepley       suffix: wxy_0
269d21efd2eSMatthew G. Knepley       args: -func constant
270d21efd2eSMatthew G. Knepley 
271d21efd2eSMatthew G. Knepley     test:
272d21efd2eSMatthew G. Knepley       suffix: wxy_1
273d21efd2eSMatthew G. Knepley       args: -func linear
274d21efd2eSMatthew G. Knepley 
275e239af90SMatthew G. Knepley     test:
276e239af90SMatthew G. Knepley       suffix: wxy_2
277e239af90SMatthew G. Knepley       args: -func prime
278e239af90SMatthew G. Knepley 
279e239af90SMatthew G. Knepley     test:
280e239af90SMatthew G. Knepley       suffix: wxy_3
281e239af90SMatthew G. Knepley       args: -func linear -shear 1 -flatten 1e-5
282e239af90SMatthew G. Knepley 
283d21efd2eSMatthew G. Knepley TEST*/
284