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