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