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