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