xref: /petsc/src/snes/tutorials/ex12.c (revision 4e278199b78715991f5c71ebbd945c1489263e6c)
1 static char help[] = "Poisson Problem in 2d and 3d with simplicial 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 discretized auxiliary fields (conductivity) as well as\n\
5 multilevel nonlinear solvers.\n\n\n";
6 
7 /*
8 A visualization of the adaptation can be accomplished using:
9 
10   -dm_adapt_view hdf5:$PWD/adapt.h5 -sol_adapt_view hdf5:$PWD/adapt.h5::append -dm_adapt_pre_view hdf5:$PWD/orig.h5 -sol_adapt_pre_view hdf5:$PWD/orig.h5::append
11 
12 Information on refinement:
13 
14    -info :~sys,vec,is,mat,ksp,snes,ts
15 */
16 
17 #include <petscdmplex.h>
18 #include <petscdmadaptor.h>
19 #include <petscsnes.h>
20 #include <petscds.h>
21 #include <petscviewerhdf5.h>
22 
23 typedef enum {NEUMANN, DIRICHLET, NONE} BCType;
24 typedef enum {RUN_FULL, RUN_EXACT, RUN_TEST, RUN_PERF} RunType;
25 typedef enum {COEFF_NONE, COEFF_ANALYTIC, COEFF_FIELD, COEFF_NONLINEAR, COEFF_CIRCLE, COEFF_CROSS, COEFF_CHECKERBOARD_0, COEFF_CHECKERBOARD_1} CoeffType;
26 
27 typedef struct {
28   RunType        runType;           /* Whether to run tests, or solve the full problem */
29   PetscBool      jacobianMF;        /* Whether to calculate the Jacobian action on the fly */
30   PetscBool      showInitial, showSolution, restart, quiet, nonzInit;
31   /* Problem definition */
32   BCType         bcType;
33   CoeffType      variableCoefficient;
34   PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx);
35   PetscBool      fieldBC;
36   void           (**exactFields)(PetscInt, PetscInt, PetscInt,
37                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
38                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
39                                  PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);
40   PetscBool      bdIntegral;        /* Compute the integral of the solution on the boundary */
41   /* Reproducing tests from SISC 40(3), pp. A1473-A1493, 2018 */
42   PetscInt       div;               /* Number of divisions */
43   PetscInt       k;                 /* Parameter for checkerboard coefficient */
44   PetscInt      *kgrid;             /* Random parameter grid */
45   PetscBool      rand;              /* Make random assignments */
46   /* Solver */
47   PC             pcmg;              /* This is needed for error monitoring */
48   PetscBool      checkksp;          /* Whether to check the KSPSolve for runType == RUN_TEST */
49 } AppCtx;
50 
51 static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
52 {
53   u[0] = 0.0;
54   return 0;
55 }
56 
57 static PetscErrorCode ecks(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
58 {
59   u[0] = x[0];
60   return 0;
61 }
62 
63 /*
64   In 2D for Dirichlet conditions, we use exact solution:
65 
66     u = x^2 + y^2
67     f = 4
68 
69   so that
70 
71     -\Delta u + f = -4 + 4 = 0
72 
73   For Neumann conditions, we have
74 
75     -\nabla u \cdot -\hat y |_{y=0} =  (2y)|_{y=0} =  0 (bottom)
76     -\nabla u \cdot  \hat y |_{y=1} = -(2y)|_{y=1} = -2 (top)
77     -\nabla u \cdot -\hat x |_{x=0} =  (2x)|_{x=0} =  0 (left)
78     -\nabla u \cdot  \hat x |_{x=1} = -(2x)|_{x=1} = -2 (right)
79 
80   Which we can express as
81 
82     \nabla u \cdot  \hat n|_\Gamma = {2 x, 2 y} \cdot \hat n = 2 (x + y)
83 
84   The boundary integral of this solution is (assuming we are not orienting the edges)
85 
86     \int^1_0 x^2 dx + \int^1_0 (1 + y^2) dy + \int^1_0 (x^2 + 1) dx + \int^1_0 y^2 dy = 1/3 + 4/3 + 4/3 + 1/3 = 3 1/3
87 */
88 static PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
89 {
90   *u = x[0]*x[0] + x[1]*x[1];
91   return 0;
92 }
93 
94 static void quadratic_u_field_2d(PetscInt dim, PetscInt Nf, PetscInt NfAux,
95                                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
96                                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
97                                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar uexact[])
98 {
99   uexact[0] = a[0];
100 }
101 
102 static PetscErrorCode circle_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
103 {
104   const PetscReal alpha   = 500.;
105   const PetscReal radius2 = PetscSqr(0.15);
106   const PetscReal r2      = PetscSqr(x[0] - 0.5) + PetscSqr(x[1] - 0.5);
107   const PetscReal xi      = alpha*(radius2 - r2);
108 
109   *u = PetscTanhScalar(xi) + 1.0;
110   return 0;
111 }
112 
113 static PetscErrorCode cross_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
114 {
115   const PetscReal alpha = 50*4;
116   const PetscReal xy    = (x[0]-0.5)*(x[1]-0.5);
117 
118   *u = PetscSinReal(alpha*xy) * (alpha*PetscAbsReal(xy) < 2*PETSC_PI ? 1 : 0.01);
119   return 0;
120 }
121 
122 static void f0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
123                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
124                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
125                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
126 {
127   f0[0] = 4.0;
128 }
129 
130 static void f0_circle_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
131                         const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
132                         const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
133                         PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
134 {
135   const PetscReal alpha   = 500.;
136   const PetscReal radius2 = PetscSqr(0.15);
137   const PetscReal r2      = PetscSqr(x[0] - 0.5) + PetscSqr(x[1] - 0.5);
138   const PetscReal xi      = alpha*(radius2 - r2);
139 
140   f0[0] = (-4.0*alpha - 8.0*PetscSqr(alpha)*r2*PetscTanhReal(xi)) * PetscSqr(1.0/PetscCoshReal(xi));
141 }
142 
143 static void f0_cross_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
144                        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
145                        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
146                        PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
147 {
148   const PetscReal alpha = 50*4;
149   const PetscReal xy    = (x[0]-0.5)*(x[1]-0.5);
150 
151   f0[0] = PetscSinReal(alpha*xy) * (alpha*PetscAbsReal(xy) < 2*PETSC_PI ? 1 : 0.01);
152 }
153 
154 static void f0_checkerboard_0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
155                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
156                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
157                                 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
158 {
159   f0[0] = -20.0*PetscExpReal(-(PetscSqr(x[0] - 0.5) + PetscSqr(x[1] - 0.5)));
160 }
161 
162 static void f0_bd_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
163                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
164                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
165                     PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
166 {
167   PetscInt d;
168   for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += -n[d]*2.0*x[d];
169 }
170 
171 /* gradU[comp*dim+d] = {u_x, u_y} or {u_x, u_y, u_z} */
172 static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
173                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
174                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
175                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
176 {
177   PetscInt d;
178   for (d = 0; d < dim; ++d) f1[d] = u_x[d];
179 }
180 
181 /* < \nabla v, \nabla u + {\nabla u}^T >
182    This just gives \nabla u, give the perdiagonal for the transpose */
183 static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
184                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
185                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
186                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
187 {
188   PetscInt d;
189   for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0;
190 }
191 
192 /*
193   In 2D for x periodicity and y Dirichlet conditions, we use exact solution:
194 
195     u = sin(2 pi x)
196     f = -4 pi^2 sin(2 pi x)
197 
198   so that
199 
200     -\Delta u + f = 4 pi^2 sin(2 pi x) - 4 pi^2 sin(2 pi x) = 0
201 */
202 static PetscErrorCode xtrig_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
203 {
204   *u = PetscSinReal(2.0*PETSC_PI*x[0]);
205   return 0;
206 }
207 
208 static void f0_xtrig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
209                        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
210                        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
211                        PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
212 {
213   f0[0] = -4.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[0]);
214 }
215 
216 /*
217   In 2D for x-y periodicity, we use exact solution:
218 
219     u = sin(2 pi x) sin(2 pi y)
220     f = -8 pi^2 sin(2 pi x)
221 
222   so that
223 
224     -\Delta u + f = 4 pi^2 sin(2 pi x) sin(2 pi y) + 4 pi^2 sin(2 pi x) sin(2 pi y) - 8 pi^2 sin(2 pi x) = 0
225 */
226 static PetscErrorCode xytrig_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
227 {
228   *u = PetscSinReal(2.0*PETSC_PI*x[0])*PetscSinReal(2.0*PETSC_PI*x[1]);
229   return 0;
230 }
231 
232 static void f0_xytrig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
233                         const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
234                         const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
235                         PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
236 {
237   f0[0] = -8.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[0]);
238 }
239 
240 /*
241   In 2D for Dirichlet conditions with a variable coefficient, we use exact solution:
242 
243     u  = x^2 + y^2
244     f  = 6 (x + y)
245     nu = (x + y)
246 
247   so that
248 
249     -\div \nu \grad u + f = -6 (x + y) + 6 (x + y) = 0
250 */
251 static PetscErrorCode nu_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
252 {
253   *u = x[0] + x[1];
254   return 0;
255 }
256 
257 static PetscErrorCode checkerboardCoeff(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
258 {
259   AppCtx  *user = (AppCtx *) ctx;
260   PetscInt div  = user->div;
261   PetscInt k    = user->k;
262   PetscInt mask = 0, ind = 0, d;
263 
264   PetscFunctionBeginUser;
265   for (d = 0; d < dim; ++d) mask = (mask + (PetscInt) (x[d]*div)) % 2;
266   if (user->kgrid) {
267     for (d = 0; d < dim; ++d) {
268       if (d > 0) ind *= dim;
269       ind += (PetscInt) (x[d]*div);
270     }
271     k = user->kgrid[ind];
272   }
273   u[0] = mask ? 1.0 : PetscPowRealInt(10.0, -k);
274   PetscFunctionReturn(0);
275 }
276 
277 void f0_analytic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
278                    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
279                    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
280                    PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
281 {
282   f0[0] = 6.0*(x[0] + x[1]);
283 }
284 
285 /* gradU[comp*dim+d] = {u_x, u_y} or {u_x, u_y, u_z} */
286 void f1_analytic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
287                    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
288                    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
289                    PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
290 {
291   PetscInt d;
292   for (d = 0; d < dim; ++d) f1[d] = (x[0] + x[1])*u_x[d];
293 }
294 
295 void f1_field_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
296                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
297                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
298                 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
299 {
300   PetscInt d;
301   for (d = 0; d < dim; ++d) f1[d] = a[0]*u_x[d];
302 }
303 
304 /* < \nabla v, \nabla u + {\nabla u}^T >
305    This just gives \nabla u, give the perdiagonal for the transpose */
306 void g3_analytic_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
307                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
308                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
309                     PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
310 {
311   PetscInt d;
312   for (d = 0; d < dim; ++d) g3[d*dim+d] = x[0] + x[1];
313 }
314 
315 void g3_field_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
316                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
317                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
318                  PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
319 {
320   PetscInt d;
321   for (d = 0; d < dim; ++d) g3[d*dim+d] = a[0];
322 }
323 
324 /*
325   In 2D for Dirichlet conditions with a nonlinear coefficient (p-Laplacian with p = 4), we use exact solution:
326 
327     u  = x^2 + y^2
328     f  = 16 (x^2 + y^2)
329     nu = 1/2 |grad u|^2
330 
331   so that
332 
333     -\div \nu \grad u + f = -16 (x^2 + y^2) + 16 (x^2 + y^2) = 0
334 */
335 void f0_analytic_nonlinear_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
336                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
337                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
338                              PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
339 {
340   f0[0] = 16.0*(x[0]*x[0] + x[1]*x[1]);
341 }
342 
343 /* gradU[comp*dim+d] = {u_x, u_y} or {u_x, u_y, u_z} */
344 void f1_analytic_nonlinear_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
345                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
346                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
347                              PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
348 {
349   PetscScalar nu = 0.0;
350   PetscInt    d;
351   for (d = 0; d < dim; ++d) nu += u_x[d]*u_x[d];
352   for (d = 0; d < dim; ++d) f1[d] = 0.5*nu*u_x[d];
353 }
354 
355 /*
356   grad (u + eps w) - grad u = eps grad w
357 
358   1/2 |grad (u + eps w)|^2 grad (u + eps w) - 1/2 |grad u|^2 grad u
359 = 1/2 (|grad u|^2 + 2 eps <grad u,grad w>) (grad u + eps grad w) - 1/2 |grad u|^2 grad u
360 = 1/2 (eps |grad u|^2 grad w + 2 eps <grad u,grad w> grad u)
361 = eps (1/2 |grad u|^2 grad w + grad u <grad u,grad w>)
362 */
363 void g3_analytic_nonlinear_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
364                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
365                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
366                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
367 {
368   PetscScalar nu = 0.0;
369   PetscInt    d, e;
370   for (d = 0; d < dim; ++d) nu += u_x[d]*u_x[d];
371   for (d = 0; d < dim; ++d) {
372     g3[d*dim+d] = 0.5*nu;
373     for (e = 0; e < dim; ++e) {
374       g3[d*dim+e] += u_x[d]*u_x[e];
375     }
376   }
377 }
378 
379 /*
380   In 3D for Dirichlet conditions we use exact solution:
381 
382     u = 2/3 (x^2 + y^2 + z^2)
383     f = 4
384 
385   so that
386 
387     -\Delta u + f = -2/3 * 6 + 4 = 0
388 
389   For Neumann conditions, we have
390 
391     -\nabla u \cdot -\hat z |_{z=0} =  (2z)|_{z=0} =  0 (bottom)
392     -\nabla u \cdot  \hat z |_{z=1} = -(2z)|_{z=1} = -2 (top)
393     -\nabla u \cdot -\hat y |_{y=0} =  (2y)|_{y=0} =  0 (front)
394     -\nabla u \cdot  \hat y |_{y=1} = -(2y)|_{y=1} = -2 (back)
395     -\nabla u \cdot -\hat x |_{x=0} =  (2x)|_{x=0} =  0 (left)
396     -\nabla u \cdot  \hat x |_{x=1} = -(2x)|_{x=1} = -2 (right)
397 
398   Which we can express as
399 
400     \nabla u \cdot  \hat n|_\Gamma = {2 x, 2 y, 2z} \cdot \hat n = 2 (x + y + z)
401 */
402 static PetscErrorCode quadratic_u_3d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
403 {
404   *u = 2.0*(x[0]*x[0] + x[1]*x[1] + x[2]*x[2])/3.0;
405   return 0;
406 }
407 
408 static void quadratic_u_field_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux,
409                                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
410                                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
411                                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar uexact[])
412 {
413   uexact[0] = a[0];
414 }
415 
416 static void bd_integral_2d(PetscInt dim, PetscInt Nf, PetscInt NfAux,
417                            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
418                            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
419                            PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *uint)
420 {
421   uint[0] = u[0];
422 }
423 
424 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
425 {
426   const char    *bcTypes[3]  = {"neumann", "dirichlet", "none"};
427   const char    *runTypes[4] = {"full", "exact", "test", "perf"};
428   const char    *coeffTypes[8] = {"none", "analytic", "field", "nonlinear", "circle", "cross", "checkerboard_0", "checkerboard_1"};
429   PetscInt       bc, run, coeff;
430   PetscErrorCode ierr;
431 
432   PetscFunctionBeginUser;
433   options->runType             = RUN_FULL;
434   options->bcType              = DIRICHLET;
435   options->variableCoefficient = COEFF_NONE;
436   options->fieldBC             = PETSC_FALSE;
437   options->jacobianMF          = PETSC_FALSE;
438   options->showInitial         = PETSC_FALSE;
439   options->showSolution        = PETSC_FALSE;
440   options->restart             = PETSC_FALSE;
441   options->quiet               = PETSC_FALSE;
442   options->nonzInit            = PETSC_FALSE;
443   options->bdIntegral          = PETSC_FALSE;
444   options->checkksp            = PETSC_FALSE;
445   options->div                 = 4;
446   options->k                   = 1;
447   options->kgrid               = NULL;
448   options->rand                = PETSC_FALSE;
449 
450   ierr = PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX");CHKERRQ(ierr);
451   run  = options->runType;
452   ierr = PetscOptionsEList("-run_type", "The run type", "ex12.c", runTypes, 4, runTypes[options->runType], &run, NULL);CHKERRQ(ierr);
453   options->runType = (RunType) run;
454   bc   = options->bcType;
455   ierr = PetscOptionsEList("-bc_type","Type of boundary condition","ex12.c",bcTypes,3,bcTypes[options->bcType],&bc,NULL);CHKERRQ(ierr);
456   options->bcType = (BCType) bc;
457   coeff = options->variableCoefficient;
458   ierr = PetscOptionsEList("-variable_coefficient","Type of variable coefficent","ex12.c",coeffTypes,8,coeffTypes[options->variableCoefficient],&coeff,NULL);CHKERRQ(ierr);
459   options->variableCoefficient = (CoeffType) coeff;
460 
461   ierr = PetscOptionsBool("-field_bc", "Use a field representation for the BC", "ex12.c", options->fieldBC, &options->fieldBC, NULL);CHKERRQ(ierr);
462   ierr = PetscOptionsBool("-jacobian_mf", "Calculate the action of the Jacobian on the fly", "ex12.c", options->jacobianMF, &options->jacobianMF, NULL);CHKERRQ(ierr);
463   ierr = PetscOptionsBool("-show_initial", "Output the initial guess for verification", "ex12.c", options->showInitial, &options->showInitial, NULL);CHKERRQ(ierr);
464   ierr = PetscOptionsBool("-show_solution", "Output the solution for verification", "ex12.c", options->showSolution, &options->showSolution, NULL);CHKERRQ(ierr);
465   ierr = PetscOptionsBool("-restart", "Read in the mesh and solution from a file", "ex12.c", options->restart, &options->restart, NULL);CHKERRQ(ierr);
466   ierr = PetscOptionsBool("-quiet", "Don't print any vecs", "ex12.c", options->quiet, &options->quiet, NULL);CHKERRQ(ierr);
467   ierr = PetscOptionsBool("-nonzero_initial_guess", "nonzero initial guess", "ex12.c", options->nonzInit, &options->nonzInit, NULL);CHKERRQ(ierr);
468   ierr = PetscOptionsBool("-bd_integral", "Compute the integral of the solution on the boundary", "ex12.c", options->bdIntegral, &options->bdIntegral, NULL);CHKERRQ(ierr);
469   if (options->runType == RUN_TEST) {
470     ierr = PetscOptionsBool("-run_test_check_ksp", "Check solution of KSP", "ex12.c", options->checkksp, &options->checkksp, NULL);CHKERRQ(ierr);
471   }
472   ierr = PetscOptionsInt("-div", "The number of division for the checkerboard coefficient", "ex12.c", options->div, &options->div, NULL);CHKERRQ(ierr);
473   ierr = PetscOptionsInt("-k", "The exponent for the checkerboard coefficient", "ex12.c", options->k, &options->k, NULL);CHKERRQ(ierr);
474   ierr = PetscOptionsBool("-k_random", "Assign random k values to checkerboard", "ex12.c", options->rand, &options->rand, NULL);CHKERRQ(ierr);
475   ierr = PetscOptionsEnd();
476   PetscFunctionReturn(0);
477 }
478 
479 static PetscErrorCode CreateBCLabel(DM dm, const char name[])
480 {
481   DM             plex;
482   DMLabel        label;
483   PetscErrorCode ierr;
484 
485   PetscFunctionBeginUser;
486   ierr = DMCreateLabel(dm, name);CHKERRQ(ierr);
487   ierr = DMGetLabel(dm, name, &label);CHKERRQ(ierr);
488   ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr);
489   ierr = DMPlexMarkBoundaryFaces(plex, 1, label);CHKERRQ(ierr);
490   ierr = DMDestroy(&plex);CHKERRQ(ierr);
491   PetscFunctionReturn(0);
492 }
493 
494 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
495 {
496   PetscErrorCode ierr;
497 
498   PetscFunctionBeginUser;
499   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
500   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
501   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
502   {
503     char      convType[256];
504     PetscBool flg;
505 
506     ierr = PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX");CHKERRQ(ierr);
507     ierr = PetscOptionsFList("-dm_plex_convert_type","Convert DMPlex to another format","ex12",DMList,DMPLEX,convType,256,&flg);CHKERRQ(ierr);
508     ierr = PetscOptionsEnd();
509     if (flg) {
510       DM dmConv;
511 
512       ierr = DMConvert(*dm,convType,&dmConv);CHKERRQ(ierr);
513       if (dmConv) {
514         ierr = DMDestroy(dm);CHKERRQ(ierr);
515         *dm  = dmConv;
516       }
517       ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
518       ierr = DMSetUp(*dm);CHKERRQ(ierr);
519     }
520   }
521   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
522   if (user->rand) {
523     PetscRandom r;
524     PetscReal   val;
525     PetscInt    dim, N, i;
526 
527     ierr = DMGetDimension(*dm, &dim);CHKERRQ(ierr);
528     N    = PetscPowInt(user->div, dim);
529     ierr = PetscMalloc1(N, &user->kgrid);CHKERRQ(ierr);
530     ierr = PetscRandomCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr);
531     ierr = PetscRandomSetFromOptions(r);CHKERRQ(ierr);
532     ierr = PetscRandomSetInterval(r, 0.0, user->k);CHKERRQ(ierr);
533     ierr = PetscRandomSetSeed(r, 1973);CHKERRQ(ierr);
534     ierr = PetscRandomSeed(r);CHKERRQ(ierr);
535     for (i = 0; i < N; ++i) {
536       ierr = PetscRandomGetValueReal(r, &val);CHKERRQ(ierr);
537       user->kgrid[i] = 1 + (PetscInt) val;
538     }
539     ierr = PetscRandomDestroy(&r);CHKERRQ(ierr);
540   }
541   PetscFunctionReturn(0);
542 }
543 
544 static PetscErrorCode SetupProblem(DM dm, AppCtx *user)
545 {
546   PetscDS         ds;
547   DMLabel         label;
548   PetscWeakForm   wf;
549   const DMBoundaryType *periodicity;
550   const PetscInt  id = 1;
551   PetscInt        bd, dim;
552   PetscErrorCode  ierr;
553 
554   PetscFunctionBeginUser;
555   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
556   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
557   ierr = DMGetPeriodicity(dm, NULL, NULL, NULL, &periodicity);CHKERRQ(ierr);
558   switch (user->variableCoefficient) {
559   case COEFF_NONE:
560     if (periodicity && periodicity[0]) {
561       if (periodicity && periodicity[1]) {
562         ierr = PetscDSSetResidual(ds, 0, f0_xytrig_u, f1_u);CHKERRQ(ierr);
563         ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu);CHKERRQ(ierr);
564       } else {
565         ierr = PetscDSSetResidual(ds, 0, f0_xtrig_u,  f1_u);CHKERRQ(ierr);
566         ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu);CHKERRQ(ierr);
567       }
568     } else {
569       ierr = PetscDSSetResidual(ds, 0, f0_u, f1_u);CHKERRQ(ierr);
570       ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu);CHKERRQ(ierr);
571     }
572     break;
573   case COEFF_ANALYTIC:
574     ierr = PetscDSSetResidual(ds, 0, f0_analytic_u, f1_analytic_u);CHKERRQ(ierr);
575     ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_analytic_uu);CHKERRQ(ierr);
576     break;
577   case COEFF_FIELD:
578     ierr = PetscDSSetResidual(ds, 0, f0_analytic_u, f1_field_u);CHKERRQ(ierr);
579     ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_field_uu);CHKERRQ(ierr);
580     break;
581   case COEFF_NONLINEAR:
582     ierr = PetscDSSetResidual(ds, 0, f0_analytic_nonlinear_u, f1_analytic_nonlinear_u);CHKERRQ(ierr);
583     ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_analytic_nonlinear_uu);CHKERRQ(ierr);
584     break;
585   case COEFF_CIRCLE:
586     ierr = PetscDSSetResidual(ds, 0, f0_circle_u, f1_u);CHKERRQ(ierr);
587     ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu);CHKERRQ(ierr);
588     break;
589   case COEFF_CROSS:
590     ierr = PetscDSSetResidual(ds, 0, f0_cross_u, f1_u);CHKERRQ(ierr);
591     ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu);CHKERRQ(ierr);
592     break;
593   case COEFF_CHECKERBOARD_0:
594     ierr = PetscDSSetResidual(ds, 0, f0_checkerboard_0_u, f1_field_u);CHKERRQ(ierr);
595     ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_field_uu);CHKERRQ(ierr);
596     break;
597   default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid variable coefficient type %d", user->variableCoefficient);
598   }
599   switch (dim) {
600   case 2:
601     switch (user->variableCoefficient) {
602     case COEFF_CIRCLE:
603       user->exactFuncs[0]  = circle_u_2d;break;
604     case COEFF_CROSS:
605       user->exactFuncs[0]  = cross_u_2d;break;
606     case COEFF_CHECKERBOARD_0:
607       user->exactFuncs[0]  = zero;break;
608     default:
609       if (periodicity && periodicity[0]) {
610         if (periodicity && periodicity[1]) {
611           user->exactFuncs[0] = xytrig_u_2d;
612         } else {
613           user->exactFuncs[0] = xtrig_u_2d;
614         }
615       } else {
616         user->exactFuncs[0]  = quadratic_u_2d;
617         user->exactFields[0] = quadratic_u_field_2d;
618       }
619     }
620     if (user->bcType == NEUMANN) {
621       ierr = DMGetLabel(dm, "boundary", &label);CHKERRQ(ierr);
622       ierr = DMAddBoundary(dm, DM_BC_NATURAL, "wall", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd);CHKERRQ(ierr);
623       ierr = PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
624       ierr = PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_bd_u, 0, NULL);CHKERRQ(ierr);
625     }
626     break;
627   case 3:
628     user->exactFuncs[0]  = quadratic_u_3d;
629     user->exactFields[0] = quadratic_u_field_3d;
630     if (user->bcType == NEUMANN) {
631       ierr = DMGetLabel(dm, "boundary", &label);CHKERRQ(ierr);
632       ierr = DMAddBoundary(dm, DM_BC_NATURAL, "wall", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd);CHKERRQ(ierr);
633       ierr = PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
634       ierr = PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_bd_u, 0, NULL);CHKERRQ(ierr);
635     }
636     break;
637   default:
638     SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %d", dim);
639   }
640   /* Setup constants */
641   switch (user->variableCoefficient) {
642   case COEFF_CHECKERBOARD_0:
643   {
644     PetscScalar constants[2];
645 
646     constants[0] = user->div;
647     constants[1] = user->k;
648     ierr = PetscDSSetConstants(ds, 2, constants);CHKERRQ(ierr);
649   }
650   break;
651   default: break;
652   }
653   ierr = PetscDSSetExactSolution(ds, 0, user->exactFuncs[0], user);CHKERRQ(ierr);
654   /* Setup Boundary Conditions */
655   if (user->bcType == DIRICHLET) {
656     ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr);
657     if (!label) {
658       /* Right now, p4est cannot create labels immediately */
659       ierr = PetscDSAddBoundaryByName(ds, user->fieldBC ? DM_BC_ESSENTIAL_FIELD : DM_BC_ESSENTIAL, "wall", "marker", 1, &id, 0, 0, NULL, user->fieldBC ? (void (*)(void)) user->exactFields[0] : (void (*)(void)) user->exactFuncs[0], NULL, user, NULL);CHKERRQ(ierr);
660     } else {
661       ierr = DMAddBoundary(dm, user->fieldBC ? DM_BC_ESSENTIAL_FIELD : DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, user->fieldBC ? (void (*)(void)) user->exactFields[0] : (void (*)(void)) user->exactFuncs[0], NULL, user, NULL);CHKERRQ(ierr);
662     }
663   }
664   PetscFunctionReturn(0);
665 }
666 
667 static PetscErrorCode SetupMaterial(DM dm, DM dmAux, AppCtx *user)
668 {
669   PetscErrorCode (*matFuncs[1])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx) = {nu_2d};
670   void            *ctx[1];
671   Vec              nu;
672   PetscErrorCode   ierr;
673 
674   PetscFunctionBegin;
675   ctx[0] = user;
676   if (user->variableCoefficient == COEFF_CHECKERBOARD_0) {matFuncs[0] = checkerboardCoeff;}
677   ierr = DMCreateLocalVector(dmAux, &nu);CHKERRQ(ierr);
678   ierr = PetscObjectSetName((PetscObject) nu, "Coefficient");CHKERRQ(ierr);
679   ierr = DMProjectFunctionLocal(dmAux, 0.0, matFuncs, ctx, INSERT_ALL_VALUES, nu);CHKERRQ(ierr);
680   ierr = DMSetAuxiliaryVec(dm, NULL, 0, nu);CHKERRQ(ierr);
681   ierr = VecDestroy(&nu);CHKERRQ(ierr);
682   PetscFunctionReturn(0);
683 }
684 
685 static PetscErrorCode SetupBC(DM dm, DM dmAux, AppCtx *user)
686 {
687   PetscErrorCode (*bcFuncs[1])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx);
688   Vec            uexact;
689   PetscInt       dim;
690   PetscErrorCode ierr;
691 
692   PetscFunctionBegin;
693   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
694   if (dim == 2) bcFuncs[0] = quadratic_u_2d;
695   else          bcFuncs[0] = quadratic_u_3d;
696   ierr = DMCreateLocalVector(dmAux, &uexact);CHKERRQ(ierr);
697   ierr = DMProjectFunctionLocal(dmAux, 0.0, bcFuncs, NULL, INSERT_ALL_VALUES, uexact);CHKERRQ(ierr);
698   ierr = DMSetAuxiliaryVec(dm, NULL, 0, uexact);CHKERRQ(ierr);
699   ierr = VecDestroy(&uexact);CHKERRQ(ierr);
700   PetscFunctionReturn(0);
701 }
702 
703 static PetscErrorCode SetupAuxDM(DM dm, PetscFE feAux, AppCtx *user)
704 {
705   DM             dmAux, coordDM;
706   PetscErrorCode ierr;
707 
708   PetscFunctionBegin;
709   /* MUST call DMGetCoordinateDM() in order to get p4est setup if present */
710   ierr = DMGetCoordinateDM(dm, &coordDM);CHKERRQ(ierr);
711   if (!feAux) PetscFunctionReturn(0);
712   ierr = DMClone(dm, &dmAux);CHKERRQ(ierr);
713   ierr = DMSetCoordinateDM(dmAux, coordDM);CHKERRQ(ierr);
714   ierr = DMSetField(dmAux, 0, NULL, (PetscObject) feAux);CHKERRQ(ierr);
715   ierr = DMCreateDS(dmAux);CHKERRQ(ierr);
716   if (user->fieldBC) {ierr = SetupBC(dm, dmAux, user);CHKERRQ(ierr);}
717   else               {ierr = SetupMaterial(dm, dmAux, user);CHKERRQ(ierr);}
718   ierr = DMDestroy(&dmAux);CHKERRQ(ierr);
719   PetscFunctionReturn(0);
720 }
721 
722 static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
723 {
724   DM             plex, cdm = dm;
725   PetscFE        fe, feAux = NULL;
726   PetscBool      simplex;
727   PetscInt       dim;
728   MPI_Comm       comm;
729   PetscErrorCode ierr;
730 
731   PetscFunctionBeginUser;
732   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
733   ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr);
734   ierr = DMPlexIsSimplex(plex, &simplex);CHKERRQ(ierr);
735   ierr = DMDestroy(&plex);CHKERRQ(ierr);
736   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
737   ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, NULL, -1, &fe);CHKERRQ(ierr);
738   ierr = PetscObjectSetName((PetscObject) fe, "potential");CHKERRQ(ierr);
739   if (user->variableCoefficient == COEFF_FIELD || user->variableCoefficient == COEFF_CHECKERBOARD_0) {
740     ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "mat_", -1, &feAux);CHKERRQ(ierr);
741     ierr = PetscObjectSetName((PetscObject) feAux, "coefficient");CHKERRQ(ierr);
742     ierr = PetscFECopyQuadrature(fe, feAux);CHKERRQ(ierr);
743   } else if (user->fieldBC) {
744     ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "bc_", -1, &feAux);CHKERRQ(ierr);
745     ierr = PetscFECopyQuadrature(fe, feAux);CHKERRQ(ierr);
746   }
747   /* Set discretization and boundary conditions for each mesh */
748   ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr);
749   ierr = DMCreateDS(dm);CHKERRQ(ierr);
750   ierr = SetupProblem(dm, user);CHKERRQ(ierr);
751   while (cdm) {
752     ierr = SetupAuxDM(cdm, feAux, user);CHKERRQ(ierr);
753     if (user->bcType == DIRICHLET) {
754       PetscBool hasLabel;
755 
756       ierr = DMHasLabel(cdm, "marker", &hasLabel);CHKERRQ(ierr);
757       if (!hasLabel) {ierr = CreateBCLabel(cdm, "marker");CHKERRQ(ierr);}
758     }
759     ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr);
760     ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
761   }
762   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
763   ierr = PetscFEDestroy(&feAux);CHKERRQ(ierr);
764   PetscFunctionReturn(0);
765 }
766 
767 int main(int argc, char **argv)
768 {
769   DM             dm;          /* Problem specification */
770   SNES           snes;        /* nonlinear solver */
771   Vec            u;           /* solution vector */
772   Mat            A,J;         /* Jacobian matrix */
773   MatNullSpace   nullSpace;   /* May be necessary for Neumann conditions */
774   AppCtx         user;        /* user-defined work context */
775   JacActionCtx   userJ;       /* context for Jacobian MF action */
776   PetscReal      error = 0.0; /* L_2 error in the solution */
777   PetscErrorCode ierr;
778 
779   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
780   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
781   ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr);
782   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
783   ierr = SNESSetDM(snes, dm);CHKERRQ(ierr);
784   ierr = DMSetApplicationContext(dm, &user);CHKERRQ(ierr);
785 
786   ierr = PetscMalloc2(1, &user.exactFuncs, 1, &user.exactFields);CHKERRQ(ierr);
787   ierr = SetupDiscretization(dm, &user);CHKERRQ(ierr);
788 
789   ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr);
790   ierr = PetscObjectSetName((PetscObject) u, "potential");CHKERRQ(ierr);
791 
792   ierr = DMCreateMatrix(dm, &J);CHKERRQ(ierr);
793   if (user.jacobianMF) {
794     PetscInt M, m, N, n;
795 
796     ierr = MatGetSize(J, &M, &N);CHKERRQ(ierr);
797     ierr = MatGetLocalSize(J, &m, &n);CHKERRQ(ierr);
798     ierr = MatCreate(PETSC_COMM_WORLD, &A);CHKERRQ(ierr);
799     ierr = MatSetSizes(A, m, n, M, N);CHKERRQ(ierr);
800     ierr = MatSetType(A, MATSHELL);CHKERRQ(ierr);
801     ierr = MatSetUp(A);CHKERRQ(ierr);
802 #if 0
803     ierr = MatShellSetOperation(A, MATOP_MULT, (void (*)(void))FormJacobianAction);CHKERRQ(ierr);
804 #endif
805 
806     userJ.dm   = dm;
807     userJ.J    = J;
808     userJ.user = &user;
809 
810     ierr = DMCreateLocalVector(dm, &userJ.u);CHKERRQ(ierr);
811     if (user.fieldBC) {ierr = DMProjectFieldLocal(dm, 0.0, userJ.u, user.exactFields, INSERT_BC_VALUES, userJ.u);CHKERRQ(ierr);}
812     else              {ierr = DMProjectFunctionLocal(dm, 0.0, user.exactFuncs, NULL, INSERT_BC_VALUES, userJ.u);CHKERRQ(ierr);}
813     ierr = MatShellSetContext(A, &userJ);CHKERRQ(ierr);
814   } else {
815     A = J;
816   }
817 
818   nullSpace = NULL;
819   if (user.bcType != DIRICHLET) {
820     ierr = MatNullSpaceCreate(PetscObjectComm((PetscObject) dm), PETSC_TRUE, 0, NULL, &nullSpace);CHKERRQ(ierr);
821     ierr = MatSetNullSpace(A, nullSpace);CHKERRQ(ierr);
822   }
823 
824   ierr = DMPlexSetSNESLocalFEM(dm,&user,&user,&user);CHKERRQ(ierr);
825   ierr = SNESSetJacobian(snes, A, J, NULL, NULL);CHKERRQ(ierr);
826 
827   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
828 
829   if (user.fieldBC) {ierr = DMProjectField(dm, 0.0, u, user.exactFields, INSERT_ALL_VALUES, u);CHKERRQ(ierr);}
830   else              {ierr = DMProjectFunction(dm, 0.0, user.exactFuncs, NULL, INSERT_ALL_VALUES, u);CHKERRQ(ierr);}
831   if (user.restart) {
832 #if defined(PETSC_HAVE_HDF5)
833     PetscViewer viewer;
834     char        filename[PETSC_MAX_PATH_LEN];
835 
836     ierr = PetscOptionsGetString(NULL, NULL, "-dm_plex_filename", filename, sizeof(filename), NULL);CHKERRQ(ierr);
837     ierr = PetscViewerCreate(PETSC_COMM_WORLD, &viewer);CHKERRQ(ierr);
838     ierr = PetscViewerSetType(viewer, PETSCVIEWERHDF5);CHKERRQ(ierr);
839     ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr);
840     ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr);
841     ierr = PetscViewerHDF5PushGroup(viewer, "/fields");CHKERRQ(ierr);
842     ierr = VecLoad(u, viewer);CHKERRQ(ierr);
843     ierr = PetscViewerHDF5PopGroup(viewer);CHKERRQ(ierr);
844     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
845 #endif
846   }
847   if (user.showInitial) {
848     Vec lv;
849     ierr = DMGetLocalVector(dm, &lv);CHKERRQ(ierr);
850     ierr = DMGlobalToLocalBegin(dm, u, INSERT_VALUES, lv);CHKERRQ(ierr);
851     ierr = DMGlobalToLocalEnd(dm, u, INSERT_VALUES, lv);CHKERRQ(ierr);
852     ierr = DMPrintLocalVec(dm, "Local function", 1.0e-10, lv);CHKERRQ(ierr);
853     ierr = DMRestoreLocalVector(dm, &lv);CHKERRQ(ierr);
854   }
855   if (user.runType == RUN_FULL || user.runType == RUN_EXACT) {
856     PetscErrorCode (*initialGuess[1])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx) = {zero};
857 
858     if (user.nonzInit) initialGuess[0] = ecks;
859     if (user.runType == RUN_FULL) {
860       ierr = DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u);CHKERRQ(ierr);
861     }
862     ierr = VecViewFromOptions(u, NULL, "-guess_vec_view");CHKERRQ(ierr);
863     ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr);
864     ierr = SNESGetSolution(snes, &u);CHKERRQ(ierr);
865     ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr);
866 
867     if (user.showSolution) {
868       ierr = PetscPrintf(PETSC_COMM_WORLD, "Solution\n");CHKERRQ(ierr);
869       ierr = VecChop(u, 3.0e-9);CHKERRQ(ierr);
870       ierr = VecView(u, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
871     }
872   } else if (user.runType == RUN_PERF) {
873     Vec       r;
874     PetscReal res = 0.0;
875 
876     ierr = SNESGetFunction(snes, &r, NULL, NULL);CHKERRQ(ierr);
877     ierr = SNESComputeFunction(snes, u, r);CHKERRQ(ierr);
878     ierr = PetscPrintf(PETSC_COMM_WORLD, "Initial Residual\n");CHKERRQ(ierr);
879     ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr);
880     ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr);
881     ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Residual: %g\n", (double)res);CHKERRQ(ierr);
882   } else {
883     Vec       r;
884     PetscReal res = 0.0, tol = 1.0e-11;
885 
886     /* Check discretization error */
887     ierr = SNESGetFunction(snes, &r, NULL, NULL);CHKERRQ(ierr);
888     ierr = PetscPrintf(PETSC_COMM_WORLD, "Initial guess\n");CHKERRQ(ierr);
889     if (!user.quiet) {ierr = VecView(u, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);}
890     ierr = DMComputeL2Diff(dm, 0.0, user.exactFuncs, NULL, u, &error);CHKERRQ(ierr);
891     if (error < tol) {ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: < %2.1e\n", (double)tol);CHKERRQ(ierr);}
892     else             {ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: %g\n", (double)error);CHKERRQ(ierr);}
893     /* Check residual */
894     ierr = SNESComputeFunction(snes, u, r);CHKERRQ(ierr);
895     ierr = PetscPrintf(PETSC_COMM_WORLD, "Initial Residual\n");CHKERRQ(ierr);
896     ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr);
897     if (!user.quiet) {ierr = VecView(r, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);}
898     ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr);
899     ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Residual: %g\n", (double)res);CHKERRQ(ierr);
900     /* Check Jacobian */
901     {
902       Vec b;
903 
904       ierr = SNESComputeJacobian(snes, u, A, A);CHKERRQ(ierr);
905       ierr = VecDuplicate(u, &b);CHKERRQ(ierr);
906       ierr = VecSet(r, 0.0);CHKERRQ(ierr);
907       ierr = SNESComputeFunction(snes, r, b);CHKERRQ(ierr);
908       ierr = MatMult(A, u, r);CHKERRQ(ierr);
909       ierr = VecAXPY(r, 1.0, b);CHKERRQ(ierr);
910       ierr = PetscPrintf(PETSC_COMM_WORLD, "Au - b = Au + F(0)\n");CHKERRQ(ierr);
911       ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr);
912       if (!user.quiet) {ierr = VecView(r, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);}
913       ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr);
914       ierr = PetscPrintf(PETSC_COMM_WORLD, "Linear L_2 Residual: %g\n", (double)res);CHKERRQ(ierr);
915       /* check solver */
916       if (user.checkksp) {
917         KSP ksp;
918 
919         if (nullSpace) {
920           ierr = MatNullSpaceRemove(nullSpace, u);CHKERRQ(ierr);
921         }
922         ierr = SNESComputeJacobian(snes, u, A, J);CHKERRQ(ierr);
923         ierr = MatMult(A, u, b);CHKERRQ(ierr);
924         ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr);
925         ierr = KSPSetOperators(ksp, A, J);CHKERRQ(ierr);
926         ierr = KSPSolve(ksp, b, r);CHKERRQ(ierr);
927         ierr = VecAXPY(r, -1.0, u);CHKERRQ(ierr);
928         ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr);
929         ierr = PetscPrintf(PETSC_COMM_WORLD, "KSP Error: %g\n", (double)res);CHKERRQ(ierr);
930       }
931       ierr = VecDestroy(&b);CHKERRQ(ierr);
932     }
933   }
934   ierr = VecViewFromOptions(u, NULL, "-vec_view");CHKERRQ(ierr);
935   {
936     Vec nu;
937 
938     ierr = DMGetAuxiliaryVec(dm, NULL, 0, &nu);CHKERRQ(ierr);
939     if (nu) {ierr = VecViewFromOptions(nu, NULL, "-coeff_view");CHKERRQ(ierr);}
940   }
941 
942   if (user.bdIntegral) {
943     DMLabel   label;
944     PetscInt  id = 1;
945     PetscScalar bdInt = 0.0;
946     PetscReal   exact = 3.3333333333;
947 
948     ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr);
949     ierr = DMPlexComputeBdIntegral(dm, u, label, 1, &id, bd_integral_2d, &bdInt, NULL);CHKERRQ(ierr);
950     ierr = PetscPrintf(PETSC_COMM_WORLD, "Solution boundary integral: %.4g\n", (double) PetscAbsScalar(bdInt));CHKERRQ(ierr);
951     if (PetscAbsReal(PetscAbsScalar(bdInt) - exact) > PETSC_SQRT_MACHINE_EPSILON) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Invalid boundary integral %g != %g", (double) PetscAbsScalar(bdInt), (double)exact);
952   }
953 
954   ierr = MatNullSpaceDestroy(&nullSpace);CHKERRQ(ierr);
955   if (user.jacobianMF) {ierr = VecDestroy(&userJ.u);CHKERRQ(ierr);}
956   if (A != J) {ierr = MatDestroy(&A);CHKERRQ(ierr);}
957   ierr = MatDestroy(&J);CHKERRQ(ierr);
958   ierr = VecDestroy(&u);CHKERRQ(ierr);
959   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
960   ierr = DMDestroy(&dm);CHKERRQ(ierr);
961   ierr = PetscFree2(user.exactFuncs, user.exactFields);CHKERRQ(ierr);
962   ierr = PetscFree(user.kgrid);CHKERRQ(ierr);
963   ierr = PetscFinalize();
964   return ierr;
965 }
966 
967 /*TEST
968   # 2D serial P1 test 0-4
969   test:
970     suffix: 2d_p1_0
971     requires: triangle
972     args: -run_type test -bc_type dirichlet -dm_plex_interpolate 0 -petscspace_degree 1 -show_initial -dm_plex_print_fem 1
973 
974   test:
975     suffix: 2d_p1_1
976     requires: triangle
977     args: -run_type test -bc_type dirichlet -petscspace_degree 1 -show_initial -dm_plex_print_fem 1
978 
979   test:
980     suffix: 2d_p1_2
981     requires: triangle
982     args: -run_type test -dm_refine_volume_limit_pre 0.0625 -bc_type dirichlet -petscspace_degree 1 -show_initial -dm_plex_print_fem 1
983 
984   test:
985     suffix: 2d_p1_neumann_0
986     requires: triangle
987     args: -dm_coord_space 0 -run_type test -bc_type neumann -dm_plex_boundary_label boundary -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 -dm_view ascii::ascii_info_detail
988 
989   test:
990     suffix: 2d_p1_neumann_1
991     requires: triangle
992     args: -run_type test -dm_refine_volume_limit_pre 0.0625 -bc_type neumann -dm_plex_boundary_label boundary -petscspace_degree 1 -show_initial -dm_plex_print_fem 1
993 
994   # 2D serial P2 test 5-8
995   test:
996     suffix: 2d_p2_0
997     requires: triangle
998     args: -run_type test -bc_type dirichlet -petscspace_degree 2 -show_initial -dm_plex_print_fem 1
999 
1000   test:
1001     suffix: 2d_p2_1
1002     requires: triangle
1003     args: -run_type test -dm_refine_volume_limit_pre 0.0625 -bc_type dirichlet -petscspace_degree 2 -show_initial -dm_plex_print_fem 1
1004 
1005   test:
1006     suffix: 2d_p2_neumann_0
1007     requires: triangle
1008     args: -dm_coord_space 0 -run_type test -bc_type neumann -dm_plex_boundary_label boundary -petscspace_degree 2 -show_initial -dm_plex_print_fem 1 -dm_view ascii::ascii_info_detail
1009 
1010   test:
1011     suffix: 2d_p2_neumann_1
1012     requires: triangle
1013     args: -dm_coord_space 0 -run_type test -dm_refine_volume_limit_pre 0.0625 -bc_type neumann -dm_plex_boundary_label boundary -petscspace_degree 2 -show_initial -dm_plex_print_fem 1 -dm_view ascii::ascii_info_detail
1014 
1015   test:
1016     suffix: bd_int_0
1017     requires: triangle
1018     args: -run_type test -bc_type dirichlet -petscspace_degree 2 -bd_integral -dm_view -quiet
1019 
1020   test:
1021     suffix: bd_int_1
1022     requires: triangle
1023     args: -run_type test -dm_refine 2 -bc_type dirichlet -petscspace_degree 2 -bd_integral -dm_view -quiet
1024 
1025   # 3D serial P1 test 9-12
1026   test:
1027     suffix: 3d_p1_0
1028     requires: ctetgen
1029     args: -run_type test -dm_plex_dim 3 -bc_type dirichlet -dm_plex_interpolate 0 -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 -dm_view
1030 
1031   test:
1032     suffix: 3d_p1_1
1033     requires: ctetgen
1034     args: -run_type test -dm_plex_dim 3 -bc_type dirichlet -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 -dm_view
1035 
1036   test:
1037     suffix: 3d_p1_2
1038     requires: ctetgen
1039     args: -run_type test -dm_plex_dim 3 -dm_refine_volume_limit_pre 0.0125 -bc_type dirichlet -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 -dm_view
1040 
1041   test:
1042     suffix: 3d_p1_neumann_0
1043     requires: ctetgen
1044     args: -run_type test -dm_plex_dim 3 -bc_type neumann -dm_plex_boundary_label boundary -petscspace_degree 1 -snes_fd -show_initial -dm_plex_print_fem 1 -dm_view
1045 
1046   # Analytic variable coefficient 13-20
1047   test:
1048     suffix: 13
1049     requires: triangle
1050     args: -run_type test -variable_coefficient analytic -petscspace_degree 1 -show_initial -dm_plex_print_fem 1
1051   test:
1052     suffix: 14
1053     requires: triangle
1054     args: -run_type test -dm_refine_volume_limit_pre 0.0625 -variable_coefficient analytic -petscspace_degree 1 -show_initial -dm_plex_print_fem 1
1055   test:
1056     suffix: 15
1057     requires: triangle
1058     args: -run_type test -variable_coefficient analytic -petscspace_degree 2 -show_initial -dm_plex_print_fem 1
1059   test:
1060     suffix: 16
1061     requires: triangle
1062     args: -run_type test -dm_refine_volume_limit_pre 0.0625 -variable_coefficient analytic -petscspace_degree 2 -show_initial -dm_plex_print_fem 1
1063   test:
1064     suffix: 17
1065     requires: ctetgen
1066     args: -run_type test -dm_plex_dim 3 -variable_coefficient analytic -petscspace_degree 1 -show_initial -dm_plex_print_fem 1
1067 
1068   test:
1069     suffix: 18
1070     requires: ctetgen
1071     args: -run_type test -dm_plex_dim 3 -dm_refine_volume_limit_pre 0.0125 -variable_coefficient analytic -petscspace_degree 1 -show_initial -dm_plex_print_fem 1
1072 
1073   test:
1074     suffix: 19
1075     requires: ctetgen
1076     args: -run_type test -dm_plex_dim 3 -variable_coefficient analytic -petscspace_degree 2 -show_initial -dm_plex_print_fem 1
1077 
1078   test:
1079     suffix: 20
1080     requires: ctetgen
1081     args: -run_type test -dm_plex_dim 3 -dm_refine_volume_limit_pre 0.0125 -variable_coefficient analytic -petscspace_degree 2 -show_initial -dm_plex_print_fem 1
1082 
1083   # P1 variable coefficient 21-28
1084   test:
1085     suffix: 21
1086     requires: triangle
1087     args: -run_type test -variable_coefficient field -petscspace_degree 1 -mat_petscspace_degree 1 -show_initial -dm_plex_print_fem 1
1088 
1089   test:
1090     suffix: 22
1091     requires: triangle
1092     args: -run_type test -dm_refine_volume_limit_pre 0.0625 -variable_coefficient field -petscspace_degree 1 -mat_petscspace_degree 1 -show_initial -dm_plex_print_fem 1
1093 
1094   test:
1095     suffix: 23
1096     requires: triangle
1097     args: -run_type test -variable_coefficient field -petscspace_degree 2 -mat_petscspace_degree 1 -show_initial -dm_plex_print_fem 1
1098 
1099   test:
1100     suffix: 24
1101     requires: triangle
1102     args: -run_type test -dm_refine_volume_limit_pre 0.0625 -variable_coefficient field -petscspace_degree 2 -mat_petscspace_degree 1 -show_initial -dm_plex_print_fem 1
1103 
1104   test:
1105     suffix: 25
1106     requires: ctetgen
1107     args: -run_type test -dm_plex_dim 3 -variable_coefficient field -petscspace_degree 1 -mat_petscspace_degree 1 -show_initial -dm_plex_print_fem 1
1108 
1109   test:
1110     suffix: 26
1111     requires: ctetgen
1112     args: -run_type test -dm_plex_dim 3 -dm_refine_volume_limit_pre 0.0125 -variable_coefficient field -petscspace_degree 1 -mat_petscspace_degree 1 -show_initial -dm_plex_print_fem 1
1113 
1114   test:
1115     suffix: 27
1116     requires: ctetgen
1117     args: -run_type test -dm_plex_dim 3 -variable_coefficient field -petscspace_degree 2 -mat_petscspace_degree 1 -show_initial -dm_plex_print_fem 1
1118 
1119   test:
1120     suffix: 28
1121     requires: ctetgen
1122     args: -run_type test -dm_plex_dim 3 -dm_refine_volume_limit_pre 0.0125 -variable_coefficient field -petscspace_degree 2 -mat_petscspace_degree 1 -show_initial -dm_plex_print_fem 1
1123 
1124   # P0 variable coefficient 29-36
1125   test:
1126     suffix: 29
1127     requires: triangle
1128     args: -run_type test -variable_coefficient field -petscspace_degree 1 -show_initial -dm_plex_print_fem 1
1129 
1130   test:
1131     suffix: 30
1132     requires: triangle
1133     args: -run_type test -dm_refine_volume_limit_pre 0.0625 -variable_coefficient field -petscspace_degree 1 -show_initial -dm_plex_print_fem 1
1134 
1135   test:
1136     suffix: 31
1137     requires: triangle
1138     args: -run_type test -variable_coefficient field -petscspace_degree 2 -show_initial -dm_plex_print_fem 1
1139 
1140   test:
1141     requires: triangle
1142     suffix: 32
1143     args: -run_type test -dm_refine_volume_limit_pre 0.0625 -variable_coefficient field -petscspace_degree 2 -show_initial -dm_plex_print_fem 1
1144 
1145   test:
1146     requires: ctetgen
1147     suffix: 33
1148     args: -run_type test -dm_plex_dim 3 -variable_coefficient field -petscspace_degree 1 -show_initial -dm_plex_print_fem 1
1149 
1150   test:
1151     suffix: 34
1152     requires: ctetgen
1153     args: -run_type test -dm_plex_dim 3 -dm_refine_volume_limit_pre 0.0125 -variable_coefficient field -petscspace_degree 1 -show_initial -dm_plex_print_fem 1
1154 
1155   test:
1156     suffix: 35
1157     requires: ctetgen
1158     args: -run_type test -dm_plex_dim 3 -variable_coefficient field -petscspace_degree 2 -show_initial -dm_plex_print_fem 1
1159 
1160   test:
1161     suffix: 36
1162     requires: ctetgen
1163     args: -run_type test -dm_plex_dim 3 -dm_refine_volume_limit_pre 0.0125 -variable_coefficient field -petscspace_degree 2 -show_initial -dm_plex_print_fem 1
1164 
1165   # Full solve 39-44
1166   test:
1167     suffix: 39
1168     requires: triangle !single
1169     args: -run_type full -dm_refine_volume_limit_pre 0.015625 -petscspace_degree 2 -pc_type gamg -ksp_rtol 1.0e-10 -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason ::ascii_info_detail
1170   test:
1171     suffix: 40
1172     requires: triangle !single
1173     args: -run_type full -dm_refine_volume_limit_pre 0.015625 -variable_coefficient nonlinear -petscspace_degree 2 -pc_type svd -ksp_rtol 1.0e-10 -snes_monitor_short -snes_converged_reason ::ascii_info_detail
1174   test:
1175     suffix: 41
1176     requires: triangle !single
1177     args: -run_type full -dm_refine_volume_limit_pre 0.03125 -variable_coefficient nonlinear -petscspace_degree 1 -snes_type fas -snes_fas_levels 2 -fas_coarse_pc_type svd -fas_coarse_ksp_rtol 1.0e-10 -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -dm_refine_hierarchy 1 -snes_view -fas_levels_1_snes_type newtonls -fas_levels_1_pc_type svd -fas_levels_1_ksp_rtol 1.0e-10 -fas_levels_1_snes_monitor_short
1178   test:
1179     suffix: 42
1180     requires: triangle !single
1181     args: -run_type full -dm_refine_volume_limit_pre 0.0625 -variable_coefficient nonlinear -petscspace_degree 1 -snes_type fas -snes_fas_levels 3 -fas_coarse_pc_type svd -fas_coarse_ksp_rtol 1.0e-10 -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -dm_refine_hierarchy 2 -snes_view -fas_levels_1_snes_type newtonls -fas_levels_1_pc_type svd -fas_levels_1_ksp_rtol 1.0e-10 -fas_levels_1_snes_monitor_short -fas_levels_2_snes_type newtonls -fas_levels_2_pc_type svd -fas_levels_2_ksp_rtol 1.0e-10 -fas_levels_2_snes_atol 1.0e-11 -fas_levels_2_snes_monitor_short
1182   test:
1183     suffix: 43
1184     requires: triangle !single
1185     nsize: 2
1186     args: -run_type full -dm_distribute -dm_refine_volume_limit_pre 0.03125 -variable_coefficient nonlinear -petscspace_degree 1 -snes_type fas -snes_fas_levels 2 -fas_coarse_pc_type svd -fas_coarse_ksp_rtol 1.0e-10 -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -dm_refine_hierarchy 1 -snes_view -fas_levels_1_snes_type newtonls -fas_levels_1_pc_type svd -fas_levels_1_ksp_rtol 1.0e-10 -fas_levels_1_snes_monitor_short
1187 
1188   test:
1189     suffix: 44
1190     requires: triangle !single
1191     nsize: 2
1192     args: -run_type full -dm_distribute -dm_refine_volume_limit_pre 0.0625 -variable_coefficient nonlinear -petscspace_degree 1 -snes_type fas -snes_fas_levels 3 -fas_coarse_pc_type svd -fas_coarse_ksp_rtol 1.0e-10 -fas_coarse_snes_monitor_short -snes_monitor_short  -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -dm_refine_hierarchy 2 -dm_plex_print_fem 0 -snes_view -fas_levels_1_snes_type newtonls -fas_levels_1_pc_type svd -fas_levels_1_ksp_rtol 1.0e-10 -fas_levels_1_snes_monitor_short -fas_levels_2_snes_type newtonls -fas_levels_2_pc_type svd -fas_levels_2_ksp_rtol 1.0e-10 -fas_levels_2_snes_atol 1.0e-11 -fas_levels_2_snes_monitor_short
1193 
1194   # These tests use a loose tolerance just to exercise the PtAP operations for MATIS and multiple PCBDDC setup calls inside PCMG
1195   testset:
1196     requires: triangle !single
1197     nsize: 3
1198     args: -run_type full -dm_distribute -petscspace_degree 1 -dm_mat_type is -pc_type mg -pc_mg_levels 2 -mg_coarse_pc_type bddc -pc_mg_galerkin pmat -ksp_rtol 1.0e-2 -snes_converged_reason -dm_refine_hierarchy 2 -snes_max_it 4
1199     test:
1200       suffix: gmg_bddc
1201       filter: sed -e "s/CONVERGED_FNORM_RELATIVE iterations 3/CONVERGED_FNORM_RELATIVE iterations 4/g"
1202       args: -mg_levels_pc_type jacobi
1203     test:
1204       filter: sed -e "s/iterations [0-4]/iterations 4/g"
1205       suffix: gmg_bddc_lev
1206       args: -mg_levels_pc_type bddc
1207 
1208   # Restarting
1209   testset:
1210     suffix: restart
1211     requires: hdf5 triangle !complex
1212     args: -run_type test -bc_type dirichlet -petscspace_degree 1
1213     test:
1214       args: -dm_view hdf5:sol.h5 -vec_view hdf5:sol.h5::append
1215     test:
1216       args: -dm_plex_filename sol.h5 -restart
1217 
1218   # Periodicity
1219   test:
1220     suffix: periodic_0
1221     requires: triangle
1222     args: -run_type full -bc_type dirichlet -petscspace_degree 1 -snes_converged_reason ::ascii_info_detail
1223 
1224   test:
1225     requires: !complex
1226     suffix: periodic_1
1227     args: -quiet -run_type test -dm_plex_simplex 0 -dm_plex_box_faces 3,3 -dm_plex_box_bd periodic,periodic -vec_view vtk:test.vtu:vtk_vtu -petscspace_degree 1 -dm_refine 1
1228 
1229   # 2D serial P1 test with field bc
1230   test:
1231     suffix: field_bc_2d_p1_0
1232     requires: triangle
1233     args: -run_type test -bc_type dirichlet -field_bc -petscspace_degree 1 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1
1234 
1235   test:
1236     suffix: field_bc_2d_p1_1
1237     requires: triangle
1238     args: -run_type test -dm_refine 1 -bc_type dirichlet -field_bc -petscspace_degree 1 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1
1239 
1240   test:
1241     suffix: field_bc_2d_p1_neumann_0
1242     requires: triangle
1243     args: -run_type test -bc_type neumann -dm_plex_boundary_label boundary -field_bc -petscspace_degree 1 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1
1244 
1245   test:
1246     suffix: field_bc_2d_p1_neumann_1
1247     requires: triangle
1248     args: -run_type test -dm_refine 1 -bc_type neumann -dm_plex_boundary_label boundary -field_bc -petscspace_degree 1 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1
1249 
1250   # 3D serial P1 test with field bc
1251   test:
1252     suffix: field_bc_3d_p1_0
1253     requires: ctetgen
1254     args: -run_type test -dm_plex_dim 3 -bc_type dirichlet -field_bc -petscspace_degree 1 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1
1255 
1256   test:
1257     suffix: field_bc_3d_p1_1
1258     requires: ctetgen
1259     args: -run_type test -dm_plex_dim 3 -dm_refine 1 -bc_type dirichlet -field_bc -petscspace_degree 1 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1
1260 
1261   test:
1262     suffix: field_bc_3d_p1_neumann_0
1263     requires: ctetgen
1264     args: -run_type test -dm_plex_dim 3 -bc_type neumann -dm_plex_boundary_label boundary -field_bc -petscspace_degree 1 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1
1265 
1266   test:
1267     suffix: field_bc_3d_p1_neumann_1
1268     requires: ctetgen
1269     args: -run_type test -dm_plex_dim 3 -dm_refine 1 -bc_type neumann -dm_plex_boundary_label boundary -field_bc -petscspace_degree 1 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1
1270 
1271   # 2D serial P2 test with field bc
1272   test:
1273     suffix: field_bc_2d_p2_0
1274     requires: triangle
1275     args: -run_type test -bc_type dirichlet -field_bc -petscspace_degree 2 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1
1276 
1277   test:
1278     suffix: field_bc_2d_p2_1
1279     requires: triangle
1280     args: -run_type test -dm_refine 1 -bc_type dirichlet -field_bc -petscspace_degree 2 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1
1281 
1282   test:
1283     suffix: field_bc_2d_p2_neumann_0
1284     requires: triangle
1285     args: -run_type test -bc_type neumann -dm_plex_boundary_label boundary -field_bc -petscspace_degree 2 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1
1286 
1287   test:
1288     suffix: field_bc_2d_p2_neumann_1
1289     requires: triangle
1290     args: -run_type test -dm_refine 1 -bc_type neumann -dm_plex_boundary_label boundary -field_bc -petscspace_degree 2 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1
1291 
1292   # 3D serial P2 test with field bc
1293   test:
1294     suffix: field_bc_3d_p2_0
1295     requires: ctetgen
1296     args: -run_type test -dm_plex_dim 3 -bc_type dirichlet -field_bc -petscspace_degree 2 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1
1297 
1298   test:
1299     suffix: field_bc_3d_p2_1
1300     requires: ctetgen
1301     args: -run_type test -dm_plex_dim 3 -dm_refine 1 -bc_type dirichlet -field_bc -petscspace_degree 2 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1
1302 
1303   test:
1304     suffix: field_bc_3d_p2_neumann_0
1305     requires: ctetgen
1306     args: -run_type test -dm_plex_dim 3 -bc_type neumann -dm_plex_boundary_label boundary -field_bc -petscspace_degree 2 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1
1307 
1308   test:
1309     suffix: field_bc_3d_p2_neumann_1
1310     requires: ctetgen
1311     args: -run_type test -dm_plex_dim 3 -dm_refine 1 -bc_type neumann -dm_plex_boundary_label boundary -field_bc -petscspace_degree 2 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1
1312 
1313   # Full solve simplex: Convergence
1314   test:
1315     suffix: 3d_p1_conv
1316     requires: ctetgen
1317     args: -run_type full -dm_plex_dim 3 -dm_refine 1 -bc_type dirichlet -petscspace_degree 1 \
1318       -snes_convergence_estimate -convest_num_refine 1 -pc_type lu
1319 
1320   # Full solve simplex: PCBDDC
1321   test:
1322     suffix: tri_bddc
1323     requires: triangle !single
1324     nsize: 5
1325     args: -run_type full -dm_distribute -petscpartitioner_type simple -dm_refine 2 -bc_type dirichlet -petscspace_degree 1 -ksp_type gmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -dm_mat_type is -pc_type bddc -snes_monitor_short -ksp_monitor_short -snes_converged_reason ::ascii_info_detail -ksp_converged_reason -snes_view -show_solution 0
1326 
1327   # Full solve simplex: PCBDDC
1328   test:
1329     suffix: tri_parmetis_bddc
1330     requires: triangle !single parmetis
1331     nsize: 4
1332     args: -run_type full -dm_distribute -petscpartitioner_type parmetis -dm_refine 2 -bc_type dirichlet -petscspace_degree 1 -ksp_type gmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -dm_mat_type is -pc_type bddc -snes_monitor_short -ksp_monitor_short -snes_converged_reason ::ascii_info_detail -ksp_converged_reason -snes_view -show_solution 0
1333 
1334   testset:
1335     args: -run_type full -dm_distribute -dm_plex_simplex 0 -dm_plex_box_faces 3,3 -petscpartitioner_type simple -dm_refine 2 -bc_type dirichlet -petscspace_degree 2 -dm_mat_type is -pc_type bddc -ksp_type gmres -snes_monitor_short -ksp_monitor_short -snes_view -petscspace_poly_tensor -pc_bddc_corner_selection -ksp_rtol 1.e-9 -pc_bddc_use_edges 0
1336     nsize: 5
1337     output_file: output/ex12_quad_bddc.out
1338     filter: sed -e "s/aijcusparse/aij/g" -e "s/aijviennacl/aij/g" -e "s/factorization: cusparse/factorization: petsc/g"
1339     test:
1340       requires: !single
1341       suffix: quad_bddc
1342     test:
1343       requires: !single cuda
1344       suffix: quad_bddc_cuda
1345       args: -matis_localmat_type aijcusparse -pc_bddc_dirichlet_pc_factor_mat_solver_type cusparse -pc_bddc_neumann_pc_factor_mat_solver_type cusparse
1346     test:
1347       requires: !single viennacl
1348       suffix: quad_bddc_viennacl
1349       args: -matis_localmat_type aijviennacl
1350 
1351   # Full solve simplex: ASM
1352   test:
1353     suffix: tri_q2q1_asm_lu
1354     requires: triangle !single
1355     args: -run_type full -dm_refine 3 -bc_type dirichlet -petscspace_degree 1 -ksp_type gmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -pc_type asm -pc_asm_type restrict -pc_asm_blocks 4 -sub_pc_type lu -snes_monitor_short -ksp_monitor_short -snes_converged_reason ::ascii_info_detail -ksp_converged_reason -snes_view -show_solution 0
1356 
1357   test:
1358     suffix: tri_q2q1_msm_lu
1359     requires: triangle !single
1360     args: -run_type full -dm_refine 3 -bc_type dirichlet -petscspace_degree 1 -ksp_type gmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -pc_type asm -pc_asm_type restrict -pc_asm_local_type multiplicative -pc_asm_blocks 4 -sub_pc_type lu -snes_monitor_short -ksp_monitor_short -snes_converged_reason ::ascii_info_detail -ksp_converged_reason -snes_view -show_solution 0
1361 
1362   test:
1363     suffix: tri_q2q1_asm_sor
1364     requires: triangle !single
1365     args: -run_type full -dm_refine 3 -bc_type dirichlet -petscspace_degree 1 -ksp_type gmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -pc_type asm -pc_asm_type restrict -pc_asm_blocks 4 -sub_pc_type sor -snes_monitor_short -ksp_monitor_short -snes_converged_reason ::ascii_info_detail -ksp_converged_reason -snes_view -show_solution 0
1366 
1367   test:
1368     suffix: tri_q2q1_msm_sor
1369     requires: triangle !single
1370     args: -run_type full -dm_refine 3 -bc_type dirichlet -petscspace_degree 1 -ksp_type gmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -pc_type asm -pc_asm_type restrict -pc_asm_local_type multiplicative -pc_asm_blocks 4 -sub_pc_type sor -snes_monitor_short -ksp_monitor_short -snes_converged_reason ::ascii_info_detail -ksp_converged_reason -snes_view -show_solution 0
1371 
1372   # Full solve simplex: FAS
1373   test:
1374     suffix: fas_newton_0
1375     requires: triangle !single
1376     args: -run_type full -variable_coefficient nonlinear -petscspace_degree 1 -snes_type fas -snes_fas_levels 2 -fas_coarse_pc_type svd -fas_coarse_ksp_rtol 1.0e-10 -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -dm_refine_hierarchy 1 -snes_view -fas_levels_1_snes_type newtonls -fas_levels_1_pc_type svd -fas_levels_1_ksp_rtol 1.0e-10 -fas_levels_1_snes_monitor_short
1377 
1378   test:
1379     suffix: fas_newton_1
1380     requires: triangle !single
1381     args: -run_type full -dm_refine_hierarchy 3 -petscspace_degree 1 -snes_type fas -snes_fas_levels 3 -fas_coarse_pc_type lu -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_snes_linesearch_type basic -fas_levels_ksp_rtol 1.0e-10 -fas_levels_snes_monitor_short
1382     filter: sed -e "s/total number of linear solver iterations=14/total number of linear solver iterations=15/g"
1383 
1384   test:
1385     suffix: fas_ngs_0
1386     requires: triangle !single
1387     args: -run_type full -variable_coefficient nonlinear -petscspace_degree 1 -snes_type fas -snes_fas_levels 2 -fas_coarse_pc_type svd -fas_coarse_ksp_rtol 1.0e-10 -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -dm_refine_hierarchy 1 -snes_view -fas_levels_1_snes_type ngs -fas_levels_1_snes_monitor_short
1388 
1389   test:
1390     suffix: fas_newton_coarse_0
1391     requires: pragmatic triangle
1392     TODO: broken
1393     args: -run_type full -dm_refine 2 -dm_plex_hash_location -variable_coefficient nonlinear -petscspace_degree 1 -snes_type fas -snes_fas_levels 2 -pc_type svd -ksp_rtol 1.0e-10 -fas_coarse_pc_type svd -fas_coarse_ksp_rtol 1.0e-10 -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -dm_coarsen_hierarchy 1 -snes_view -fas_levels_1_snes_type newtonls -fas_levels_1_pc_type svd -fas_levels_1_ksp_rtol 1.0e-10 -fas_levels_1_snes_monitor_short
1394 
1395   test:
1396     suffix: mg_newton_coarse_0
1397     requires: triangle pragmatic
1398     TODO: broken
1399     args: -run_type full -dm_refine 3 -petscspace_degree 1 -snes_monitor_short -ksp_monitor_true_residual -snes_converged_reason ::ascii_info_detail -dm_coarsen_hierarchy 3 -dm_plex_hash_location -snes_view -dm_view -ksp_type richardson -pc_type mg  -pc_mg_levels 4 -snes_atol 1.0e-8 -ksp_atol 1.0e-8 -snes_rtol 0.0 -ksp_rtol 0.0 -ksp_norm_type unpreconditioned -mg_levels_ksp_type gmres -mg_levels_pc_type ilu -mg_levels_ksp_max_it 10
1400 
1401   test:
1402     suffix: mg_newton_coarse_1
1403     requires: triangle pragmatic
1404     TODO: broken
1405     args: -run_type full -dm_refine 5 -petscspace_degree 1 -dm_coarsen_hierarchy 5 -dm_plex_hash_location -dm_plex_separate_marker -dm_plex_coarsen_bd_label marker -dm_plex_remesh_bd -ksp_type richardson -ksp_rtol 1.0e-12 -pc_type mg -pc_mg_levels 3 -mg_levels_ksp_max_it 2 -snes_converged_reason ::ascii_info_detail -snes_monitor -ksp_monitor_true_residual -mg_levels_ksp_monitor_true_residual -dm_view -ksp_view
1406 
1407   test:
1408     suffix: mg_newton_coarse_2
1409     requires: triangle pragmatic
1410     TODO: broken
1411     args: -run_type full -dm_refine 5 -petscspace_degree 1 -dm_coarsen_hierarchy 5 -dm_plex_hash_location -dm_plex_separate_marker -dm_plex_remesh_bd -ksp_type richardson -ksp_rtol 1.0e-12 -pc_type mg -pc_mg_levels 3 -mg_levels_ksp_max_it 2 -snes_converged_reason ::ascii_info_detail -snes_monitor -ksp_monitor_true_residual -mg_levels_ksp_monitor_true_residual -dm_view -ksp_view
1412 
1413   # Full solve tensor
1414   test:
1415     suffix: tensor_plex_2d
1416     args: -run_type test -dm_plex_simplex 0 -bc_type dirichlet -petscspace_degree 1 -dm_refine_hierarchy 2
1417 
1418   test:
1419     suffix: tensor_p4est_2d
1420     requires: p4est
1421     args: -run_type test -dm_plex_simplex 0 -bc_type dirichlet -petscspace_degree 1 -dm_forest_initial_refinement 2 -dm_forest_minimum_refinement 0 -dm_plex_convert_type p4est
1422 
1423   test:
1424     suffix: tensor_plex_3d
1425     args: -run_type test -dm_plex_simplex 0 -bc_type dirichlet -petscspace_degree 1 -dm_plex_dim 3 -dm_refine_hierarchy 1 -dm_plex_box_faces 2,2,2
1426 
1427   test:
1428     suffix: tensor_p4est_3d
1429     requires: p4est
1430     args: -run_type test -dm_plex_simplex 0 -bc_type dirichlet -petscspace_degree 1 -dm_forest_initial_refinement 1 -dm_forest_minimum_refinement 0 -dm_plex_dim 3 -dm_plex_convert_type p8est -dm_plex_box_faces 2,2,2
1431 
1432   test:
1433     suffix: p4est_test_q2_conformal_serial
1434     requires: p4est
1435     args: -run_type test -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2
1436 
1437   test:
1438     suffix: p4est_test_q2_conformal_parallel
1439     requires: p4est
1440     nsize: 7
1441     args: -run_type test -dm_distribute -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -petscpartitioner_type simple
1442 
1443   test:
1444     suffix: p4est_test_q2_conformal_parallel_parmetis
1445     requires: parmetis p4est
1446     nsize: 4
1447     args: -run_type test -dm_distribute -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -petscpartitioner_type parmetis
1448 
1449   test:
1450     suffix: p4est_test_q2_nonconformal_serial
1451     requires: p4est
1452     filter: grep -v "CG or CGNE: variant"
1453     args: -run_type test -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash
1454 
1455   test:
1456     suffix: p4est_test_q2_nonconformal_parallel
1457     requires: p4est
1458     filter: grep -v "CG or CGNE: variant"
1459     nsize: 7
1460     args: -run_type test -dm_distribute -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type simple
1461 
1462   test:
1463     suffix: p4est_test_q2_nonconformal_parallel_parmetis
1464     requires: parmetis p4est
1465     nsize: 4
1466     args: -run_type test -dm_distribute -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type parmetis
1467 
1468   test:
1469     suffix: p4est_exact_q2_conformal_serial
1470     requires: p4est !single !complex !__float128
1471     args: -run_type exact -petscspace_degree 2 -fas_levels_snes_atol 1.e-10 -snes_max_it 1 -snes_type fas -snes_fas_levels 3 -fas_coarse_pc_type none -fas_coarse_ksp_type preonly -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type none -fas_levels_ksp_type preonly -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2
1472 
1473   test:
1474     suffix: p4est_exact_q2_conformal_parallel
1475     requires: p4est !single !complex !__float128
1476     nsize: 4
1477     args: -run_type exact -dm_distribute -petscspace_degree 2 -fas_levels_snes_atol 1.e-10 -snes_max_it 1 -snes_type fas -snes_fas_levels 3 -fas_coarse_pc_type none -fas_coarse_ksp_type preonly -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type none -fas_levels_ksp_type preonly -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2
1478 
1479   test:
1480     suffix: p4est_exact_q2_conformal_parallel_parmetis
1481     requires: parmetis p4est !single
1482     nsize: 4
1483     args: -run_type exact -dm_distribute -petscspace_degree 2 -fas_levels_snes_atol 1.e-10 -snes_max_it 1 -snes_type fas -snes_fas_levels 3 -fas_coarse_pc_type none -fas_coarse_ksp_type preonly -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type none -fas_levels_ksp_type preonly -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -petscpartitioner_type parmetis
1484 
1485   test:
1486     suffix: p4est_exact_q2_nonconformal_serial
1487     requires: p4est
1488     args: -run_type exact -petscspace_degree 2 -fas_levels_snes_atol 1.e-10 -snes_max_it 1 -snes_type fas -snes_fas_levels 3 -fas_coarse_pc_type none -fas_coarse_ksp_type preonly -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type none -fas_levels_ksp_type preonly -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash
1489 
1490   test:
1491     suffix: p4est_exact_q2_nonconformal_parallel
1492     requires: p4est
1493     nsize: 7
1494     args: -run_type exact -dm_distribute -petscspace_degree 2 -fas_levels_snes_atol 1.e-10 -snes_max_it 1 -snes_type fas -snes_fas_levels 3 -fas_coarse_pc_type none -fas_coarse_ksp_type preonly -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type none -fas_levels_ksp_type preonly -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type simple
1495 
1496   test:
1497     suffix: p4est_exact_q2_nonconformal_parallel_parmetis
1498     requires: parmetis p4est
1499     nsize: 4
1500     args: -run_type exact -dm_distribute -petscspace_degree 2 -fas_levels_snes_atol 1.e-10 -snes_max_it 1 -snes_type fas -snes_fas_levels 3 -fas_coarse_pc_type none -fas_coarse_ksp_type preonly -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type none -fas_levels_ksp_type preonly -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type parmetis
1501 
1502   test:
1503     suffix: p4est_full_q2_nonconformal_serial
1504     requires: p4est !single
1505     filter: grep -v "variant HERMITIAN"
1506     args: -run_type full -petscspace_degree 2 -snes_max_it 20 -snes_type fas -snes_fas_levels 3 -fas_coarse_pc_type jacobi -fas_coarse_ksp_type cg -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type jacobi -fas_levels_ksp_type cg -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash
1507 
1508   test:
1509     suffix: p4est_full_q2_nonconformal_parallel
1510     requires: p4est !single
1511     filter: grep -v "variant HERMITIAN"
1512     nsize: 7
1513     args: -run_type full -dm_distribute -petscspace_degree 2 -snes_max_it 20 -snes_type fas -snes_fas_levels 3 -fas_coarse_pc_type jacobi -fas_coarse_ksp_type cg -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type jacobi -fas_levels_ksp_type cg -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type simple
1514 
1515   test:
1516     suffix: p4est_full_q2_nonconformal_parallel_bddcfas
1517     requires: p4est !single
1518     filter: grep -v "variant HERMITIAN"
1519     nsize: 7
1520     args: -run_type full -dm_distribute -petscspace_degree 2 -snes_max_it 20 -snes_type fas -snes_fas_levels 3 -dm_mat_type is -fas_coarse_pc_type bddc -fas_coarse_ksp_type cg -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type bddc -fas_levels_ksp_type cg -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type simple
1521 
1522   test:
1523     suffix: p4est_full_q2_nonconformal_parallel_bddc
1524     requires: p4est !single
1525     filter: grep -v "variant HERMITIAN"
1526     nsize: 7
1527     args: -run_type full -dm_distribute -petscspace_degree 2 -snes_max_it 20 -snes_type newtonls -dm_mat_type is -pc_type bddc -ksp_type cg -snes_monitor_short -snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type simple
1528 
1529   test:
1530     TODO: broken
1531     suffix: p4est_fas_q2_conformal_serial
1532     requires: p4est !complex !__float128
1533     args: -run_type full -variable_coefficient nonlinear -petscspace_degree 2 -snes_max_it 20 -snes_type fas -snes_fas_levels 3 -pc_type jacobi -ksp_type gmres -fas_coarse_pc_type svd -fas_coarse_ksp_type gmres -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type svd -fas_levels_ksp_type gmres -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_refine_hierarchy 3
1534 
1535   test:
1536     TODO: broken
1537     suffix: p4est_fas_q2_nonconformal_serial
1538     requires: p4est
1539     args: -run_type full -variable_coefficient nonlinear -petscspace_degree 2 -snes_max_it 20 -snes_type fas -snes_fas_levels 3 -pc_type jacobi -ksp_type gmres -fas_coarse_pc_type jacobi -fas_coarse_ksp_type gmres -fas_coarse_ksp_monitor_true_residual -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type jacobi -fas_levels_ksp_type gmres -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash
1540 
1541   test:
1542     suffix: fas_newton_0_p4est
1543     requires: p4est !single !__float128
1544     args: -run_type full -variable_coefficient nonlinear -petscspace_degree 1 -snes_type fas -snes_fas_levels 2 -fas_coarse_pc_type svd -fas_coarse_ksp_rtol 1.0e-10 -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_1_snes_type newtonls -fas_levels_1_pc_type svd -fas_levels_1_ksp_rtol 1.0e-10 -fas_levels_1_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash
1545 
1546   # Full solve simplicial AMR
1547   test:
1548     suffix: tri_p1_adapt_0
1549     requires: pragmatic
1550     TODO: broken
1551     args: -run_type exact -dm_refine 5 -bc_type dirichlet -petscspace_degree 1 -variable_coefficient circle -snes_converged_reason ::ascii_info_detail -pc_type lu -adaptor_refinement_factor 1.0 -dm_view -dm_adapt_view -snes_adapt_initial 1
1552 
1553   test:
1554     suffix: tri_p1_adapt_1
1555     requires: pragmatic
1556     TODO: broken
1557     args: -run_type exact -dm_refine 5 -bc_type dirichlet -petscspace_degree 1 -variable_coefficient circle -snes_converged_reason ::ascii_info_detail -pc_type lu -adaptor_refinement_factor 1.0 -dm_view -dm_adapt_iter_view -dm_adapt_view -snes_adapt_sequence 2
1558 
1559   test:
1560     suffix: tri_p1_adapt_analytic_0
1561     requires: pragmatic
1562     TODO: broken
1563     args: -run_type exact -dm_refine 3 -bc_type dirichlet -petscspace_degree 1 -variable_coefficient cross -snes_adapt_initial 4 -adaptor_target_num 500 -adaptor_monitor -dm_view -dm_adapt_iter_view
1564 
1565   # Full solve tensor AMR
1566   test:
1567     suffix: quad_q1_adapt_0
1568     requires: p4est
1569     args: -run_type exact -dm_plex_simplex 0 -dm_plex_convert_type p4est -bc_type dirichlet -petscspace_degree 1 -variable_coefficient circle -snes_converged_reason ::ascii_info_detail -pc_type lu -dm_forest_initial_refinement 4 -snes_adapt_initial 1 -dm_view
1570     filter: grep -v DM_
1571 
1572   test:
1573     suffix: amr_0
1574     nsize: 5
1575     args: -run_type test -dm_distribute -petscpartitioner_type simple -dm_plex_simplex 0 -bc_type dirichlet -petscspace_degree 1 -dm_refine 1
1576 
1577   test:
1578     suffix: amr_1
1579     requires: p4est !complex
1580     args: -run_type test -dm_plex_simplex 0 -bc_type dirichlet -petscspace_degree 1 -dm_plex_convert_type p4est -dm_p4est_refine_pattern center -dm_forest_maximum_refinement 5 -dm_view vtk:amr.vtu:vtk_vtu -vec_view vtk:amr.vtu:vtk_vtu:append
1581 
1582   test:
1583     suffix: p4est_solve_bddc
1584     requires: p4est !complex
1585     args: -run_type full -dm_distribute -variable_coefficient nonlinear -nonzero_initial_guess 1 -petscspace_degree 2 -snes_max_it 20 -snes_type newtonls -dm_mat_type is -pc_type bddc -ksp_type cg -snes_monitor_short -ksp_monitor -snes_linesearch_type bt -snes_converged_reason -snes_view -dm_plex_simplex 0 -petscspace_poly_tensor -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type simple -pc_bddc_detect_disconnected
1586     nsize: 4
1587 
1588   test:
1589     suffix: p4est_solve_fas
1590     requires: p4est
1591     args: -run_type full -dm_distribute -variable_coefficient nonlinear -nonzero_initial_guess 1 -petscspace_degree 2 -snes_max_it 10 -snes_type fas -snes_linesearch_type bt -snes_fas_levels 3 -fas_coarse_snes_type newtonls -fas_coarse_snes_linesearch_type basic -fas_coarse_ksp_type cg -fas_coarse_pc_type jacobi -fas_coarse_snes_monitor_short -fas_levels_snes_max_it 4 -fas_levels_snes_type newtonls -fas_levels_snes_linesearch_type bt -fas_levels_ksp_type cg -fas_levels_pc_type jacobi -fas_levels_snes_monitor_short -fas_levels_cycle_snes_linesearch_type bt -snes_monitor_short -snes_converged_reason -snes_view -dm_plex_simplex 0 -petscspace_poly_tensor -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash
1592     nsize: 4
1593     TODO: identical machine two runs produce slightly different solver trackers
1594 
1595   test:
1596     suffix: p4est_convergence_test_1
1597     requires: p4est
1598     args:  -quiet -run_type test -dm_distribute -petscspace_degree 1 -dm_plex_simplex 0 -petscspace_poly_tensor -dm_plex_convert_type p4est -dm_forest_minimum_refinement 2 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash
1599     nsize: 4
1600 
1601   test:
1602     suffix: p4est_convergence_test_2
1603     requires: p4est
1604     args: -quiet -run_type test -petscspace_degree 1 -dm_plex_simplex 0 -petscspace_poly_tensor -dm_plex_convert_type p4est -dm_forest_minimum_refinement 3 -dm_forest_initial_refinement 3 -dm_forest_maximum_refinement 5 -dm_p4est_refine_pattern hash
1605 
1606   test:
1607     suffix: p4est_convergence_test_3
1608     requires: p4est
1609     args: -quiet -run_type test -petscspace_degree 1 -dm_plex_simplex 0 -petscspace_poly_tensor -dm_plex_convert_type p4est -dm_forest_minimum_refinement 4 -dm_forest_initial_refinement 4 -dm_forest_maximum_refinement 6 -dm_p4est_refine_pattern hash
1610 
1611   test:
1612     suffix: p4est_convergence_test_4
1613     requires: p4est
1614     args: -quiet -run_type test -petscspace_degree 1 -dm_plex_simplex 0 -petscspace_poly_tensor -dm_plex_convert_type p4est -dm_forest_minimum_refinement 5 -dm_forest_initial_refinement 5 -dm_forest_maximum_refinement 7 -dm_p4est_refine_pattern hash
1615     timeoutfactor: 5
1616 
1617   # Serial tests with GLVis visualization
1618   test:
1619     suffix: glvis_2d_tet_p1
1620     args: -quiet -run_type test -bc_type dirichlet -petscspace_degree 1 -vec_view glvis: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -dm_plex_boundary_label marker -dm_plex_gmsh_periodic 0 -dm_coord_space 0
1621   test:
1622     suffix: glvis_2d_tet_p2
1623     args: -quiet -run_type test -bc_type dirichlet -petscspace_degree 2 -vec_view glvis: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -dm_plex_boundary_label marker -dm_plex_gmsh_periodic 0 -dm_coord_space 0
1624   test:
1625     suffix: glvis_2d_hex_p1
1626     args: -quiet -run_type test -bc_type dirichlet -petscspace_degree 1 -vec_view glvis: -dm_plex_simplex 0 -dm_refine 1 -dm_coord_space 0
1627   test:
1628     suffix: glvis_2d_hex_p2
1629     args: -quiet -run_type test -bc_type dirichlet -petscspace_degree 2 -vec_view glvis: -dm_plex_simplex 0 -dm_refine 1 -dm_coord_space 0
1630   test:
1631     suffix: glvis_2d_hex_p2_p4est
1632     requires: p4est
1633     args: -quiet -run_type test -bc_type dirichlet -petscspace_degree 2 -vec_view glvis: -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 1 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -viewer_glvis_dm_plex_enable_ncmesh
1634   test:
1635     suffix: glvis_2d_tet_p0
1636     args: -run_type exact -guess_vec_view glvis: -nonzero_initial_guess 1 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -dm_plex_boundary_label marker -petscspace_degree 0 -dm_coord_space 0
1637   test:
1638     suffix: glvis_2d_hex_p0
1639     args: -run_type exact -guess_vec_view glvis: -nonzero_initial_guess 1 -dm_plex_box_faces 5,7 -dm_plex_simplex 0 -petscspace_degree 0 -dm_coord_space 0
1640 
1641   # PCHPDDM tests
1642   testset:
1643     nsize: 4
1644     requires: hpddm slepc !single define(PETSC_HAVE_DYNAMIC_LIBRARIES) define(PETSC_USE_SHARED_LIBRARIES)
1645     args: -run_type test -run_test_check_ksp -quiet -petscspace_degree 1 -dm_distribute -petscpartitioner_type simple -bc_type none -dm_plex_simplex 0 -pc_type hpddm -pc_hpddm_levels_1_sub_pc_type lu -pc_hpddm_levels_1_eps_nev 2 -pc_hpddm_coarse_p 1 -pc_hpddm_coarse_pc_type svd -ksp_rtol 1.e-10 -pc_hpddm_levels_1_st_pc_factor_shift_type INBLOCKS -ksp_converged_reason
1646     test:
1647       suffix: quad_singular_hpddm
1648       args: -dm_plex_box_faces 6,7
1649     test:
1650       requires: p4est
1651       suffix: p4est_singular_2d_hpddm
1652       args: -dm_plex_convert_type p4est -dm_forest_minimum_refinement 1 -dm_forest_initial_refinement 3 -dm_forest_maximum_refinement 3
1653     test:
1654       requires: p4est
1655       suffix: p4est_nc_singular_2d_hpddm
1656       args: -dm_plex_convert_type p4est -dm_forest_minimum_refinement 1 -dm_forest_initial_refinement 1 -dm_forest_maximum_refinement 3 -dm_p4est_refine_pattern hash
1657   testset:
1658     nsize: 4
1659     requires: hpddm slepc triangle !single define(PETSC_HAVE_DYNAMIC_LIBRARIES) define(PETSC_USE_SHARED_LIBRARIES)
1660     args: -run_type full -dm_distribute -petscpartitioner_type simple -dm_refine 2 -bc_type dirichlet -petscspace_degree 2 -ksp_type gmres -ksp_gmres_restart 100 -pc_type hpddm -snes_monitor_short -ksp_monitor_short -snes_converged_reason ::ascii_info_detail -ksp_converged_reason -snes_view -show_solution 0 -pc_type hpddm -pc_hpddm_levels_1_sub_pc_type lu -pc_hpddm_levels_1_eps_nev 4 -pc_hpddm_coarse_p 2 -pc_hpddm_coarse_pc_type redundant -ksp_rtol 1.e-1
1661     test:
1662       args: -pc_hpddm_coarse_mat_type baij -options_left no
1663       suffix: tri_hpddm_reuse_baij
1664     test:
1665       requires: !complex
1666       suffix: tri_hpddm_reuse
1667   testset:
1668     nsize: 4
1669     requires: hpddm slepc !single define(PETSC_HAVE_DYNAMIC_LIBRARIES) define(PETSC_USE_SHARED_LIBRARIES)
1670     args: -run_type full -dm_distribute -petscpartitioner_type simple -dm_plex_box_faces 7,5 -dm_refine 2 -dm_plex_simplex 0 -bc_type dirichlet -petscspace_degree 2 -ksp_type gmres -ksp_gmres_restart 100 -pc_type hpddm -snes_monitor_short -ksp_monitor_short -snes_converged_reason ::ascii_info_detail -ksp_converged_reason -snes_view -show_solution 0 -pc_type hpddm -pc_hpddm_levels_1_sub_pc_type lu -pc_hpddm_levels_1_eps_nev 4 -pc_hpddm_coarse_p 2 -pc_hpddm_coarse_pc_type redundant -ksp_rtol 1.e-1
1671     test:
1672       args: -pc_hpddm_coarse_mat_type baij -options_left no
1673       suffix: quad_hpddm_reuse_baij
1674     test:
1675       requires: !complex
1676       suffix: quad_hpddm_reuse
1677   testset:
1678     nsize: 4
1679     requires: hpddm slepc !single define(PETSC_HAVE_DYNAMIC_LIBRARIES) define(PETSC_USE_SHARED_LIBRARIES)
1680     args: -run_type full -dm_distribute -petscpartitioner_type simple -dm_plex_box_faces 7,5 -dm_refine 2 -dm_plex_simplex 0 -bc_type dirichlet -petscspace_degree 1 -ksp_type gmres -ksp_gmres_restart 100 -pc_type hpddm -snes_monitor_short -ksp_monitor_short -snes_converged_reason ::ascii_info_detail -ksp_converged_reason -snes_view -show_solution 0 -pc_type hpddm -pc_hpddm_levels_1_sub_pc_type lu -pc_hpddm_levels_1_eps_threshold 0.1 -pc_hpddm_coarse_p 2 -pc_hpddm_coarse_pc_type redundant -ksp_rtol 1.e-1
1681     test:
1682       args: -pc_hpddm_coarse_mat_type baij -options_left no
1683       suffix: quad_hpddm_reuse_threshold_baij
1684     test:
1685       requires: !complex
1686       suffix: quad_hpddm_reuse_threshold
1687   testset:
1688     nsize: 4
1689     requires: hpddm slepc parmetis !single define(PETSC_HAVE_DYNAMIC_LIBRARIES) define(PETSC_USE_SHARED_LIBRARIES)
1690     filter: sed -e "s/linear solver iterations=17/linear solver iterations=16/g"
1691     args: -run_type full -dm_distribute -petscpartitioner_type parmetis -dm_refine 3 -bc_type dirichlet -petscspace_degree 1 -ksp_type gmres -ksp_gmres_restart 100 -pc_type hpddm -snes_monitor_short -snes_converged_reason ::ascii_info_detail -snes_view -show_solution 0 -pc_type hpddm -pc_hpddm_levels_1_sub_pc_type icc -pc_hpddm_levels_1_eps_nev 20 -pc_hpddm_coarse_p 2 -pc_hpddm_coarse_pc_type redundant -ksp_rtol 1.e-10 -dm_plex_filename ${PETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -dm_plex_boundary_label marker -pc_hpddm_levels_1_sub_pc_factor_levels 3 -variable_coefficient circle -dm_plex_gmsh_periodic 0
1692     test:
1693       args: -pc_hpddm_coarse_mat_type baij -options_left no
1694       filter: grep -v "      total: nonzeros=" | grep -v "      rows=" | sed -e "s/total number of linear solver iterations=1[5-7]/total number of linear solver iterations=16/g"
1695       suffix: tri_parmetis_hpddm_baij
1696     test:
1697       filter: grep -v "      total: nonzeros=" | grep -v "      rows=" | sed -e "s/total number of linear solver iterations=1[5-7]/total number of linear solver iterations=16/g"
1698       requires: !complex
1699       suffix: tri_parmetis_hpddm
1700 
1701   # 2D serial P1 tests for adaptive MG
1702   test:
1703     suffix: 2d_p1_adaptmg_0
1704     requires: triangle bamg
1705     args: -dm_refine_hierarchy 3 -dm_plex_box_faces 4,4 -bc_type dirichlet -petscspace_degree 1 \
1706           -variable_coefficient checkerboard_0 -mat_petscspace_degree 0 -div 16 -k 3 \
1707           -snes_max_it 1 -ksp_converged_reason \
1708           -ksp_rtol 1e-8 -pc_type mg
1709   # -ksp_monitor_true_residual -ksp_converged_reason -mg_levels_ksp_monitor_true_residual -pc_mg_mesp_monitor -dm_adapt_interp_view_fine draw -dm_adapt_interp_view_coarse draw -draw_pause 1
1710   test:
1711     suffix: 2d_p1_adaptmg_1
1712     requires: triangle bamg
1713     args: -dm_refine_hierarchy 3 -dm_plex_box_faces 4,4 -bc_type dirichlet -petscspace_degree 1 \
1714           -variable_coefficient checkerboard_0 -mat_petscspace_degree 0 -div 16 -k 3 \
1715           -snes_max_it 1 -ksp_converged_reason \
1716           -ksp_rtol 1e-8 -pc_type mg -pc_mg_galerkin -pc_mg_adapt_interp -pc_mg_adapt_interp_coarse_space eigenvector -pc_mg_adapt_interp_n 1 \
1717             -pc_mg_mesp_ksp_type richardson -pc_mg_mesp_ksp_richardson_self_scale -pc_mg_mesp_ksp_max_it 100 -pc_mg_mesp_pc_type none
1718 
1719 TEST*/
1720