xref: /petsc/src/snes/tutorials/ex24.c (revision 45480ffeba491aca5d3f3b8c78679069fb317d32)
1 static char help[] = "Poisson Problem in mixed form with 2d and 3d with finite elements.\n\
2 We solve the Poisson problem in a rectangular\n\
3 domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\
4 This example supports automatic convergence estimation\n\
5 and Hdiv elements.\n\n\n";
6 
7 /*
8 
9 The mixed form of Poisson's equation, e.g. https://firedrakeproject.org/demos/poisson_mixed.py.html, is given
10 in the strong form by
11 \begin{align}
12   \vb{q} - \nabla u   &= 0 \\
13   \nabla \cdot \vb{q} &= f
14 \end{align}
15 where $u$ is the potential, as in the original problem, but we also discretize the gradient of potential $\vb{q}$,
16 or flux, directly. The weak form is then
17 \begin{align}
18   <t, \vb{q}> + <\nabla \cdot t, u> - <t_n, u>_\Gamma &= 0 \\
19   <v, \nabla \cdot \vb{q}> &= <v, f>
20 \end{align}
21 For the original Poisson problem, the Dirichlet boundary forces an essential boundary condition on the potential space,
22 and the Neumann boundary gives a natural boundary condition in the weak form. In the mixed formulation, the Neumann
23 boundary gives an essential boundary condition on the flux space, $\vb{q} \cdot \vu{n} = h$, and the Dirichlet condition
24 becomes a natural condition in the weak form, <t_n, g>_\Gamma.
25 */
26 
27 #include <petscdmplex.h>
28 #include <petscsnes.h>
29 #include <petscds.h>
30 #include <petscconvest.h>
31 
32 typedef enum {SOL_LINEAR, SOL_QUADRATIC, SOL_QUARTIC, SOL_QUARTIC_NEUMANN, SOL_UNKNOWN, NUM_SOLTYPE} SolType;
33 const char *SolTypeNames[NUM_SOLTYPE+3] = {"linear", "quadratic", "quartic", "quartic_neumann", "unknown", "SolType", "SOL_", NULL};
34 
35 typedef struct {
36   SolType solType; /* The type of exact solution */
37 } AppCtx;
38 
39 static void f0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
40                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
41                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
42                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
43 {
44   PetscInt d;
45   for (d = 0; d < dim; ++d) f0[0] += u_x[uOff_x[0]+d*dim+d];
46 }
47 
48 /* 2D Dirichlet potential example
49 
50   u = x
51   q = <1, 0>
52   f = 0
53 
54   We will need a boundary integral of u over \Gamma.
55 */
56 static PetscErrorCode linear_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
57 {
58   u[0] = x[0];
59   return 0;
60 }
61 static PetscErrorCode linear_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
62 {
63   PetscInt c;
64   for (c = 0; c < Nc; ++c) u[c] = c ? 0.0 : 1.0;
65   return 0;
66 }
67 
68 static void f0_linear_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
69                         const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
70                         const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
71                         PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
72 {
73   f0[0] = 0.0;
74   f0_u(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, a_t, a_x, t, x, numConstants, constants, f0);
75 }
76 static void f0_bd_linear_q(PetscInt dim, PetscInt Nf, PetscInt NfAux,
77                            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
78                            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
79                            PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
80 {
81   PetscScalar potential;
82   PetscInt    d;
83 
84   linear_u(dim, t, x, dim, &potential, NULL);
85   for (d = 0; d < dim; ++d) f0[d] = -potential*n[d];
86 }
87 
88 /* 2D Dirichlet potential example
89 
90   u = x^2 + y^2
91   q = <2x, 2y>
92   f = 4
93 
94   We will need a boundary integral of u over \Gamma.
95 */
96 static PetscErrorCode quadratic_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
97 {
98   PetscInt d;
99 
100   u[0] = 0.0;
101   for (d = 0; d < dim; ++d) u[0] += x[d]*x[d];
102   return 0;
103 }
104 static PetscErrorCode quadratic_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
105 {
106   PetscInt c;
107   for (c = 0; c < Nc; ++c) u[c] = 2.0*x[c];
108   return 0;
109 }
110 
111 static void f0_quadratic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
112                            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
113                            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
114                            PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
115 {
116   f0[0] = -4.0;
117   f0_u(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, a_t, a_x, t, x, numConstants, constants, f0);
118 }
119 static void f0_bd_quadratic_q(PetscInt dim, PetscInt Nf, PetscInt NfAux,
120                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
121                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
122                               PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
123 {
124   PetscScalar potential;
125   PetscInt    d;
126 
127   quadratic_u(dim, t, x, dim, &potential, NULL);
128   for (d = 0; d < dim; ++d) f0[d] = -potential*n[d];
129 }
130 
131 /* 2D Dirichlet potential example
132 
133   u = x (1-x) y (1-y)
134   q = <(1-2x) y (1-y), x (1-x) (1-2y)>
135   f = -y (1-y) - x (1-x)
136 
137   u|_\Gamma = 0 so that the boundary integral vanishes
138 */
139 static PetscErrorCode quartic_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
140 {
141   PetscInt d;
142 
143   u[0] = 1.0;
144   for (d = 0; d < dim; ++d) u[0] *= x[d]*(1.0 - x[d]);
145   return 0;
146 }
147 static PetscErrorCode quartic_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
148 {
149   PetscInt c, d;
150 
151   for (c = 0; c < Nc; ++c) {
152     u[c] = 1.0;
153     for (d = 0; d < dim; ++d) {
154       if (c == d) u[c] *= 1 - 2.0*x[d];
155       else        u[c] *= x[d]*(1.0 - x[d]);
156     }
157   }
158   return 0;
159 }
160 
161 static void f0_quartic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
162                         const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
163                         const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
164                         PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
165 {
166   PetscInt d;
167   f0[0] = 0.0;
168   for (d = 0; d < dim; ++d) f0[0] += 2.0*x[d]*(1.0 - x[d]);
169   f0_u(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, a_t, a_x, t, x, numConstants, constants, f0);
170 }
171 
172 /* 2D Dirichlet potential example
173 
174   u = x (1-x) y (1-y)
175   q = <(1-2x) y (1-y), x (1-x) (1-2y)>
176   f = -y (1-y) - x (1-x)
177 
178   du/dn_\Gamma = {(1-2x) y (1-y), x (1-x) (1-2y)} . n produces an essential condition on q
179 */
180 
181 static void f0_q(PetscInt dim, PetscInt Nf, PetscInt NfAux,
182                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
183                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
184                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
185 {
186   PetscInt c;
187   for (c = 0; c < dim; ++c) f0[c] = u[uOff[0]+c];
188 }
189 
190 /* <\nabla\cdot w, u> = <\nabla w, Iu> */
191 static void f1_q(PetscInt dim, PetscInt Nf, PetscInt NfAux,
192                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
193                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
194                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
195 {
196   PetscInt c;
197   for (c = 0; c < dim; ++c) f1[c*dim+c] = u[uOff[1]];
198 }
199 
200 /* <t, q> */
201 static void g0_qq(PetscInt dim, PetscInt Nf, PetscInt NfAux,
202                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
203                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
204                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
205 {
206   PetscInt c;
207   for (c = 0; c < dim; ++c) g0[c*dim+c] = 1.0;
208 }
209 
210 /* <\nabla\cdot t, u> = <\nabla t, Iu> */
211 static void g2_qu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
212                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
213                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
214                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
215 {
216   PetscInt d;
217   for (d = 0; d < dim; ++d) g2[d*dim+d] = 1.0;
218 }
219 
220 /* <v, \nabla\cdot q> */
221 static void g1_uq(PetscInt dim, PetscInt Nf, PetscInt NfAux,
222                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
223                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
224                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
225 {
226   PetscInt d;
227   for (d = 0; d < dim; ++d) g1[d*dim+d] = 1.0;
228 }
229 
230 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
231 {
232   PetscErrorCode ierr;
233 
234   PetscFunctionBeginUser;
235   options->solType = SOL_LINEAR;
236 
237   ierr = PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX");CHKERRQ(ierr);
238   ierr = PetscOptionsEnum("-sol_type", "Type of exact solution", "ex24.c", SolTypeNames, (PetscEnum) options->solType, (PetscEnum *) &options->solType, NULL);CHKERRQ(ierr);
239   ierr = PetscOptionsEnd();
240   PetscFunctionReturn(0);
241 }
242 
243 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
244 {
245   PetscErrorCode ierr;
246 
247   PetscFunctionBeginUser;
248   if (0) {
249     DMLabel     label;
250     const char *name = "marker";
251 
252     ierr = DMPlexCreateReferenceCell(comm, 2, PETSC_TRUE, dm);CHKERRQ(ierr);
253     ierr = DMCreateLabel(*dm, name);CHKERRQ(ierr);
254     ierr = DMGetLabel(*dm, name, &label);CHKERRQ(ierr);
255     ierr = DMPlexMarkBoundaryFaces(*dm, 1, label);CHKERRQ(ierr);
256     ierr = DMPlexLabelComplete(*dm, label);CHKERRQ(ierr);
257   } else {
258     ierr = DMPlexCreateBoxMesh(comm, 2, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);
259   }
260   ierr = DMLocalizeCoordinates(*dm);CHKERRQ(ierr);
261   ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr);
262   ierr = DMSetApplicationContext(*dm, user);CHKERRQ(ierr);
263   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
264   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
265   PetscFunctionReturn(0);
266 }
267 
268 static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user)
269 {
270   PetscDS        ds;
271   DMLabel        label;
272   PetscWeakForm  wf;
273   const PetscInt id = 1;
274   PetscInt       bd;
275   PetscErrorCode ierr;
276 
277   PetscFunctionBeginUser;
278   ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr);
279   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
280   ierr = PetscDSSetResidual(ds, 0, f0_q, f1_q);CHKERRQ(ierr);
281   ierr = PetscDSSetJacobian(ds, 0, 0, g0_qq, NULL, NULL, NULL);CHKERRQ(ierr);
282   ierr = PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_qu, NULL);CHKERRQ(ierr);
283   ierr = PetscDSSetJacobian(ds, 1, 0, NULL, g1_uq, NULL, NULL);CHKERRQ(ierr);
284   switch (user->solType)
285   {
286     case SOL_LINEAR:
287       ierr = PetscDSSetResidual(ds, 1, f0_linear_u, NULL);CHKERRQ(ierr);
288       ierr = DMAddBoundary(dm, DM_BC_NATURAL, "Dirichlet Bd Integral", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd);CHKERRQ(ierr);
289       ierr = PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
290       ierr = PetscWeakFormSetIndexBdResidual(wf, label, 1, 0, 0, f0_bd_linear_q, 0, NULL);CHKERRQ(ierr);
291       ierr = PetscDSSetExactSolution(ds, 0, linear_q, user);CHKERRQ(ierr);
292       ierr = PetscDSSetExactSolution(ds, 1, linear_u, user);CHKERRQ(ierr);
293       break;
294     case SOL_QUADRATIC:
295       ierr = PetscDSSetResidual(ds, 1, f0_quadratic_u, NULL);CHKERRQ(ierr);
296       ierr = DMAddBoundary(dm, DM_BC_NATURAL, "Dirichlet Bd Integral", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd);CHKERRQ(ierr);
297       ierr = PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
298       ierr = PetscWeakFormSetIndexBdResidual(wf, label, 1, 0, 0, f0_bd_quadratic_q, 0, NULL);CHKERRQ(ierr);
299       ierr = PetscDSSetExactSolution(ds, 0, quadratic_q, user);CHKERRQ(ierr);
300       ierr = PetscDSSetExactSolution(ds, 1, quadratic_u, user);CHKERRQ(ierr);
301       break;
302     case SOL_QUARTIC:
303       ierr = PetscDSSetResidual(ds, 1, f0_quartic_u, NULL);CHKERRQ(ierr);
304       ierr = PetscDSSetExactSolution(ds, 0, quartic_q, user);CHKERRQ(ierr);
305       ierr = PetscDSSetExactSolution(ds, 1, quartic_u, user);CHKERRQ(ierr);
306       break;
307     case SOL_QUARTIC_NEUMANN:
308       ierr = PetscDSSetResidual(ds, 1, f0_quartic_u, NULL);CHKERRQ(ierr);
309       ierr = PetscDSSetExactSolution(ds, 0, quartic_q, user);CHKERRQ(ierr);
310       ierr = PetscDSSetExactSolution(ds, 1, quartic_u, user);CHKERRQ(ierr);
311       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "Flux condition", label, 1, &id, 0, 0, NULL, (void (*)(void)) quartic_q, NULL, user, NULL);CHKERRQ(ierr);
312       break;
313     default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Invalid exact solution type %s", SolTypeNames[PetscMin(user->solType, SOL_UNKNOWN)]);
314   }
315   PetscFunctionReturn(0);
316 }
317 
318 static PetscErrorCode SetupDiscretization(DM dm, PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user)
319 {
320   DM             cdm = dm;
321   PetscFE        feq, feu;
322   DMPolytopeType ct;
323   PetscBool      simplex;
324   PetscInt       dim, cStart;
325   PetscErrorCode ierr;
326 
327   PetscFunctionBeginUser;
328   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
329   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr);
330   ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr);
331   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
332   /* Create finite element */
333   ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, dim, simplex, "field_",     -1, &feq);CHKERRQ(ierr);
334   ierr = PetscObjectSetName((PetscObject) feq, "field");CHKERRQ(ierr);
335   ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, 1,   simplex, "potential_", -1, &feu);CHKERRQ(ierr);
336   ierr = PetscObjectSetName((PetscObject) feu, "potential");CHKERRQ(ierr);
337   ierr = PetscFECopyQuadrature(feq, feu);CHKERRQ(ierr);
338   /* Set discretization and boundary conditions for each mesh */
339   ierr = DMSetField(dm, 0, NULL, (PetscObject) feq);CHKERRQ(ierr);
340   ierr = DMSetField(dm, 1, NULL, (PetscObject) feu);CHKERRQ(ierr);
341   ierr = DMCreateDS(dm);CHKERRQ(ierr);
342   ierr = (*setup)(dm, user);CHKERRQ(ierr);
343   while (cdm) {
344     ierr = DMCopyDisc(dm,cdm);CHKERRQ(ierr);
345     ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
346   }
347   ierr = PetscFEDestroy(&feq);CHKERRQ(ierr);
348   ierr = PetscFEDestroy(&feu);CHKERRQ(ierr);
349   PetscFunctionReturn(0);
350 }
351 
352 int main(int argc, char **argv)
353 {
354   DM             dm;   /* Problem specification */
355   SNES           snes; /* Nonlinear solver */
356   Vec            u;    /* Solutions */
357   AppCtx         user; /* User-defined work context */
358   PetscErrorCode ierr;
359 
360   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
361   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
362   /* Primal system */
363   ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr);
364   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
365   ierr = SNESSetDM(snes, dm);CHKERRQ(ierr);
366   ierr = SetupDiscretization(dm, SetupPrimalProblem, &user);CHKERRQ(ierr);
367   ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr);
368   ierr = VecSet(u, 0.0);CHKERRQ(ierr);
369   ierr = PetscObjectSetName((PetscObject) u, "potential");CHKERRQ(ierr);
370   ierr = DMPlexSetSNESLocalFEM(dm, &user, &user, &user);CHKERRQ(ierr);
371   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
372   ierr = DMSNESCheckFromOptions(snes, u);CHKERRQ(ierr);
373   ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr);
374   ierr = SNESGetSolution(snes, &u);CHKERRQ(ierr);
375   ierr = VecViewFromOptions(u, NULL, "-potential_view");CHKERRQ(ierr);
376   /* Cleanup */
377   ierr = VecDestroy(&u);CHKERRQ(ierr);
378   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
379   ierr = DMDestroy(&dm);CHKERRQ(ierr);
380   ierr = PetscFinalize();
381   return ierr;
382 }
383 
384 /*TEST
385 
386   test:
387     suffix: 2d_bdm1_p0
388     requires: triangle
389     args: -sol_type linear \
390           -field_petscspace_degree 1 -field_petscdualspace_type bdm -dm_refine 1 \
391           -dmsnes_check .001 -snes_error_if_not_converged \
392           -ksp_rtol 1e-10 -ksp_error_if_not_converged \
393           -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \
394             -fieldsplit_field_pc_type lu \
395             -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu
396   test:
397     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: [2.0, 1.0]
398     # Using -sol_type quadratic -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: [2.9, 1.0]
399     suffix: 2d_bdm1_p0_conv
400     requires: triangle
401     args: -sol_type quartic \
402           -field_petscspace_degree 1 -field_petscdualspace_type bdm -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate \
403           -snes_error_if_not_converged \
404           -ksp_rtol 1e-10 -ksp_error_if_not_converged \
405           -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \
406             -fieldsplit_field_pc_type lu \
407             -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu
408 
409   test:
410     # HDF5 output: -dm_view hdf5:${PETSC_DIR}/sol.h5 -potential_view hdf5:${PETSC_DIR}/sol.h5::append -exact_vec_view hdf5:${PETSC_DIR}/sol.h5::append
411     # VTK output: -potential_view vtk: -exact_vec_view vtk:
412     # VTU output: -potential_view vtk:multifield.vtu -exact_vec_view vtk:exact.vtu
413     suffix: 2d_q2_p0
414     requires: triangle
415     args: -sol_type linear -dm_plex_box_simplex 0 \
416           -field_petscspace_degree 2 -dm_refine 0 \
417           -dmsnes_check .001 -snes_error_if_not_converged \
418           -ksp_rtol 1e-10 -ksp_error_if_not_converged \
419           -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \
420             -fieldsplit_field_pc_type lu \
421             -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu
422   test:
423     # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: [0.0057, 1.0]
424     suffix: 2d_q2_p0_conv
425     requires: triangle
426     args: -sol_type linear -dm_plex_box_simplex 0 \
427           -field_petscspace_degree 2 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate \
428           -snes_error_if_not_converged \
429           -ksp_rtol 1e-10 -ksp_error_if_not_converged \
430           -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \
431             -fieldsplit_field_pc_type lu \
432             -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu
433   test:
434     # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: [-0.014, 0.0066]
435     suffix: 2d_q2_p0_neumann_conv
436     requires: triangle
437     args: -sol_type quartic_neumann -dm_plex_box_simplex 0 \
438           -field_petscspace_degree 2 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate \
439           -snes_error_if_not_converged \
440           -ksp_rtol 1e-10 -ksp_error_if_not_converged \
441           -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \
442             -fieldsplit_field_pc_type lu \
443             -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type svd
444 
445 TEST*/
446 
447 /* These tests will be active once tensor P^- is working
448 
449   test:
450     suffix: 2d_bdmq1_p0_0
451     requires: triangle
452     args: -dm_plex_box_simplex 0 -sol_type linear \
453           -field_petscspace_poly_type pminus_hdiv -field_petscspace_degree 1 -field_petscdualspace_type bdm -dm_refine 0 -convest_num_refine 3 -snes_convergence_estimate \
454           -dmsnes_check .001 -snes_error_if_not_converged \
455           -ksp_rtol 1e-10 -ksp_error_if_not_converged \
456           -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \
457             -fieldsplit_field_pc_type lu \
458             -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu
459   test:
460     suffix: 2d_bdmq1_p0_2
461     requires: triangle
462     args: -dm_plex_box_simplex 0 -sol_type quartic \
463           -field_petscspace_poly_type_no pminus_hdiv -field_petscspace_degree 1 -field_petscdualspace_type bdm -dm_refine 0 -convest_num_refine 3 -snes_convergence_estimate \
464           -dmsnes_check .001 -snes_error_if_not_converged \
465           -ksp_rtol 1e-10 -ksp_error_if_not_converged \
466           -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \
467             -fieldsplit_field_pc_type lu \
468             -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu
469 
470 */
471