Actual source code: ex12.c

petsc-3.12.5 2020-03-29
Report Typos and Errors
  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";

  7: /*
  8: A visualization of the adaptation can be accomplished using:

 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

 12: Information on refinement:

 14:    -info -info_exclude null,sys,vec,is,mat,ksp,snes,ts
 15: */

 17:  #include <petscdmplex.h>
 18:  #include <petscdmadaptor.h>
 19:  #include <petscsnes.h>
 20:  #include <petscds.h>
 21: #include <petscviewerhdf5.h>

 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;

 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;

 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: }

 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: }

 69: /*
 70:   In 2D for Dirichlet conditions, we use exact solution:

 72:     u = x^2 + y^2
 73:     f = 4

 75:   so that

 77:     -\Delta u + f = -4 + 4 = 0

 79:   For Neumann conditions, we have

 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)

 86:   Which we can express as

 88:     \nabla u \cdot  \hat n|_\Gamma = {2 x, 2 y} \cdot \hat n = 2 (x + y)

 90:   The boundary integral of this solution is (assuming we are not orienting the edges)

 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: }

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: }

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);

115:   *u = PetscTanhScalar(xi) + 1.0;
116:   return 0;
117: }

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);

124:   *u = PetscSinReal(alpha*xy) * (alpha*PetscAbsReal(xy) < 2*PETSC_PI ? 1 : 0.01);
125:   return 0;
126: }

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: }

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);

146:   f0[0] = (-4.0*alpha - 8.0*PetscSqr(alpha)*r2*PetscTanhReal(xi)) * PetscSqr(1.0/PetscCoshReal(xi));
147: }

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);

157:   f0[0] = PetscSinReal(alpha*xy) * (alpha*PetscAbsReal(xy) < 2*PETSC_PI ? 1 : 0.01);
158: }

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: }

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: }

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: }

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: }

199: /*
200:   In 2D for x periodicity and y Dirichlet conditions, we use exact solution:

202:     u = sin(2 pi x)
203:     f = -4 pi^2 sin(2 pi x)

205:   so that

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: }

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: }

223: /*
224:   In 2D for x-y periodicity, we use exact solution:

226:     u = sin(2 pi x) sin(2 pi y)
227:     f = -8 pi^2 sin(2 pi x)

229:   so that

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: }

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: }

247: /*
248:   In 2D for Dirichlet conditions with a variable coefficient, we use exact solution:

250:     u  = x^2 + y^2
251:     f  = 6 (x + y)
252:     nu = (x + y)

254:   so that

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: }

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: }

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: }

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: }

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: }

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: }

311: /*
312:   In 2D for Dirichlet conditions with a nonlinear coefficient (p-Laplacian with p = 4), we use exact solution:

314:     u  = x^2 + y^2
315:     f  = 16 (x^2 + y^2)
316:     nu = 1/2 |grad u|^2

318:   so that

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: }

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: }

342: /*
343:   grad (u + eps w) - grad u = eps grad w

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: }

366: /*
367:   In 3D for Dirichlet conditions we use exact solution:

369:     u = 2/3 (x^2 + y^2 + z^2)
370:     f = 4

372:   so that

374:     -\Delta u + f = -2/3 * 6 + 4 = 0

376:   For Neumann conditions, we have

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)

385:   Which we can express as

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: }

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: }

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: }

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;

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_FALSE;
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;

447:   PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX");
448:   PetscOptionsInt("-debug", "The debugging level", "ex12.c", options->debug, &options->debug, NULL);
449:   run  = options->runType;
450:   PetscOptionsEList("-run_type", "The run type", "ex12.c", runTypes, 4, runTypes[options->runType], &run, NULL);

452:   options->runType = (RunType) run;

454:   PetscOptionsInt("-dim", "The topological mesh dimension", "ex12.c", options->dim, &options->dim, NULL);
455:   bd = options->periodicity[0];
456:   PetscOptionsEList("-x_periodicity", "The x-boundary periodicity", "ex12.c", DMBoundaryTypes, 5, DMBoundaryTypes[options->periodicity[0]], &bd, NULL);
457:   options->periodicity[0] = (DMBoundaryType) bd;
458:   bd = options->periodicity[1];
459:   PetscOptionsEList("-y_periodicity", "The y-boundary periodicity", "ex12.c", DMBoundaryTypes, 5, DMBoundaryTypes[options->periodicity[1]], &bd, NULL);
460:   options->periodicity[1] = (DMBoundaryType) bd;
461:   bd = options->periodicity[2];
462:   PetscOptionsEList("-z_periodicity", "The z-boundary periodicity", "ex12.c", DMBoundaryTypes, 5, DMBoundaryTypes[options->periodicity[2]], &bd, NULL);
463:   options->periodicity[2] = (DMBoundaryType) bd;
464:   n = 3;
465:   PetscOptionsIntArray("-cells", "The initial mesh division", "ex12.c", options->cells, &n, NULL);
466:   PetscOptionsString("-f", "Mesh filename to read", "ex12.c", options->filename, options->filename, sizeof(options->filename), &flg);
467:   PetscOptionsBool("-interpolate", "Generate intermediate mesh elements", "ex12.c", options->interpolate, &options->interpolate, NULL);
468:   PetscOptionsReal("-refinement_limit", "The largest allowable cell volume", "ex12.c", options->refinementLimit, &options->refinementLimit, NULL);
469:   bc   = options->bcType;
470:   PetscOptionsEList("-bc_type","Type of boundary condition","ex12.c",bcTypes,3,bcTypes[options->bcType],&bc,NULL);
471:   options->bcType = (BCType) bc;
472:   coeff = options->variableCoefficient;
473:   PetscOptionsEList("-variable_coefficient","Type of variable coefficent","ex12.c",coeffTypes,6,coeffTypes[options->variableCoefficient],&coeff,NULL);
474:   options->variableCoefficient = (CoeffType) coeff;

476:   PetscOptionsBool("-field_bc", "Use a field representation for the BC", "ex12.c", options->fieldBC, &options->fieldBC, NULL);
477:   PetscOptionsBool("-jacobian_mf", "Calculate the action of the Jacobian on the fly", "ex12.c", options->jacobianMF, &options->jacobianMF, NULL);
478:   PetscOptionsBool("-show_initial", "Output the initial guess for verification", "ex12.c", options->showInitial, &options->showInitial, NULL);
479:   PetscOptionsBool("-show_solution", "Output the solution for verification", "ex12.c", options->showSolution, &options->showSolution, NULL);
480:   PetscOptionsBool("-restart", "Read in the mesh and solution from a file", "ex12.c", options->restart, &options->restart, NULL);
481:   PetscOptionsBool("-dm_view_hierarchy", "View the coarsened hierarchy", "ex12.c", options->viewHierarchy, &options->viewHierarchy, NULL);
482:   PetscOptionsBool("-simplex", "Simplicial (true) or tensor (false) mesh", "ex12.c", options->simplex, &options->simplex, NULL);
483:   PetscOptionsBool("-quiet", "Don't print any vecs", "ex12.c", options->quiet, &options->quiet, NULL);
484:   PetscOptionsBool("-nonzero_initial_guess", "nonzero intial guess", "ex12.c", options->nonzInit, &options->nonzInit, NULL);
485:   PetscOptionsBool("-bd_integral", "Compute the integral of the solution on the boundary", "ex12.c", options->bdIntegral, &options->bdIntegral, NULL);
486:   if (options->runType == RUN_TEST) {
487:     PetscOptionsBool("-run_test_check_ksp", "Check solution of KSP", "ex12.c", options->checkksp, &options->checkksp, NULL);
488:   }
489:   PetscOptionsEnd();
490:   PetscLogEventRegister("CreateMesh", DM_CLASSID, &options->createMeshEvent);
491:   return(0);
492: }

494: static PetscErrorCode CreateBCLabel(DM dm, const char name[])
495: {
496:   DMLabel        label;

500:   DMCreateLabel(dm, name);
501:   DMGetLabel(dm, name, &label);
502:   DMPlexMarkBoundaryFaces(dm, 1, label);
503:   DMPlexLabelComplete(dm, label);
504:   return(0);
505: }

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;

517:   PetscLogEventBegin(user->createMeshEvent,0,0,0,0);
518:   PetscStrlen(filename, &len);
519:   if (!len) {
520:     PetscInt d;

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:     DMPlexCreateBoxMesh(comm, dim, user->simplex, user->cells, NULL, NULL, user->periodicity, interpolate, dm);
524:     PetscObjectSetName((PetscObject) *dm, "Mesh");
525:   } else {
526:     DMPlexCreateFromFile(comm, filename, interpolate, dm);
527:     DMPlexSetRefinementUniform(*dm, PETSC_FALSE);
528:   }
529:   {
530:     PetscPartitioner part;
531:     DM               refinedMesh     = NULL;
532:     DM               distributedMesh = NULL;

534:     /* Refine mesh using a volume constraint */
535:     if (refinementLimit > 0.0) {
536:       DMPlexSetRefinementLimit(*dm, refinementLimit);
537:       DMRefine(*dm, comm, &refinedMesh);
538:       if (refinedMesh) {
539:         const char *name;

541:         PetscObjectGetName((PetscObject) *dm,         &name);
542:         PetscObjectSetName((PetscObject) refinedMesh,  name);
543:         DMDestroy(dm);
544:         *dm  = refinedMesh;
545:       }
546:     }
547:     /* Distribute mesh over processes */
548:     DMPlexGetPartitioner(*dm,&part);
549:     PetscPartitionerSetFromOptions(part);
550:     DMPlexDistribute(*dm, 0, NULL, &distributedMesh);
551:     if (distributedMesh) {
552:       DMDestroy(dm);
553:       *dm  = distributedMesh;
554:     }
555:   }
556:   if (user->bcType == NEUMANN) {
557:     DMLabel   label;

559:     DMCreateLabel(*dm, "boundary");
560:     DMGetLabel(*dm, "boundary", &label);
561:     DMPlexMarkBoundaryFaces(*dm, 1, label);
562:   } else if (user->bcType == DIRICHLET) {
563:     PetscBool hasLabel;

565:     DMHasLabel(*dm,"marker",&hasLabel);
566:     if (!hasLabel) {CreateBCLabel(*dm, "marker");}
567:   }
568:   {
569:     char      convType[256];
570:     PetscBool flg;

572:     PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX");
573:     PetscOptionsFList("-dm_plex_convert_type","Convert DMPlex to another format","ex12",DMList,DMPLEX,convType,256,&flg);
574:     PetscOptionsEnd();
575:     if (flg) {
576:       DM dmConv;

578:       DMConvert(*dm,convType,&dmConv);
579:       if (dmConv) {
580:         DMDestroy(dm);
581:         *dm  = dmConv;
582:       }
583:     }
584:   }
585:   DMLocalizeCoordinates(*dm); /* needed for periodic */
586:   DMSetFromOptions(*dm);
587:   DMViewFromOptions(*dm, NULL, "-dm_view");
588:   if (user->viewHierarchy) {
589:     DM       cdm = *dm;
590:     PetscInt i   = 0;
591:     char     buf[256];

593:     while (cdm) {
594:       DMSetUp(cdm);
595:       DMGetCoarseDM(cdm, &cdm);
596:       ++i;
597:     }
598:     cdm = *dm;
599:     while (cdm) {
600:       PetscViewer       viewer;
601:       PetscBool   isHDF5, isVTK;

603:       --i;
604:       PetscViewerCreate(comm,&viewer);
605:       PetscViewerSetType(viewer,PETSCVIEWERHDF5);
606:       PetscViewerSetOptionsPrefix(viewer,"hierarchy_");
607:       PetscViewerSetFromOptions(viewer);
608:       PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&isHDF5);
609:       PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isVTK);
610:       if (isHDF5) {
611:         PetscSNPrintf(buf, 256, "ex12-%d.h5", i);
612:       } else if (isVTK) {
613:         PetscSNPrintf(buf, 256, "ex12-%d.vtu", i);
614:         PetscViewerPushFormat(viewer,PETSC_VIEWER_VTK_VTU);
615:       } else {
616:         PetscSNPrintf(buf, 256, "ex12-%d", i);
617:       }
618:       PetscViewerFileSetMode(viewer,FILE_MODE_WRITE);
619:       PetscViewerFileSetName(viewer,buf);
620:       DMView(cdm, viewer);
621:       PetscViewerDestroy(&viewer);
622:       DMGetCoarseDM(cdm, &cdm);
623:     }
624:   }
625:   PetscLogEventEnd(user->createMeshEvent,0,0,0,0);
626:   return(0);
627: }

629: static PetscErrorCode SetupProblem(DM dm, AppCtx *user)
630: {
631:   PetscDS        prob;
632:   const PetscInt id = 1;

636:   DMGetDS(dm, &prob);
637:   switch (user->variableCoefficient) {
638:   case COEFF_NONE:
639:     if (user->periodicity[0]) {
640:       if (user->periodicity[1]) {
641:         PetscDSSetResidual(prob, 0, f0_xytrig_u, f1_u);
642:         PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu);
643:       } else {
644:         PetscDSSetResidual(prob, 0, f0_xtrig_u,  f1_u);
645:         PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu);
646:       }
647:     } else {
648:       PetscDSSetResidual(prob, 0, f0_u, f1_u);
649:       PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu);
650:     }
651:     break;
652:   case COEFF_ANALYTIC:
653:     PetscDSSetResidual(prob, 0, f0_analytic_u, f1_analytic_u);
654:     PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_analytic_uu);
655:     break;
656:   case COEFF_FIELD:
657:     PetscDSSetResidual(prob, 0, f0_analytic_u, f1_field_u);
658:     PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_field_uu);
659:     break;
660:   case COEFF_NONLINEAR:
661:     PetscDSSetResidual(prob, 0, f0_analytic_nonlinear_u, f1_analytic_nonlinear_u);
662:     PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_analytic_nonlinear_uu);
663:     break;
664:   case COEFF_CIRCLE:
665:     PetscDSSetResidual(prob, 0, f0_circle_u, f1_u);
666:     PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu);
667:     break;
668:   case COEFF_CROSS:
669:     PetscDSSetResidual(prob, 0, f0_cross_u, f1_u);
670:     PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu);
671:     break;
672:   default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid variable coefficient type %d", user->variableCoefficient);
673:   }
674:   switch (user->dim) {
675:   case 2:
676:     switch (user->variableCoefficient) {
677:     case COEFF_CIRCLE:
678:       user->exactFuncs[0]  = circle_u_2d;break;
679:     case COEFF_CROSS:
680:       user->exactFuncs[0]  = cross_u_2d;break;
681:     default:
682:       if (user->periodicity[0]) {
683:         if (user->periodicity[1]) {
684:           user->exactFuncs[0] = xytrig_u_2d;
685:         } else {
686:           user->exactFuncs[0] = xtrig_u_2d;
687:         }
688:       } else {
689:         user->exactFuncs[0]  = quadratic_u_2d;
690:         user->exactFields[0] = quadratic_u_field_2d;
691:       }
692:     }
693:     if (user->bcType == NEUMANN) {PetscDSSetBdResidual(prob, 0, f0_bd_u, f1_bd_zero);}
694:     break;
695:   case 3:
696:     user->exactFuncs[0]  = quadratic_u_3d;
697:     user->exactFields[0] = quadratic_u_field_3d;
698:     if (user->bcType == NEUMANN) {PetscDSSetBdResidual(prob, 0, f0_bd_u, f1_bd_zero);}
699:     break;
700:   default:
701:     SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %d", user->dim);
702:   }
703:   if (user->bcType != NONE) {
704:     PetscDSAddBoundary(prob, user->bcType == DIRICHLET ? (user->fieldBC ? DM_BC_ESSENTIAL_FIELD : DM_BC_ESSENTIAL) : DM_BC_NATURAL,
705:                               "wall", user->bcType == DIRICHLET ? "marker" : "boundary", 0, 0, NULL,
706:                               user->fieldBC ? (void (*)(void)) user->exactFields[0] : (void (*)(void)) user->exactFuncs[0], 1, &id, user);
707:   }
708:   PetscDSSetExactSolution(prob, 0, user->exactFuncs[0], user);
709:   return(0);
710: }

712: static PetscErrorCode SetupMaterial(DM dm, DM dmAux, AppCtx *user)
713: {
714:   PetscErrorCode (*matFuncs[1])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx) = {nu_2d};
715:   Vec            nu;

719:   DMCreateLocalVector(dmAux, &nu);
720:   DMProjectFunctionLocal(dmAux, 0.0, matFuncs, NULL, INSERT_ALL_VALUES, nu);
721:   PetscObjectCompose((PetscObject) dm, "A", (PetscObject) nu);
722:   VecDestroy(&nu);
723:   return(0);
724: }

726: static PetscErrorCode SetupBC(DM dm, DM dmAux, AppCtx *user)
727: {
728:   PetscErrorCode (*bcFuncs[1])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx);
729:   Vec            uexact;
730:   PetscInt       dim;

734:   DMGetDimension(dm, &dim);
735:   if (dim == 2) bcFuncs[0] = quadratic_u_2d;
736:   else          bcFuncs[0] = quadratic_u_3d;
737:   DMCreateLocalVector(dmAux, &uexact);
738:   DMProjectFunctionLocal(dmAux, 0.0, bcFuncs, NULL, INSERT_ALL_VALUES, uexact);
739:   PetscObjectCompose((PetscObject) dm, "A", (PetscObject) uexact);
740:   VecDestroy(&uexact);
741:   return(0);
742: }

744: static PetscErrorCode SetupAuxDM(DM dm, PetscFE feAux, AppCtx *user)
745: {
746:   DM             dmAux, coordDM;

750:   /* MUST call DMGetCoordinateDM() in order to get p4est setup if present */
751:   DMGetCoordinateDM(dm, &coordDM);
752:   if (!feAux) return(0);
753:   DMClone(dm, &dmAux);
754:   PetscObjectCompose((PetscObject) dm, "dmAux", (PetscObject) dmAux);
755:   DMSetCoordinateDM(dmAux, coordDM);
756:   DMSetField(dmAux, 0, NULL, (PetscObject) feAux);
757:   DMCreateDS(dmAux);
758:   if (user->fieldBC) {SetupBC(dm, dmAux, user);}
759:   else               {SetupMaterial(dm, dmAux, user);}
760:   DMDestroy(&dmAux);
761:   return(0);
762: }

764: static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
765: {
766:   DM             cdm = dm;
767:   const PetscInt dim = user->dim;
768:   PetscFE        fe, feAux = NULL;
769:   PetscBool      simplex   = user->simplex;
770:   MPI_Comm       comm;

774:   /* Create finite element for each field and auxiliary field */
775:   PetscObjectGetComm((PetscObject) dm, &comm);
776:   PetscFECreateDefault(comm, dim, 1, simplex, NULL, -1, &fe);
777:   PetscObjectSetName((PetscObject) fe, "potential");
778:   if (user->variableCoefficient == COEFF_FIELD) {
779:     PetscFECreateDefault(comm, dim, 1, simplex, "mat_", -1, &feAux);
780:     PetscFECopyQuadrature(fe, feAux);
781:   } else if (user->fieldBC) {
782:     PetscFECreateDefault(comm, dim, 1, simplex, "bc_", -1, &feAux);
783:     PetscFECopyQuadrature(fe, feAux);
784:   }
785:   /* Set discretization and boundary conditions for each mesh */
786:   DMSetField(dm, 0, NULL, (PetscObject) fe);
787:   DMCreateDS(dm);
788:   SetupProblem(dm, user);
789:   while (cdm) {
790:     DMCopyDisc(dm, cdm);
791:     SetupAuxDM(cdm, feAux, user);
792:     if (user->bcType == DIRICHLET) {
793:       PetscBool hasLabel;

795:       DMHasLabel(cdm, "marker", &hasLabel);
796:       if (!hasLabel) {CreateBCLabel(cdm, "marker");}
797:     }
798:     DMGetCoarseDM(cdm, &cdm);
799:   }
800:   PetscFEDestroy(&fe);
801:   PetscFEDestroy(&feAux);
802:   return(0);
803: }

805: #include "petsc/private/petscimpl.h"

807: /*@C
808:   KSPMonitorError - Outputs the error at each iteration of an iterative solver.

810:   Collective on KSP

812:   Input Parameters:
813: + ksp   - the KSP
814: . its   - iteration number
815: . rnorm - 2-norm, preconditioned residual value (may be estimated).
816: - ctx   - monitor context

818:   Level: intermediate

820: .seealso: KSPMonitorSet(), KSPMonitorTrueResidualNorm(), KSPMonitorDefault()
821: @*/
822: static PetscErrorCode KSPMonitorError(KSP ksp, PetscInt its, PetscReal rnorm, void *ctx)
823: {
824:   AppCtx        *user = (AppCtx *) ctx;
825:   DM             dm;
826:   Vec            du = NULL, r;
827:   PetscInt       level = 0;
828:   PetscBool      hasLevel;
829: #if defined(PETSC_HAVE_HDF5)
830:   PetscViewer    viewer;
831:   char           buf[256];
832: #endif

836:   KSPGetDM(ksp, &dm);
837:   /* Calculate solution */
838:   {
839:     PC        pc = user->pcmg; /* The MG PC */
840:     DM        fdm = NULL,  cdm = NULL;
841:     KSP       fksp, cksp;
842:     Vec       fu,   cu = NULL;
843:     PetscInt  levels, l;

845:     KSPBuildSolution(ksp, NULL, &du);
846:     PetscObjectComposedDataGetInt((PetscObject) ksp, PetscMGLevelId, level, hasLevel);
847:     PCMGGetLevels(pc, &levels);
848:     PCMGGetSmoother(pc, levels-1, &fksp);
849:     KSPBuildSolution(fksp, NULL, &fu);
850:     for (l = levels-1; l > level; --l) {
851:       Mat R;
852:       Vec s;

854:       PCMGGetSmoother(pc, l-1, &cksp);
855:       KSPGetDM(cksp, &cdm);
856:       DMGetGlobalVector(cdm, &cu);
857:       PCMGGetRestriction(pc, l, &R);
858:       PCMGGetRScale(pc, l, &s);
859:       MatRestrict(R, fu, cu);
860:       VecPointwiseMult(cu, cu, s);
861:       if (l < levels-1) {DMRestoreGlobalVector(fdm, &fu);}
862:       fdm  = cdm;
863:       fu   = cu;
864:     }
865:     if (levels-1 > level) {
866:       VecAXPY(du, 1.0, cu);
867:       DMRestoreGlobalVector(cdm, &cu);
868:     }
869:   }
870:   /* Calculate error */
871:   DMGetGlobalVector(dm, &r);
872:   DMProjectFunction(dm, 0.0, user->exactFuncs, NULL, INSERT_ALL_VALUES, r);
873:   VecAXPY(r,-1.0,du);
874:   PetscObjectSetName((PetscObject) r, "solution error");
875:   /* View error */
876: #if defined(PETSC_HAVE_HDF5)
877:   PetscSNPrintf(buf, 256, "ex12-%D.h5", level);
878:   PetscViewerHDF5Open(PETSC_COMM_WORLD, buf, FILE_MODE_APPEND, &viewer);
879:   VecView(r, viewer);
880:   PetscViewerDestroy(&viewer);
881: #endif
882:   DMRestoreGlobalVector(dm, &r);
883:   return(0);
884: }

886: /*@C
887:   SNESMonitorError - Outputs the error at each iteration of an iterative solver.

889:   Collective on SNES

891:   Input Parameters:
892: + snes  - the SNES
893: . its   - iteration number
894: . rnorm - 2-norm of residual
895: - ctx   - user context

897:   Level: intermediate

899: .seealso: SNESMonitorDefault(), SNESMonitorSet(), SNESMonitorSolution()
900: @*/
901: static PetscErrorCode SNESMonitorError(SNES snes, PetscInt its, PetscReal rnorm, void *ctx)
902: {
903:   AppCtx        *user = (AppCtx *) ctx;
904:   DM             dm;
905:   Vec            u, r;
906:   PetscInt       level = -1;
907:   PetscBool      hasLevel;
908: #if defined(PETSC_HAVE_HDF5)
909:   PetscViewer    viewer;
910: #endif
911:   char           buf[256];

915:   SNESGetDM(snes, &dm);
916:   /* Calculate error */
917:   SNESGetSolution(snes, &u);
918:   DMGetGlobalVector(dm, &r);
919:   PetscObjectSetName((PetscObject) r, "solution error");
920:   DMProjectFunction(dm, 0.0, user->exactFuncs, NULL, INSERT_ALL_VALUES, r);
921:   VecAXPY(r, -1.0, u);
922:   /* View error */
923:   PetscObjectComposedDataGetInt((PetscObject) snes, PetscMGLevelId, level, hasLevel);
924:   PetscSNPrintf(buf, 256, "ex12-%D.h5", level);
925: #if defined(PETSC_HAVE_HDF5)
926:   PetscViewerHDF5Open(PETSC_COMM_WORLD, buf, FILE_MODE_APPEND, &viewer);
927:   VecView(r, viewer);
928:   PetscViewerDestroy(&viewer);
929:   /* Cleanup */
930:   DMRestoreGlobalVector(dm, &r);
931:   return(0);
932: #else
933:   SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"You need to configure with --download-hdf5");
934: #endif
935: }

937: int main(int argc, char **argv)
938: {
939:   DM             dm;          /* Problem specification */
940:   SNES           snes;        /* nonlinear solver */
941:   Vec            u;           /* solution vector */
942:   Mat            A,J;         /* Jacobian matrix */
943:   MatNullSpace   nullSpace;   /* May be necessary for Neumann conditions */
944:   AppCtx         user;        /* user-defined work context */
945:   JacActionCtx   userJ;       /* context for Jacobian MF action */
946:   PetscReal      error = 0.0; /* L_2 error in the solution */
947:   PetscBool      isFAS;

950:   PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
951:   ProcessOptions(PETSC_COMM_WORLD, &user);
952:   SNESCreate(PETSC_COMM_WORLD, &snes);
953:   CreateMesh(PETSC_COMM_WORLD, &user, &dm);
954:   SNESSetDM(snes, dm);
955:   DMSetApplicationContext(dm, &user);

957:   PetscMalloc2(1, &user.exactFuncs, 1, &user.exactFields);
958:   SetupDiscretization(dm, &user);

960:   DMCreateGlobalVector(dm, &u);
961:   PetscObjectSetName((PetscObject) u, "potential");

963:   DMCreateMatrix(dm, &J);
964:   if (user.jacobianMF) {
965:     PetscInt M, m, N, n;

967:     MatGetSize(J, &M, &N);
968:     MatGetLocalSize(J, &m, &n);
969:     MatCreate(PETSC_COMM_WORLD, &A);
970:     MatSetSizes(A, m, n, M, N);
971:     MatSetType(A, MATSHELL);
972:     MatSetUp(A);
973: #if 0
974:     MatShellSetOperation(A, MATOP_MULT, (void (*)(void))FormJacobianAction);
975: #endif

977:     userJ.dm   = dm;
978:     userJ.J    = J;
979:     userJ.user = &user;

981:     DMCreateLocalVector(dm, &userJ.u);
982:     if (user.fieldBC) {DMProjectFieldLocal(dm, 0.0, userJ.u, user.exactFields, INSERT_BC_VALUES, userJ.u);}
983:     else              {DMProjectFunctionLocal(dm, 0.0, user.exactFuncs, NULL, INSERT_BC_VALUES, userJ.u);}
984:     MatShellSetContext(A, &userJ);
985:   } else {
986:     A = J;
987:   }

989:   nullSpace = NULL;
990:   if (user.bcType != DIRICHLET) {
991:     MatNullSpaceCreate(PetscObjectComm((PetscObject) dm), PETSC_TRUE, 0, NULL, &nullSpace);
992:     MatSetNullSpace(A, nullSpace);
993:   }

995:   DMPlexSetSNESLocalFEM(dm,&user,&user,&user);
996:   SNESSetJacobian(snes, A, J, NULL, NULL);

998:   SNESSetFromOptions(snes);

1000:   if (user.fieldBC) {DMProjectField(dm, 0.0, u, user.exactFields, INSERT_ALL_VALUES, u);}
1001:   else              {DMProjectFunction(dm, 0.0, user.exactFuncs, NULL, INSERT_ALL_VALUES, u);}
1002:   if (user.restart) {
1003: #if defined(PETSC_HAVE_HDF5)
1004:     PetscViewer viewer;

1006:     PetscViewerCreate(PETSC_COMM_WORLD, &viewer);
1007:     PetscViewerSetType(viewer, PETSCVIEWERHDF5);
1008:     PetscViewerFileSetMode(viewer, FILE_MODE_READ);
1009:     PetscViewerFileSetName(viewer, user.filename);
1010:     PetscViewerHDF5PushGroup(viewer, "/fields");
1011:     VecLoad(u, viewer);
1012:     PetscViewerHDF5PopGroup(viewer);
1013:     PetscViewerDestroy(&viewer);
1014: #endif
1015:   }
1016:   if (user.showInitial) {
1017:     Vec lv;
1018:     DMGetLocalVector(dm, &lv);
1019:     DMGlobalToLocalBegin(dm, u, INSERT_VALUES, lv);
1020:     DMGlobalToLocalEnd(dm, u, INSERT_VALUES, lv);
1021:     DMPrintLocalVec(dm, "Local function", 1.0e-10, lv);
1022:     DMRestoreLocalVector(dm, &lv);
1023:   }
1024:   if (user.viewHierarchy) {
1025:     SNES      lsnes;
1026:     KSP       ksp;
1027:     PC        pc;
1028:     PetscInt  numLevels, l;
1029:     PetscBool isMG;

1031:     PetscObjectTypeCompare((PetscObject) snes, SNESFAS, &isFAS);
1032:     if (isFAS) {
1033:       SNESFASGetLevels(snes, &numLevels);
1034:       for (l = 0; l < numLevels; ++l) {
1035:         SNESFASGetCycleSNES(snes, l, &lsnes);
1036:         SNESMonitorSet(lsnes, SNESMonitorError, &user, NULL);
1037:       }
1038:     } else {
1039:       SNESGetKSP(snes, &ksp);
1040:       KSPGetPC(ksp, &pc);
1041:       PetscObjectTypeCompare((PetscObject) pc, PCMG, &isMG);
1042:       if (isMG) {
1043:         user.pcmg = pc;
1044:         PCMGGetLevels(pc, &numLevels);
1045:         for (l = 0; l < numLevels; ++l) {
1046:           PCMGGetSmootherDown(pc, l, &ksp);
1047:           KSPMonitorSet(ksp, KSPMonitorError, &user, NULL);
1048:         }
1049:       }
1050:     }
1051:   }
1052:   if (user.runType == RUN_FULL || user.runType == RUN_EXACT) {
1053:     PetscErrorCode (*initialGuess[1])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx) = {zero};

1055:     if (user.nonzInit) initialGuess[0] = ecks;
1056:     if (user.runType == RUN_FULL) {
1057:       DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u);
1058:     }
1059:     if (user.debug) {
1060:       PetscPrintf(PETSC_COMM_WORLD, "Initial guess\n");
1061:       VecView(u, PETSC_VIEWER_STDOUT_WORLD);
1062:     }
1063:     VecViewFromOptions(u, NULL, "-guess_vec_view");
1064:     SNESSolve(snes, NULL, u);
1065:     SNESGetSolution(snes, &u);
1066:     SNESGetDM(snes, &dm);

1068:     if (user.showSolution) {
1069:       PetscPrintf(PETSC_COMM_WORLD, "Solution\n");
1070:       VecChop(u, 3.0e-9);
1071:       VecView(u, PETSC_VIEWER_STDOUT_WORLD);
1072:     }
1073:     VecViewFromOptions(u, NULL, "-vec_view");
1074:   } else if (user.runType == RUN_PERF) {
1075:     Vec       r;
1076:     PetscReal res = 0.0;

1078:     SNESGetFunction(snes, &r, NULL, NULL);
1079:     SNESComputeFunction(snes, u, r);
1080:     PetscPrintf(PETSC_COMM_WORLD, "Initial Residual\n");
1081:     VecChop(r, 1.0e-10);
1082:     VecNorm(r, NORM_2, &res);
1083:     PetscPrintf(PETSC_COMM_WORLD, "L_2 Residual: %g\n", (double)res);
1084:   } else {
1085:     Vec       r;
1086:     PetscReal res = 0.0, tol = 1.0e-11;

1088:     /* Check discretization error */
1089:     SNESGetFunction(snes, &r, NULL, NULL);
1090:     PetscPrintf(PETSC_COMM_WORLD, "Initial guess\n");
1091:     if (!user.quiet) {VecView(u, PETSC_VIEWER_STDOUT_WORLD);}
1092:     DMComputeL2Diff(dm, 0.0, user.exactFuncs, NULL, u, &error);
1093:     if (error < tol) {PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: < %2.1e\n", (double)tol);}
1094:     else             {PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: %g\n", (double)error);}
1095:     /* Check residual */
1096:     SNESComputeFunction(snes, u, r);
1097:     PetscPrintf(PETSC_COMM_WORLD, "Initial Residual\n");
1098:     VecChop(r, 1.0e-10);
1099:     if (!user.quiet) {VecView(r, PETSC_VIEWER_STDOUT_WORLD);}
1100:     VecNorm(r, NORM_2, &res);
1101:     PetscPrintf(PETSC_COMM_WORLD, "L_2 Residual: %g\n", (double)res);
1102:     /* Check Jacobian */
1103:     {
1104:       Vec b;

1106:       SNESComputeJacobian(snes, u, A, A);
1107:       VecDuplicate(u, &b);
1108:       VecSet(r, 0.0);
1109:       SNESComputeFunction(snes, r, b);
1110:       MatMult(A, u, r);
1111:       VecAXPY(r, 1.0, b);
1112:       PetscPrintf(PETSC_COMM_WORLD, "Au - b = Au + F(0)\n");
1113:       VecChop(r, 1.0e-10);
1114:       if (!user.quiet) {VecView(r, PETSC_VIEWER_STDOUT_WORLD);}
1115:       VecNorm(r, NORM_2, &res);
1116:       PetscPrintf(PETSC_COMM_WORLD, "Linear L_2 Residual: %g\n", (double)res);
1117:       /* check solver */
1118:       if (user.checkksp) {
1119:         KSP ksp;

1121:         if (nullSpace) {
1122:           MatNullSpaceRemove(nullSpace, u);
1123:         }
1124:         SNESComputeJacobian(snes, u, A, J);
1125:         MatMult(A, u, b);
1126:         SNESGetKSP(snes, &ksp);
1127:         KSPSetOperators(ksp, A, J);
1128:         KSPSolve(ksp, b, r);
1129:         VecAXPY(r, -1.0, u);
1130:         VecNorm(r, NORM_2, &res);
1131:         PetscPrintf(PETSC_COMM_WORLD, "KSP Error: %g\n", (double)res);
1132:       }
1133:       VecDestroy(&b);
1134:     }
1135:   }
1136:   VecViewFromOptions(u, NULL, "-vec_view");

1138:   if (user.bdIntegral) {
1139:     DMLabel   label;
1140:     PetscInt  id = 1;
1141:     PetscScalar bdInt = 0.0;
1142:     PetscReal   exact = 3.3333333333;

1144:     DMGetLabel(dm, "marker", &label);
1145:     DMPlexComputeBdIntegral(dm, u, label, 1, &id, bd_integral_2d, &bdInt, NULL);
1146:     PetscPrintf(PETSC_COMM_WORLD, "Solution boundary integral: %.4g\n", (double) PetscAbsScalar(bdInt));
1147:     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);
1148:   }

1150:   MatNullSpaceDestroy(&nullSpace);
1151:   if (user.jacobianMF) {VecDestroy(&userJ.u);}
1152:   if (A != J) {MatDestroy(&A);}
1153:   MatDestroy(&J);
1154:   VecDestroy(&u);
1155:   SNESDestroy(&snes);
1156:   DMDestroy(&dm);
1157:   PetscFree2(user.exactFuncs, user.exactFields);
1158:   PetscFinalize();
1159:   return ierr;
1160: }

1162: /*TEST
1163:   # 2D serial P1 test 0-4
1164:   test:
1165:     suffix: 2d_p1_0
1166:     requires: triangle
1167:     args: -run_type test -refinement_limit 0.0    -bc_type dirichlet -interpolate 0 -petscspace_degree 1 -show_initial -dm_plex_print_fem 1

1169:   test:
1170:     suffix: 2d_p1_1
1171:     requires: triangle
1172:     args: -run_type test -refinement_limit 0.0    -bc_type dirichlet -interpolate 1 -petscspace_degree 1 -show_initial -dm_plex_print_fem 1

1174:   test:
1175:     suffix: 2d_p1_2
1176:     requires: triangle
1177:     args: -run_type test -refinement_limit 0.0625 -bc_type dirichlet -interpolate 1 -petscspace_degree 1 -show_initial -dm_plex_print_fem 1

1179:   test:
1180:     suffix: 2d_p1_neumann_0
1181:     requires: triangle
1182:     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

1184:   test:
1185:     suffix: 2d_p1_neumann_1
1186:     requires: triangle
1187:     args: -run_type test -refinement_limit 0.0625 -bc_type neumann   -interpolate 1 -petscspace_degree 1 -show_initial -dm_plex_print_fem 1

1189:   # 2D serial P2 test 5-8
1190:   test:
1191:     suffix: 2d_p2_0
1192:     requires: triangle
1193:     args: -run_type test -refinement_limit 0.0    -bc_type dirichlet -interpolate 1 -petscspace_degree 2 -show_initial -dm_plex_print_fem 1

1195:   test:
1196:     suffix: 2d_p2_1
1197:     requires: triangle
1198:     args: -run_type test -refinement_limit 0.0625 -bc_type dirichlet -interpolate 1 -petscspace_degree 2 -show_initial -dm_plex_print_fem 1

1200:   test:
1201:     suffix: 2d_p2_neumann_0
1202:     requires: triangle
1203:     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

1205:   test:
1206:     suffix: 2d_p2_neumann_1
1207:     requires: triangle
1208:     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

1210:   test:
1211:     suffix: bd_int_0
1212:     requires: triangle
1213:     args: -run_type test -bc_type dirichlet -interpolate 1 -petscspace_degree 2 -bd_integral -dm_view -quiet

1215:   test:
1216:     suffix: bd_int_1
1217:     requires: triangle
1218:     args: -run_type test -dm_refine 2 -bc_type dirichlet -interpolate 1 -petscspace_degree 2 -bd_integral -dm_view -quiet

1220:   # 3D serial P1 test 9-12
1221:   test:
1222:     suffix: 3d_p1_0
1223:     requires: ctetgen
1224:     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

1226:   test:
1227:     suffix: 3d_p1_1
1228:     requires: ctetgen
1229:     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

1231:   test:
1232:     suffix: 3d_p1_2
1233:     requires: ctetgen
1234:     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

1236:   test:
1237:     suffix: 3d_p1_neumann_0
1238:     requires: ctetgen
1239:     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

1241:   # Analytic variable coefficient 13-20
1242:   test:
1243:     suffix: 13
1244:     requires: triangle
1245:     args: -run_type test -refinement_limit 0.0    -variable_coefficient analytic -interpolate 1 -petscspace_degree 1 -show_initial -dm_plex_print_fem 1
1246:   test:
1247:     suffix: 14
1248:     requires: triangle
1249:     args: -run_type test -refinement_limit 0.0625 -variable_coefficient analytic -interpolate 1 -petscspace_degree 1 -show_initial -dm_plex_print_fem 1
1250:   test:
1251:     suffix: 15
1252:     requires: triangle
1253:     args: -run_type test -refinement_limit 0.0    -variable_coefficient analytic -interpolate 1 -petscspace_degree 2 -show_initial -dm_plex_print_fem 1
1254:   test:
1255:     suffix: 16
1256:     requires: triangle
1257:     args: -run_type test -refinement_limit 0.0625 -variable_coefficient analytic -interpolate 1 -petscspace_degree 2 -show_initial -dm_plex_print_fem 1
1258:   test:
1259:     suffix: 17
1260:     requires: ctetgen
1261:     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

1263:   test:
1264:     suffix: 18
1265:     requires: ctetgen
1266:     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

1268:   test:
1269:     suffix: 19
1270:     requires: ctetgen
1271:     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

1273:   test:
1274:     suffix: 20
1275:     requires: ctetgen
1276:     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

1278:   # P1 variable coefficient 21-28
1279:   test:
1280:     suffix: 21
1281:     requires: triangle
1282:     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

1284:   test:
1285:     suffix: 22
1286:     requires: triangle
1287:     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

1289:   test:
1290:     suffix: 23
1291:     requires: triangle
1292:     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

1294:   test:
1295:     suffix: 24
1296:     requires: triangle
1297:     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

1299:   test:
1300:     suffix: 25
1301:     requires: ctetgen
1302:     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

1304:   test:
1305:     suffix: 26
1306:     requires: ctetgen
1307:     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

1309:   test:
1310:     suffix: 27
1311:     requires: ctetgen
1312:     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

1314:   test:
1315:     suffix: 28
1316:     requires: ctetgen
1317:     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

1319:   # P0 variable coefficient 29-36
1320:   test:
1321:     suffix: 29
1322:     requires: triangle
1323:     args: -run_type test -refinement_limit 0.0    -variable_coefficient field    -interpolate 1 -petscspace_degree 1 -show_initial -dm_plex_print_fem 1

1325:   test:
1326:     suffix: 30
1327:     requires: triangle
1328:     args: -run_type test -refinement_limit 0.0625 -variable_coefficient field    -interpolate 1 -petscspace_degree 1 -show_initial -dm_plex_print_fem 1

1330:   test:
1331:     suffix: 31
1332:     requires: triangle
1333:     args: -run_type test -refinement_limit 0.0    -variable_coefficient field    -interpolate 1 -petscspace_degree 2 -show_initial -dm_plex_print_fem 1

1335:   test:
1336:     requires: triangle
1337:     suffix: 32
1338:     args: -run_type test -refinement_limit 0.0625 -variable_coefficient field    -interpolate 1 -petscspace_degree 2 -show_initial -dm_plex_print_fem 1

1340:   test:
1341:     requires: ctetgen
1342:     suffix: 33
1343:     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

1345:   test:
1346:     suffix: 34
1347:     requires: ctetgen
1348:     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

1350:   test:
1351:     suffix: 35
1352:     requires: ctetgen
1353:     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

1355:   test:
1356:     suffix: 36
1357:     requires: ctetgen
1358:     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

1360:   # Full solve 39-44
1361:   test:
1362:     suffix: 39
1363:     requires: triangle !single
1364:     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
1365:   test:
1366:     suffix: 40
1367:     requires: triangle !single
1368:     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
1369:   test:
1370:     suffix: 41
1371:     requires: triangle !single
1372:     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
1373:   test:
1374:     suffix: 42
1375:     requires: triangle !single
1376:     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
1377:   test:
1378:     suffix: 43
1379:     requires: triangle !single
1380:     nsize: 2
1381:     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

1383:   test:
1384:     suffix: 44
1385:     requires: triangle !single
1386:     nsize: 2
1387:     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

1389:   # These tests use a loose tolerance just to exercise the PtAP operations for MATIS and multiple PCBDDC setup calls inside PCMG
1390:   testset:
1391:     requires: triangle !single
1392:     nsize: 3
1393:     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
1394:     test:
1395:       suffix: gmg_bddc
1396:       filter: sed -e "s/CONVERGED_FNORM_RELATIVE iterations 3/CONVERGED_FNORM_RELATIVE iterations 4/g"
1397:       args: -mg_levels_pc_type jacobi
1398:     test:
1399:       filter: sed -e "s/iterations [0-4]/iterations 4/g"
1400:       suffix: gmg_bddc_lev
1401:       args: -mg_levels_pc_type bddc

1403:   # Restarting
1404:   testset:
1405:     suffix: restart
1406:     requires: hdf5 triangle !complex
1407:     args: -run_type test -refinement_limit 0.0    -bc_type dirichlet -interpolate 1 -petscspace_degree 1
1408:     test:
1409:       args: -dm_view hdf5:sol.h5 -vec_view hdf5:sol.h5::append
1410:     test:
1411:       args: -f sol.h5 -restart

1413:   # Periodicity
1414:   test:
1415:     suffix: periodic_0
1416:     requires: triangle
1417:     args: -run_type full -refinement_limit 0.0    -bc_type dirichlet -interpolate 1 -petscspace_degree 1 -snes_converged_reason ::ascii_info_detail

1419:   test:
1420:     requires: !complex
1421:     suffix: periodic_1
1422:     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

1424:   # 2D serial P1 test with field bc
1425:   test:
1426:     suffix: field_bc_2d_p1_0
1427:     requires: triangle
1428:     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

1430:   test:
1431:     suffix: field_bc_2d_p1_1
1432:     requires: triangle
1433:     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

1435:   test:
1436:     suffix: field_bc_2d_p1_neumann_0
1437:     requires: triangle
1438:     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

1440:   test:
1441:     suffix: field_bc_2d_p1_neumann_1
1442:     requires: triangle
1443:     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

1445:   # 3D serial P1 test with field bc
1446:   test:
1447:     suffix: field_bc_3d_p1_0
1448:     requires: ctetgen
1449:     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

1451:   test:
1452:     suffix: field_bc_3d_p1_1
1453:     requires: ctetgen
1454:     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

1456:   test:
1457:     suffix: field_bc_3d_p1_neumann_0
1458:     requires: ctetgen
1459:     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

1461:   test:
1462:     suffix: field_bc_3d_p1_neumann_1
1463:     requires: ctetgen
1464:     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

1466:   # 2D serial P2 test with field bc
1467:   test:
1468:     suffix: field_bc_2d_p2_0
1469:     requires: triangle
1470:     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

1472:   test:
1473:     suffix: field_bc_2d_p2_1
1474:     requires: triangle
1475:     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

1477:   test:
1478:     suffix: field_bc_2d_p2_neumann_0
1479:     requires: triangle
1480:     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

1482:   test:
1483:     suffix: field_bc_2d_p2_neumann_1
1484:     requires: triangle
1485:     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

1487:   # 3D serial P2 test with field bc
1488:   test:
1489:     suffix: field_bc_3d_p2_0
1490:     requires: ctetgen
1491:     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

1493:   test:
1494:     suffix: field_bc_3d_p2_1
1495:     requires: ctetgen
1496:     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

1498:   test:
1499:     suffix: field_bc_3d_p2_neumann_0
1500:     requires: ctetgen
1501:     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

1503:   test:
1504:     suffix: field_bc_3d_p2_neumann_1
1505:     requires: ctetgen
1506:     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

1508:   # Full solve simplex: Convergence
1509:   test:
1510:     suffix: tet_conv_p1_r0
1511:     requires: ctetgen
1512:     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
1513:   test:
1514:     suffix: tet_conv_p1_r2
1515:     requires: ctetgen
1516:     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
1517:   test:
1518:     suffix: tet_conv_p1_r3
1519:     requires: ctetgen
1520:     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
1521:   test:
1522:     suffix: tet_conv_p2_r0
1523:     requires: ctetgen
1524:     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
1525:   test:
1526:     suffix: tet_conv_p2_r2
1527:     requires: ctetgen
1528:     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

1530:   # Full solve simplex: PCBDDC
1531:   test:
1532:     suffix: tri_bddc
1533:     requires: triangle !single
1534:     nsize: 5
1535:     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

1537:   # Full solve simplex: PCBDDC
1538:   test:
1539:     suffix: tri_parmetis_bddc
1540:     requires: triangle !single parmetis
1541:     nsize: 4
1542:     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

1544:   testset:
1545:     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
1546:     nsize: 5
1547:     output_file: output/ex12_quad_bddc.out
1548:     filter: sed -e "s/aijcusparse/aij/g" -e "s/aijviennacl/aij/g" -e "s/factorization: cusparse/factorization: petsc/g"
1549:     test:
1550:       requires: !single
1551:       suffix: quad_bddc
1552:     test:
1553:       requires: !single cuda
1554:       suffix: quad_bddc_cuda
1555:       args: -matis_localmat_type aijcusparse -pc_bddc_dirichlet_pc_factor_mat_solver_type cusparse -pc_bddc_neumann_pc_factor_mat_solver_type cusparse
1556:     test:
1557:       requires: !single viennacl
1558:       suffix: quad_bddc_viennacl
1559:       args: -matis_localmat_type aijviennacl

1561:   # Full solve simplex: ASM
1562:   test:
1563:     suffix: tri_q2q1_asm_lu
1564:     requires: triangle !single
1565:     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

1567:   test:
1568:     suffix: tri_q2q1_msm_lu
1569:     requires: triangle !single
1570:     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

1572:   test:
1573:     suffix: tri_q2q1_asm_sor
1574:     requires: triangle !single
1575:     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

1577:   test:
1578:     suffix: tri_q2q1_msm_sor
1579:     requires: triangle !single
1580:     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

1582:   # Full solve simplex: FAS
1583:   test:
1584:     suffix: fas_newton_0
1585:     requires: triangle !single
1586:     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

1588:   test:
1589:     suffix: fas_newton_1
1590:     requires: triangle !single
1591:     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

1593:   test:
1594:     suffix: fas_ngs_0
1595:     requires: triangle !single
1596:     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

1598:   test:
1599:     suffix: fas_newton_coarse_0
1600:     requires: pragmatic triangle
1601:     TODO: broken
1602:     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

1604:   test:
1605:     suffix: mg_newton_coarse_0
1606:     requires: triangle pragmatic
1607:     TODO: broken
1608:     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

1610:   test:
1611:     suffix: mg_newton_coarse_1
1612:     requires: triangle pragmatic
1613:     TODO: broken
1614:     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

1616:   test:
1617:     suffix: mg_newton_coarse_2
1618:     requires: triangle pragmatic
1619:     TODO: broken
1620:     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

1622:   # Full solve tensor
1623:   test:
1624:     suffix: tensor_plex_2d
1625:     args: -run_type test -refinement_limit 0.0 -simplex 0 -interpolate -bc_type dirichlet -petscspace_degree 1 -dm_refine_hierarchy 2 -cells 2,2

1627:   test:
1628:     suffix: tensor_p4est_2d
1629:     requires: p4est
1630:     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

1632:   test:
1633:     suffix: tensor_plex_3d
1634:     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

1636:   test:
1637:     suffix: tensor_p4est_3d
1638:     requires: p4est
1639:     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

1641:   test:
1642:     suffix: p4est_test_q2_conformal_serial
1643:     requires: p4est
1644:     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

1646:   test:
1647:     suffix: p4est_test_q2_conformal_parallel
1648:     requires: p4est
1649:     nsize: 7
1650:     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

1652:   test:
1653:     suffix: p4est_test_q2_conformal_parallel_parmetis
1654:     requires: parmetis p4est
1655:     nsize: 4
1656:     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

1658:   test:
1659:     suffix: p4est_test_q2_nonconformal_serial
1660:     requires: p4est
1661:     filter: grep -v "CG or CGNE: variant"
1662:     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

1664:   test:
1665:     suffix: p4est_test_q2_nonconformal_parallel
1666:     requires: p4est
1667:     filter: grep -v "CG or CGNE: variant"
1668:     nsize: 7
1669:     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

1671:   test:
1672:     suffix: p4est_test_q2_nonconformal_parallel_parmetis
1673:     requires: parmetis p4est
1674:     nsize: 4
1675:     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

1677:   test:
1678:     suffix: p4est_exact_q2_conformal_serial
1679:     requires: p4est !single !complex !__float128
1680:     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

1682:   test:
1683:     suffix: p4est_exact_q2_conformal_parallel
1684:     requires: p4est !single !complex !__float128
1685:     nsize: 4
1686:     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

1688:   test:
1689:     suffix: p4est_exact_q2_conformal_parallel_parmetis
1690:     requires: parmetis p4est !single
1691:     nsize: 4
1692:     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

1694:   test:
1695:     suffix: p4est_exact_q2_nonconformal_serial
1696:     requires: p4est
1697:     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

1699:   test:
1700:     suffix: p4est_exact_q2_nonconformal_parallel
1701:     requires: p4est
1702:     nsize: 7
1703:     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

1705:   test:
1706:     suffix: p4est_exact_q2_nonconformal_parallel_parmetis
1707:     requires: parmetis p4est
1708:     nsize: 4
1709:     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

1711:   test:
1712:     suffix: p4est_full_q2_nonconformal_serial
1713:     requires: p4est !single
1714:     filter: grep -v "variant HERMITIAN"
1715:     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

1717:   test:
1718:     suffix: p4est_full_q2_nonconformal_parallel
1719:     requires: p4est !single
1720:     filter: grep -v "variant HERMITIAN"
1721:     nsize: 7
1722:     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

1724:   test:
1725:     suffix: p4est_full_q2_nonconformal_parallel_bddcfas
1726:     requires: p4est !single
1727:     filter: grep -v "variant HERMITIAN"
1728:     nsize: 7
1729:     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

1731:   test:
1732:     suffix: p4est_full_q2_nonconformal_parallel_bddc
1733:     requires: p4est !single
1734:     filter: grep -v "variant HERMITIAN"
1735:     nsize: 7
1736:     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

1738:   test:
1739:     TODO: broken
1740:     suffix: p4est_fas_q2_conformal_serial
1741:     requires: p4est !complex !__float128
1742:     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

1744:   test:
1745:     TODO: broken
1746:     suffix: p4est_fas_q2_nonconformal_serial
1747:     requires: p4est
1748:     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

1750:   test:
1751:     suffix: fas_newton_0_p4est
1752:     requires: p4est !single !__float128
1753:     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

1755:   # Full solve simplicial AMR
1756:   test:
1757:     suffix: tri_p1_adapt_0
1758:     requires: pragmatic
1759:     TODO: broken
1760:     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

1762:   test:
1763:     suffix: tri_p1_adapt_1
1764:     requires: pragmatic
1765:     TODO: broken
1766:     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

1768:   test:
1769:     suffix: tri_p1_adapt_analytic_0
1770:     requires: pragmatic
1771:     TODO: broken
1772:     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

1774:   # Full solve tensor AMR
1775:   test:
1776:     suffix: quad_q1_adapt_0
1777:     requires: p4est
1778:     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
1779:     filter: grep -v DM_

1781:   test:
1782:     suffix: amr_0
1783:     nsize: 5
1784:     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

1786:   test:
1787:     suffix: amr_1
1788:     requires: p4est !complex
1789:     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

1791:   test:
1792:     suffix: p4est_solve_bddc
1793:     requires: p4est !complex
1794:     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
1795:     nsize: 4

1797:   test:
1798:     suffix: p4est_solve_fas
1799:     requires: p4est
1800:     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
1801:     nsize: 4
1802:     TODO: identical machine two runs produce slightly different solver trackers

1804:   test:
1805:     suffix: p4est_convergence_test_1
1806:     requires: p4est
1807:     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
1808:     nsize: 4

1810:   test:
1811:     suffix: p4est_convergence_test_2
1812:     requires: p4est
1813:     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

1815:   test:
1816:     suffix: p4est_convergence_test_3
1817:     requires: p4est
1818:     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

1820:   test:
1821:     suffix: p4est_convergence_test_4
1822:     requires: p4est
1823:     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
1824:     timeoutfactor: 5

1826:   # Serial tests with GLVis visualization
1827:   test:
1828:     suffix: glvis_2d_tet_p1
1829:     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
1830:   test:
1831:     suffix: glvis_2d_tet_p2
1832:     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
1833:   test:
1834:     suffix: glvis_2d_hex_p1
1835:     args: -quiet -run_type test -interpolate 1 -bc_type dirichlet -petscspace_degree 1 -vec_view glvis: -simplex 0 -dm_refine 1
1836:   test:
1837:     suffix: glvis_2d_hex_p2
1838:     args: -quiet -run_type test -interpolate 1 -bc_type dirichlet -petscspace_degree 2 -vec_view glvis: -simplex 0 -dm_refine 1
1839:   test:
1840:     suffix: glvis_2d_hex_p2_p4est
1841:     requires: p4est
1842:     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
1843:   test:
1844:     suffix: glvis_2d_tet_p0
1845:     args: -run_type exact  -guess_vec_view glvis: -nonzero_initial_guess 1 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -petscspace_degree 0
1846:   test:
1847:     suffix: glvis_2d_hex_p0
1848:     args: -run_type exact  -guess_vec_view glvis: -nonzero_initial_guess 1 -cells 5,7  -simplex 0 -petscspace_degree 0

1850:   # PCHPDDM tests
1851:   testset:
1852:     nsize: 4
1853:     requires: hpddm !single
1854:     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
1855:     test:
1856:       suffix: quad_singular_hpddm
1857:       args: -cells 6,7
1858:     test:
1859:       requires: p4est
1860:       suffix: p4est_singular_2d_hpddm
1861:       args: -dm_plex_convert_type p4est -dm_forest_minimum_refinement 1 -dm_forest_initial_refinement 3 -dm_forest_maximum_refinement 3
1862:     test:
1863:       requires: p4est
1864:       suffix: p4est_nc_singular_2d_hpddm
1865:       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
1866:   testset:
1867:     nsize: 4
1868:     requires: hpddm triangle !single
1869:     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
1870:     test:
1871:       args: -pc_hpddm_coarse_mat_type baij -options_left no
1872:       suffix: tri_hpddm_reuse_baij
1873:     test:
1874:       requires: !complex
1875:       suffix: tri_hpddm_reuse
1876:   testset:
1877:     nsize: 4
1878:     requires: hpddm !single
1879:     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
1880:     test:
1881:       args: -pc_hpddm_coarse_mat_type baij -options_left no
1882:       suffix: quad_hpddm_reuse_baij
1883:     test:
1884:       requires: !complex
1885:       suffix: quad_hpddm_reuse
1886:   testset:
1887:     nsize: 4
1888:     requires: hpddm !single
1889:     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
1890:     test:
1891:       args: -pc_hpddm_coarse_mat_type baij -options_left no
1892:       suffix: quad_hpddm_reuse_threshold_baij
1893:     test:
1894:       requires: !complex
1895:       suffix: quad_hpddm_reuse_threshold
1896:   testset:
1897:     nsize: 4
1898:     requires: hpddm parmetis !single
1899:     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
1900:     test:
1901:       args: -pc_hpddm_coarse_mat_type baij -options_left no
1902:       suffix: tri_parmetis_hpddm_baij
1903:     test:
1904:       requires: !complex
1905:       suffix: tri_parmetis_hpddm
1906: TEST*/