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