Actual source code: ex12.c

  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 :~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_BALL, COEFF_CROSS, COEFF_CHECKERBOARD_0, COEFF_CHECKERBOARD_1} CoeffType;

 27: typedef struct {
 28:   RunType        runType;           /* Whether to run tests, or solve the full problem */
 29:   PetscBool      jacobianMF;        /* Whether to calculate the Jacobian action on the fly */
 30:   PetscBool      showInitial, showSolution, restart, quiet, nonzInit;
 31:   /* Problem definition */
 32:   BCType         bcType;
 33:   CoeffType      variableCoefficient;
 34:   PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx);
 35:   PetscBool      fieldBC;
 36:   void           (**exactFields)(PetscInt, PetscInt, PetscInt,
 37:                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
 38:                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
 39:                                  PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);
 40:   PetscBool      bdIntegral;        /* Compute the integral of the solution on the boundary */
 41:   /* Reproducing tests from SISC 40(3), pp. A1473-A1493, 2018 */
 42:   PetscInt       div;               /* Number of divisions */
 43:   PetscInt       k;                 /* Parameter for checkerboard coefficient */
 44:   PetscInt      *kgrid;             /* Random parameter grid */
 45:   PetscBool      rand;              /* Make random assignments */
 46:   /* Solver */
 47:   PC             pcmg;              /* This is needed for error monitoring */
 48:   PetscBool      checkksp;          /* Whether to check the KSPSolve for runType == RUN_TEST */
 49: } AppCtx;

 51: static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
 52: {
 53:   u[0] = 0.0;
 54:   return 0;
 55: }

 57: static PetscErrorCode ecks(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
 58: {
 59:   u[0] = x[0];
 60:   return 0;
 61: }

 63: /*
 64:   In 2D for Dirichlet conditions, we use exact solution:

 66:     u = x^2 + y^2
 67:     f = 4

 69:   so that

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

 73:   For Neumann conditions, we have

 75:     -\nabla u \cdot -\hat y |_{y=0} =  (2y)|_{y=0} =  0 (bottom)
 76:     -\nabla u \cdot  \hat y |_{y=1} = -(2y)|_{y=1} = -2 (top)
 77:     -\nabla u \cdot -\hat x |_{x=0} =  (2x)|_{x=0} =  0 (left)
 78:     -\nabla u \cdot  \hat x |_{x=1} = -(2x)|_{x=1} = -2 (right)

 80:   Which we can express as

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

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

 86:     \int^1_0 x^2 dx + \int^1_0 (1 + y^2) dy + \int^1_0 (x^2 + 1) dx + \int^1_0 y^2 dy = 1/3 + 4/3 + 4/3 + 1/3 = 3 1/3
 87: */
 88: static PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
 89: {
 90:   *u = x[0]*x[0] + x[1]*x[1];
 91:   return 0;
 92: }

 94: static void quadratic_u_field_2d(PetscInt dim, PetscInt Nf, PetscInt NfAux,
 95:                                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
 96:                                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
 97:                                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar uexact[])
 98: {
 99:   uexact[0] = a[0];
100: }

102: static PetscErrorCode ball_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
103: {
104:   const PetscReal alpha   = 500.;
105:   const PetscReal radius2 = PetscSqr(0.15);
106:   const PetscReal r2      = PetscSqr(x[0] - 0.5) + PetscSqr(x[1] - 0.5);
107:   const PetscReal xi      = alpha*(radius2 - r2);

109:   *u = PetscTanhScalar(xi) + 1.0;
110:   return 0;
111: }

113: static PetscErrorCode cross_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
114: {
115:   const PetscReal alpha = 50*4;
116:   const PetscReal xy    = (x[0]-0.5)*(x[1]-0.5);

118:   *u = PetscSinReal(alpha*xy) * (alpha*PetscAbsReal(xy) < 2*PETSC_PI ? 1 : 0.01);
119:   return 0;
120: }

122: static void f0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
123:                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
124:                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
125:                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
126: {
127:   f0[0] = 4.0;
128: }

130: static void f0_ball_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
131:                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
132:                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
133:                       PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
134: {
135:   PetscInt        d;
136:   const PetscReal alpha = 500., radius2 = PetscSqr(0.15);
137:   PetscReal       r2, xi;

139:   for (d = 0, r2 = 0.0; d < dim; ++d) r2 += PetscSqr(x[d] - 0.5);
140:   xi = alpha*(radius2 - r2);
141:   f0[0] = (-2.0*dim*alpha - 8.0*PetscSqr(alpha)*r2*PetscTanhReal(xi)) * PetscSqr(1.0/PetscCoshReal(xi));
142: }

144: static void f0_cross_u_2d(PetscInt dim, PetscInt Nf, PetscInt NfAux,
145:                           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
146:                           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
147:                           PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
148: {
149:   const PetscReal alpha = 50*4;
150:   const PetscReal xy    = (x[0]-0.5)*(x[1]-0.5);

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

155: static void f0_checkerboard_0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
156:                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
157:                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
158:                                 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
159: {
160:   f0[0] = -20.0*PetscExpReal(-(PetscSqr(x[0] - 0.5) + PetscSqr(x[1] - 0.5)));
161: }

163: static void f0_bd_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
164:                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
165:                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
166:                     PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
167: {
168:   PetscInt d;
169:   for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += -n[d]*2.0*x[d];
170: }

172: /* gradU[comp*dim+d] = {u_x, u_y} or {u_x, u_y, u_z} */
173: static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
174:                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
175:                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
176:                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
177: {
178:   PetscInt d;
179:   for (d = 0; d < dim; ++d) f1[d] = u_x[d];
180: }

182: /* < \nabla v, \nabla u + {\nabla u}^T >
183:    This just gives \nabla u, give the perdiagonal for the transpose */
184: static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
185:                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
186:                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
187:                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
188: {
189:   PetscInt d;
190:   for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0;
191: }

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

196:     u = sin(2 pi x)
197:     f = -4 pi^2 sin(2 pi x)

199:   so that

201:     -\Delta u + f = 4 pi^2 sin(2 pi x) - 4 pi^2 sin(2 pi x) = 0
202: */
203: static PetscErrorCode xtrig_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
204: {
205:   *u = PetscSinReal(2.0*PETSC_PI*x[0]);
206:   return 0;
207: }

209: static void f0_xtrig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
210:                        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
211:                        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
212:                        PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
213: {
214:   f0[0] = -4.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[0]);
215: }

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

220:     u = sin(2 pi x) sin(2 pi y)
221:     f = -8 pi^2 sin(2 pi x)

223:   so that

225:     -\Delta u + f = 4 pi^2 sin(2 pi x) sin(2 pi y) + 4 pi^2 sin(2 pi x) sin(2 pi y) - 8 pi^2 sin(2 pi x) = 0
226: */
227: static PetscErrorCode xytrig_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
228: {
229:   *u = PetscSinReal(2.0*PETSC_PI*x[0])*PetscSinReal(2.0*PETSC_PI*x[1]);
230:   return 0;
231: }

233: static void f0_xytrig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
234:                         const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
235:                         const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
236:                         PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
237: {
238:   f0[0] = -8.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[0]);
239: }

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

244:     u  = x^2 + y^2
245:     f  = 6 (x + y)
246:     nu = (x + y)

248:   so that

250:     -\div \nu \grad u + f = -6 (x + y) + 6 (x + y) = 0
251: */
252: static PetscErrorCode nu_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
253: {
254:   *u = x[0] + x[1];
255:   return 0;
256: }

258: static PetscErrorCode checkerboardCoeff(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
259: {
260:   AppCtx  *user = (AppCtx *) ctx;
261:   PetscInt div  = user->div;
262:   PetscInt k    = user->k;
263:   PetscInt mask = 0, ind = 0, d;

266:   for (d = 0; d < dim; ++d) mask = (mask + (PetscInt) (x[d]*div)) % 2;
267:   if (user->kgrid) {
268:     for (d = 0; d < dim; ++d) {
269:       if (d > 0) ind *= dim;
270:       ind += (PetscInt) (x[d]*div);
271:     }
272:     k = user->kgrid[ind];
273:   }
274:   u[0] = mask ? 1.0 : PetscPowRealInt(10.0, -k);
275:   return 0;
276: }

278: void f0_analytic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
279:                    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
280:                    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
281:                    PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
282: {
283:   f0[0] = 6.0*(x[0] + x[1]);
284: }

286: /* gradU[comp*dim+d] = {u_x, u_y} or {u_x, u_y, u_z} */
287: void f1_analytic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
288:                    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
289:                    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
290:                    PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
291: {
292:   PetscInt d;
293:   for (d = 0; d < dim; ++d) f1[d] = (x[0] + x[1])*u_x[d];
294: }

296: void f1_field_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
297:                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
298:                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
299:                 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
300: {
301:   PetscInt d;
302:   for (d = 0; d < dim; ++d) f1[d] = a[0]*u_x[d];
303: }

305: /* < \nabla v, \nabla u + {\nabla u}^T >
306:    This just gives \nabla u, give the perdiagonal for the transpose */
307: void g3_analytic_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
308:                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
309:                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
310:                     PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
311: {
312:   PetscInt d;
313:   for (d = 0; d < dim; ++d) g3[d*dim+d] = x[0] + x[1];
314: }

316: void g3_field_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
317:                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
318:                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
319:                  PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
320: {
321:   PetscInt d;
322:   for (d = 0; d < dim; ++d) g3[d*dim+d] = a[0];
323: }

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

328:     u  = x^2 + y^2
329:     f  = 16 (x^2 + y^2)
330:     nu = 1/2 |grad u|^2

332:   so that

334:     -\div \nu \grad u + f = -16 (x^2 + y^2) + 16 (x^2 + y^2) = 0
335: */
336: void f0_analytic_nonlinear_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
337:                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
338:                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
339:                              PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
340: {
341:   f0[0] = 16.0*(x[0]*x[0] + x[1]*x[1]);
342: }

344: /* gradU[comp*dim+d] = {u_x, u_y} or {u_x, u_y, u_z} */
345: void f1_analytic_nonlinear_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
346:                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
347:                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
348:                              PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
349: {
350:   PetscScalar nu = 0.0;
351:   PetscInt    d;
352:   for (d = 0; d < dim; ++d) nu += u_x[d]*u_x[d];
353:   for (d = 0; d < dim; ++d) f1[d] = 0.5*nu*u_x[d];
354: }

356: /*
357:   grad (u + eps w) - grad u = eps grad w

359:   1/2 |grad (u + eps w)|^2 grad (u + eps w) - 1/2 |grad u|^2 grad u
360: = 1/2 (|grad u|^2 + 2 eps <grad u,grad w>) (grad u + eps grad w) - 1/2 |grad u|^2 grad u
361: = 1/2 (eps |grad u|^2 grad w + 2 eps <grad u,grad w> grad u)
362: = eps (1/2 |grad u|^2 grad w + grad u <grad u,grad w>)
363: */
364: void g3_analytic_nonlinear_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
365:                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
366:                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
367:                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
368: {
369:   PetscScalar nu = 0.0;
370:   PetscInt    d, e;
371:   for (d = 0; d < dim; ++d) nu += u_x[d]*u_x[d];
372:   for (d = 0; d < dim; ++d) {
373:     g3[d*dim+d] = 0.5*nu;
374:     for (e = 0; e < dim; ++e) {
375:       g3[d*dim+e] += u_x[d]*u_x[e];
376:     }
377:   }
378: }

380: /*
381:   In 3D for Dirichlet conditions we use exact solution:

383:     u = 2/3 (x^2 + y^2 + z^2)
384:     f = 4

386:   so that

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

390:   For Neumann conditions, we have

392:     -\nabla u \cdot -\hat z |_{z=0} =  (2z)|_{z=0} =  0 (bottom)
393:     -\nabla u \cdot  \hat z |_{z=1} = -(2z)|_{z=1} = -2 (top)
394:     -\nabla u \cdot -\hat y |_{y=0} =  (2y)|_{y=0} =  0 (front)
395:     -\nabla u \cdot  \hat y |_{y=1} = -(2y)|_{y=1} = -2 (back)
396:     -\nabla u \cdot -\hat x |_{x=0} =  (2x)|_{x=0} =  0 (left)
397:     -\nabla u \cdot  \hat x |_{x=1} = -(2x)|_{x=1} = -2 (right)

399:   Which we can express as

401:     \nabla u \cdot  \hat n|_\Gamma = {2 x, 2 y, 2z} \cdot \hat n = 2 (x + y + z)
402: */
403: static PetscErrorCode quadratic_u_3d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
404: {
405:   *u = 2.0*(x[0]*x[0] + x[1]*x[1] + x[2]*x[2])/3.0;
406:   return 0;
407: }

409: static PetscErrorCode ball_u_3d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
410: {
411:   const PetscReal alpha   = 500.;
412:   const PetscReal radius2 = PetscSqr(0.15);
413:   const PetscReal r2      = PetscSqr(x[0] - 0.5) + PetscSqr(x[1] - 0.5) + PetscSqr(x[2] - 0.5);
414:   const PetscReal xi      = alpha*(radius2 - r2);

416:   *u = PetscTanhScalar(xi) + 1.0;
417:   return 0;
418: }

420: static void quadratic_u_field_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux,
421:                                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
422:                                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
423:                                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar uexact[])
424: {
425:   uexact[0] = a[0];
426: }

428: static PetscErrorCode cross_u_3d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
429: {
430:   const PetscReal alpha = 50*4;
431:   const PetscReal xyz   = (x[0]-0.5)*(x[1]-0.5)*(x[2]-0.5);

433:   *u = PetscSinReal(alpha*xyz) * (alpha*PetscAbsReal(xyz) < 2*PETSC_PI ? (alpha*PetscAbsReal(xyz) > -2*PETSC_PI ? 1.0 : 0.01) : 0.01);
434:   return 0;
435: }

437: static void f0_cross_u_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux,
438:                           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
439:                           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
440:                           PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
441: {
442:   const PetscReal alpha = 50*4;
443:   const PetscReal xyz   = (x[0]-0.5)*(x[1]-0.5)*(x[2]-0.5);

445:   f0[0] = PetscSinReal(alpha*xyz) * (alpha*PetscAbsReal(xyz) < 2*PETSC_PI ? (alpha*PetscAbsReal(xyz) > -2*PETSC_PI ? 1.0 : 0.01) : 0.01);
446: }

448: static void bd_integral_2d(PetscInt dim, PetscInt Nf, PetscInt NfAux,
449:                            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
450:                            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
451:                            PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *uint)
452: {
453:   uint[0] = u[0];
454: }

456: static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
457: {
458:   const char    *bcTypes[3]  = {"neumann", "dirichlet", "none"};
459:   const char    *runTypes[4] = {"full", "exact", "test", "perf"};
460:   const char    *coeffTypes[8] = {"none", "analytic", "field", "nonlinear", "ball", "cross", "checkerboard_0", "checkerboard_1"};
461:   PetscInt       bc, run, coeff;

465:   options->runType             = RUN_FULL;
466:   options->bcType              = DIRICHLET;
467:   options->variableCoefficient = COEFF_NONE;
468:   options->fieldBC             = PETSC_FALSE;
469:   options->jacobianMF          = PETSC_FALSE;
470:   options->showInitial         = PETSC_FALSE;
471:   options->showSolution        = PETSC_FALSE;
472:   options->restart             = PETSC_FALSE;
473:   options->quiet               = PETSC_FALSE;
474:   options->nonzInit            = PETSC_FALSE;
475:   options->bdIntegral          = PETSC_FALSE;
476:   options->checkksp            = PETSC_FALSE;
477:   options->div                 = 4;
478:   options->k                   = 1;
479:   options->kgrid               = NULL;
480:   options->rand                = PETSC_FALSE;

482:   PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX");
483:   run  = options->runType;
484:   PetscOptionsEList("-run_type", "The run type", "ex12.c", runTypes, 4, runTypes[options->runType], &run, NULL);
485:   options->runType = (RunType) run;
486:   bc   = options->bcType;
487:   PetscOptionsEList("-bc_type","Type of boundary condition","ex12.c",bcTypes,3,bcTypes[options->bcType],&bc,NULL);
488:   options->bcType = (BCType) bc;
489:   coeff = options->variableCoefficient;
490:   PetscOptionsEList("-variable_coefficient","Type of variable coefficent","ex12.c",coeffTypes,8,coeffTypes[options->variableCoefficient],&coeff,NULL);
491:   options->variableCoefficient = (CoeffType) coeff;

493:   PetscOptionsBool("-field_bc", "Use a field representation for the BC", "ex12.c", options->fieldBC, &options->fieldBC, NULL);
494:   PetscOptionsBool("-jacobian_mf", "Calculate the action of the Jacobian on the fly", "ex12.c", options->jacobianMF, &options->jacobianMF, NULL);
495:   PetscOptionsBool("-show_initial", "Output the initial guess for verification", "ex12.c", options->showInitial, &options->showInitial, NULL);
496:   PetscOptionsBool("-show_solution", "Output the solution for verification", "ex12.c", options->showSolution, &options->showSolution, NULL);
497:   PetscOptionsBool("-restart", "Read in the mesh and solution from a file", "ex12.c", options->restart, &options->restart, NULL);
498:   PetscOptionsBool("-quiet", "Don't print any vecs", "ex12.c", options->quiet, &options->quiet, NULL);
499:   PetscOptionsBool("-nonzero_initial_guess", "nonzero initial guess", "ex12.c", options->nonzInit, &options->nonzInit, NULL);
500:   PetscOptionsBool("-bd_integral", "Compute the integral of the solution on the boundary", "ex12.c", options->bdIntegral, &options->bdIntegral, NULL);
501:   if (options->runType == RUN_TEST) {
502:     PetscOptionsBool("-run_test_check_ksp", "Check solution of KSP", "ex12.c", options->checkksp, &options->checkksp, NULL);
503:   }
504:   PetscOptionsInt("-div", "The number of division for the checkerboard coefficient", "ex12.c", options->div, &options->div, NULL);
505:   PetscOptionsInt("-k", "The exponent for the checkerboard coefficient", "ex12.c", options->k, &options->k, NULL);
506:   PetscOptionsBool("-k_random", "Assign random k values to checkerboard", "ex12.c", options->rand, &options->rand, NULL);
507:   PetscOptionsEnd();
508:   return 0;
509: }

511: static PetscErrorCode CreateBCLabel(DM dm, const char name[])
512: {
513:   DM             plex;
514:   DMLabel        label;

517:   DMCreateLabel(dm, name);
518:   DMGetLabel(dm, name, &label);
519:   DMConvert(dm, DMPLEX, &plex);
520:   DMPlexMarkBoundaryFaces(plex, 1, label);
521:   DMDestroy(&plex);
522:   return 0;
523: }

525: static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
526: {

530:   DMCreate(comm, dm);
531:   DMSetType(*dm, DMPLEX);
532:   DMSetFromOptions(*dm);
533:   {
534:     char      convType[256];
535:     PetscBool flg;

537:     PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX");
538:     PetscOptionsFList("-dm_plex_convert_type","Convert DMPlex to another format","ex12",DMList,DMPLEX,convType,256,&flg);
539:     PetscOptionsEnd();
540:     if (flg) {
541:       DM dmConv;

543:       DMConvert(*dm,convType,&dmConv);
544:       if (dmConv) {
545:         DMDestroy(dm);
546:         *dm  = dmConv;
547:       }
548:       DMSetFromOptions(*dm);
549:       DMSetUp(*dm);
550:     }
551:   }
552:   DMViewFromOptions(*dm, NULL, "-dm_view");
553:   if (user->rand) {
554:     PetscRandom r;
555:     PetscReal   val;
556:     PetscInt    dim, N, i;

558:     DMGetDimension(*dm, &dim);
559:     N    = PetscPowInt(user->div, dim);
560:     PetscMalloc1(N, &user->kgrid);
561:     PetscRandomCreate(PETSC_COMM_SELF, &r);
562:     PetscRandomSetFromOptions(r);
563:     PetscRandomSetInterval(r, 0.0, user->k);
564:     PetscRandomSetSeed(r, 1973);
565:     PetscRandomSeed(r);
566:     for (i = 0; i < N; ++i) {
567:       PetscRandomGetValueReal(r, &val);
568:       user->kgrid[i] = 1 + (PetscInt) val;
569:     }
570:     PetscRandomDestroy(&r);
571:   }
572:   return 0;
573: }

575: static PetscErrorCode SetupProblem(DM dm, AppCtx *user)
576: {
577:   PetscDS         ds;
578:   DMLabel         label;
579:   PetscWeakForm   wf;
580:   const DMBoundaryType *periodicity;
581:   const PetscInt  id = 1;
582:   PetscInt        bd, dim;

585:   DMGetDS(dm, &ds);
586:   DMGetDimension(dm, &dim);
587:   DMGetPeriodicity(dm, NULL, NULL, NULL, &periodicity);
588:   switch (user->variableCoefficient) {
589:   case COEFF_NONE:
590:     if (periodicity && periodicity[0]) {
591:       if (periodicity && periodicity[1]) {
592:         PetscDSSetResidual(ds, 0, f0_xytrig_u, f1_u);
593:         PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu);
594:       } else {
595:         PetscDSSetResidual(ds, 0, f0_xtrig_u,  f1_u);
596:         PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu);
597:       }
598:     } else {
599:       PetscDSSetResidual(ds, 0, f0_u, f1_u);
600:       PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu);
601:     }
602:     break;
603:   case COEFF_ANALYTIC:
604:     PetscDSSetResidual(ds, 0, f0_analytic_u, f1_analytic_u);
605:     PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_analytic_uu);
606:     break;
607:   case COEFF_FIELD:
608:     PetscDSSetResidual(ds, 0, f0_analytic_u, f1_field_u);
609:     PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_field_uu);
610:     break;
611:   case COEFF_NONLINEAR:
612:     PetscDSSetResidual(ds, 0, f0_analytic_nonlinear_u, f1_analytic_nonlinear_u);
613:     PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_analytic_nonlinear_uu);
614:     break;
615:   case COEFF_BALL:
616:     PetscDSSetResidual(ds, 0, f0_ball_u, f1_u);
617:     PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu);
618:     break;
619:   case COEFF_CROSS:
620:     switch (dim) {
621:     case 2:
622:       PetscDSSetResidual(ds, 0, f0_cross_u_2d, f1_u);
623:       break;
624:     case 3:
625:       PetscDSSetResidual(ds, 0, f0_cross_u_3d, f1_u);
626:       break;
627:     default:
628:       SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %d", dim);
629:     }
630:     PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu);
631:     break;
632:   case COEFF_CHECKERBOARD_0:
633:     PetscDSSetResidual(ds, 0, f0_checkerboard_0_u, f1_field_u);
634:     PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_field_uu);
635:     break;
636:   default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid variable coefficient type %d", user->variableCoefficient);
637:   }
638:   switch (dim) {
639:   case 2:
640:     switch (user->variableCoefficient) {
641:     case COEFF_BALL:
642:       user->exactFuncs[0]  = ball_u_2d;break;
643:     case COEFF_CROSS:
644:       user->exactFuncs[0]  = cross_u_2d;break;
645:     case COEFF_CHECKERBOARD_0:
646:       user->exactFuncs[0]  = zero;break;
647:     default:
648:       if (periodicity && periodicity[0]) {
649:         if (periodicity && periodicity[1]) {
650:           user->exactFuncs[0] = xytrig_u_2d;
651:         } else {
652:           user->exactFuncs[0] = xtrig_u_2d;
653:         }
654:       } else {
655:         user->exactFuncs[0]  = quadratic_u_2d;
656:         user->exactFields[0] = quadratic_u_field_2d;
657:       }
658:     }
659:     if (user->bcType == NEUMANN) {
660:       DMGetLabel(dm, "boundary", &label);
661:       DMAddBoundary(dm, DM_BC_NATURAL, "wall", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd);
662:       PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL);
663:       PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_bd_u, 0, NULL);
664:     }
665:     break;
666:   case 3:
667:     switch (user->variableCoefficient) {
668:     case COEFF_BALL:
669:       user->exactFuncs[0]  = ball_u_3d;break;
670:     case COEFF_CROSS:
671:       user->exactFuncs[0]  = cross_u_3d;break;
672:     default:
673:       user->exactFuncs[0]  = quadratic_u_3d;
674:       user->exactFields[0] = quadratic_u_field_3d;
675:     }
676:     if (user->bcType == NEUMANN) {
677:       DMGetLabel(dm, "boundary", &label);
678:       DMAddBoundary(dm, DM_BC_NATURAL, "wall", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd);
679:       PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL);
680:       PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_bd_u, 0, NULL);
681:     }
682:     break;
683:   default:
684:     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %d", dim);
685:   }
686:   /* Setup constants */
687:   switch (user->variableCoefficient) {
688:   case COEFF_CHECKERBOARD_0:
689:   {
690:     PetscScalar constants[2];

692:     constants[0] = user->div;
693:     constants[1] = user->k;
694:     PetscDSSetConstants(ds, 2, constants);
695:   }
696:   break;
697:   default: break;
698:   }
699:   PetscDSSetExactSolution(ds, 0, user->exactFuncs[0], user);
700:   /* Setup Boundary Conditions */
701:   if (user->bcType == DIRICHLET) {
702:     DMGetLabel(dm, "marker", &label);
703:     if (!label) {
704:       /* Right now, p4est cannot create labels immediately */
705:       PetscDSAddBoundaryByName(ds, user->fieldBC ? DM_BC_ESSENTIAL_FIELD : DM_BC_ESSENTIAL, "wall", "marker", 1, &id, 0, 0, NULL, user->fieldBC ? (void (*)(void)) user->exactFields[0] : (void (*)(void)) user->exactFuncs[0], NULL, user, NULL);
706:     } else {
707:       DMAddBoundary(dm, user->fieldBC ? DM_BC_ESSENTIAL_FIELD : DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, user->fieldBC ? (void (*)(void)) user->exactFields[0] : (void (*)(void)) user->exactFuncs[0], NULL, user, NULL);
708:     }
709:   }
710:   return 0;
711: }

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

719:   ctx[0] = user;
720:   if (user->variableCoefficient == COEFF_CHECKERBOARD_0) {matFuncs[0] = checkerboardCoeff;}
721:   DMCreateLocalVector(dmAux, &nu);
722:   PetscObjectSetName((PetscObject) nu, "Coefficient");
723:   DMProjectFunctionLocal(dmAux, 0.0, matFuncs, ctx, INSERT_ALL_VALUES, nu);
724:   DMSetAuxiliaryVec(dm, NULL, 0, 0, nu);
725:   VecDestroy(&nu);
726:   return 0;
727: }

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

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

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

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

762: static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
763: {
764:   DM             plex, cdm = dm;
765:   PetscFE        fe, feAux = NULL;
766:   PetscBool      simplex;
767:   PetscInt       dim;
768:   MPI_Comm       comm;

771:   DMGetDimension(dm, &dim);
772:   DMConvert(dm, DMPLEX, &plex);
773:   DMPlexIsSimplex(plex, &simplex);
774:   DMDestroy(&plex);
775:   PetscObjectGetComm((PetscObject) dm, &comm);
776:   PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, NULL, -1, &fe);
777:   PetscObjectSetName((PetscObject) fe, "potential");
778:   if (user->variableCoefficient == COEFF_FIELD || user->variableCoefficient == COEFF_CHECKERBOARD_0) {
779:     PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "mat_", -1, &feAux);
780:     PetscObjectSetName((PetscObject) feAux, "coefficient");
781:     PetscFECopyQuadrature(fe, feAux);
782:   } else if (user->fieldBC) {
783:     PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "bc_", -1, &feAux);
784:     PetscFECopyQuadrature(fe, feAux);
785:   }
786:   /* Set discretization and boundary conditions for each mesh */
787:   DMSetField(dm, 0, NULL, (PetscObject) fe);
788:   DMCreateDS(dm);
789:   SetupProblem(dm, user);
790:   while (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:     DMCopyDisc(dm, cdm);
799:     DMGetCoarseDM(cdm, &cdm);
800:   }
801:   PetscFEDestroy(&fe);
802:   PetscFEDestroy(&feAux);
803:   return 0;
804: }

806: int main(int argc, char **argv)
807: {
808:   DM             dm;          /* Problem specification */
809:   SNES           snes;        /* nonlinear solver */
810:   Vec            u;           /* solution vector */
811:   Mat            A,J;         /* Jacobian matrix */
812:   MatNullSpace   nullSpace;   /* May be necessary for Neumann conditions */
813:   AppCtx         user;        /* user-defined work context */
814:   JacActionCtx   userJ;       /* context for Jacobian MF action */
815:   PetscReal      error = 0.0; /* L_2 error in the solution */

817:   PetscInitialize(&argc, &argv, NULL,help);
818:   ProcessOptions(PETSC_COMM_WORLD, &user);
819:   SNESCreate(PETSC_COMM_WORLD, &snes);
820:   CreateMesh(PETSC_COMM_WORLD, &user, &dm);
821:   SNESSetDM(snes, dm);
822:   DMSetApplicationContext(dm, &user);

824:   PetscMalloc2(1, &user.exactFuncs, 1, &user.exactFields);
825:   SetupDiscretization(dm, &user);

827:   DMCreateGlobalVector(dm, &u);
828:   PetscObjectSetName((PetscObject) u, "potential");

830:   DMCreateMatrix(dm, &J);
831:   if (user.jacobianMF) {
832:     PetscInt M, m, N, n;

834:     MatGetSize(J, &M, &N);
835:     MatGetLocalSize(J, &m, &n);
836:     MatCreate(PETSC_COMM_WORLD, &A);
837:     MatSetSizes(A, m, n, M, N);
838:     MatSetType(A, MATSHELL);
839:     MatSetUp(A);
840: #if 0
841:     MatShellSetOperation(A, MATOP_MULT, (void (*)(void))FormJacobianAction);
842: #endif

844:     userJ.dm   = dm;
845:     userJ.J    = J;
846:     userJ.user = &user;

848:     DMCreateLocalVector(dm, &userJ.u);
849:     if (user.fieldBC) DMProjectFieldLocal(dm, 0.0, userJ.u, user.exactFields, INSERT_BC_VALUES, userJ.u);
850:     else              DMProjectFunctionLocal(dm, 0.0, user.exactFuncs, NULL, INSERT_BC_VALUES, userJ.u);
851:     MatShellSetContext(A, &userJ);
852:   } else {
853:     A = J;
854:   }

856:   nullSpace = NULL;
857:   if (user.bcType != DIRICHLET) {
858:     MatNullSpaceCreate(PetscObjectComm((PetscObject) dm), PETSC_TRUE, 0, NULL, &nullSpace);
859:     MatSetNullSpace(A, nullSpace);
860:   }

862:   DMPlexSetSNESLocalFEM(dm,&user,&user,&user);
863:   SNESSetJacobian(snes, A, J, NULL, NULL);

865:   SNESSetFromOptions(snes);

867:   if (user.fieldBC) DMProjectField(dm, 0.0, u, user.exactFields, INSERT_ALL_VALUES, u);
868:   else              DMProjectFunction(dm, 0.0, user.exactFuncs, NULL, INSERT_ALL_VALUES, u);
869:   if (user.restart) {
870: #if defined(PETSC_HAVE_HDF5)
871:     PetscViewer viewer;
872:     char        filename[PETSC_MAX_PATH_LEN];

874:     PetscOptionsGetString(NULL, NULL, "-dm_plex_filename", filename, sizeof(filename), NULL);
875:     PetscViewerCreate(PETSC_COMM_WORLD, &viewer);
876:     PetscViewerSetType(viewer, PETSCVIEWERHDF5);
877:     PetscViewerFileSetMode(viewer, FILE_MODE_READ);
878:     PetscViewerFileSetName(viewer, filename);
879:     PetscViewerHDF5PushGroup(viewer, "/fields");
880:     VecLoad(u, viewer);
881:     PetscViewerHDF5PopGroup(viewer);
882:     PetscViewerDestroy(&viewer);
883: #endif
884:   }
885:   if (user.showInitial) {
886:     Vec lv;
887:     DMGetLocalVector(dm, &lv);
888:     DMGlobalToLocalBegin(dm, u, INSERT_VALUES, lv);
889:     DMGlobalToLocalEnd(dm, u, INSERT_VALUES, lv);
890:     DMPrintLocalVec(dm, "Local function", 1.0e-10, lv);
891:     DMRestoreLocalVector(dm, &lv);
892:   }
893:   if (user.runType == RUN_FULL || user.runType == RUN_EXACT) {
894:     PetscErrorCode (*initialGuess[1])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx) = {zero};

896:     if (user.nonzInit) initialGuess[0] = ecks;
897:     if (user.runType == RUN_FULL) {
898:       DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u);
899:     }
900:     VecViewFromOptions(u, NULL, "-guess_vec_view");
901:     SNESSolve(snes, NULL, u);
902:     SNESGetSolution(snes, &u);
903:     SNESGetDM(snes, &dm);

905:     if (user.showSolution) {
906:       PetscPrintf(PETSC_COMM_WORLD, "Solution\n");
907:       VecChop(u, 3.0e-9);
908:       VecView(u, PETSC_VIEWER_STDOUT_WORLD);
909:     }
910:   } else if (user.runType == RUN_PERF) {
911:     Vec       r;
912:     PetscReal res = 0.0;

914:     SNESGetFunction(snes, &r, NULL, NULL);
915:     SNESComputeFunction(snes, u, r);
916:     PetscPrintf(PETSC_COMM_WORLD, "Initial Residual\n");
917:     VecChop(r, 1.0e-10);
918:     VecNorm(r, NORM_2, &res);
919:     PetscPrintf(PETSC_COMM_WORLD, "L_2 Residual: %g\n", (double)res);
920:   } else {
921:     Vec       r;
922:     PetscReal res = 0.0, tol = 1.0e-11;

924:     /* Check discretization error */
925:     SNESGetFunction(snes, &r, NULL, NULL);
926:     PetscPrintf(PETSC_COMM_WORLD, "Initial guess\n");
927:     if (!user.quiet) VecView(u, PETSC_VIEWER_STDOUT_WORLD);
928:     DMComputeL2Diff(dm, 0.0, user.exactFuncs, NULL, u, &error);
929:     if (error < tol) PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: < %2.1e\n", (double)tol);
930:     else             PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: %g\n", (double)error);
931:     /* Check residual */
932:     SNESComputeFunction(snes, u, r);
933:     PetscPrintf(PETSC_COMM_WORLD, "Initial Residual\n");
934:     VecChop(r, 1.0e-10);
935:     if (!user.quiet) VecView(r, PETSC_VIEWER_STDOUT_WORLD);
936:     VecNorm(r, NORM_2, &res);
937:     PetscPrintf(PETSC_COMM_WORLD, "L_2 Residual: %g\n", (double)res);
938:     /* Check Jacobian */
939:     {
940:       Vec b;

942:       SNESComputeJacobian(snes, u, A, A);
943:       VecDuplicate(u, &b);
944:       VecSet(r, 0.0);
945:       SNESComputeFunction(snes, r, b);
946:       MatMult(A, u, r);
947:       VecAXPY(r, 1.0, b);
948:       PetscPrintf(PETSC_COMM_WORLD, "Au - b = Au + F(0)\n");
949:       VecChop(r, 1.0e-10);
950:       if (!user.quiet) VecView(r, PETSC_VIEWER_STDOUT_WORLD);
951:       VecNorm(r, NORM_2, &res);
952:       PetscPrintf(PETSC_COMM_WORLD, "Linear L_2 Residual: %g\n", (double)res);
953:       /* check solver */
954:       if (user.checkksp) {
955:         KSP ksp;

957:         if (nullSpace) {
958:           MatNullSpaceRemove(nullSpace, u);
959:         }
960:         SNESComputeJacobian(snes, u, A, J);
961:         MatMult(A, u, b);
962:         SNESGetKSP(snes, &ksp);
963:         KSPSetOperators(ksp, A, J);
964:         KSPSolve(ksp, b, r);
965:         VecAXPY(r, -1.0, u);
966:         VecNorm(r, NORM_2, &res);
967:         PetscPrintf(PETSC_COMM_WORLD, "KSP Error: %g\n", (double)res);
968:       }
969:       VecDestroy(&b);
970:     }
971:   }
972:   VecViewFromOptions(u, NULL, "-vec_view");
973:   {
974:     Vec nu;

976:     DMGetAuxiliaryVec(dm, NULL, 0, 0, &nu);
977:     if (nu) VecViewFromOptions(nu, NULL, "-coeff_view");
978:   }

980:   if (user.bdIntegral) {
981:     DMLabel   label;
982:     PetscInt  id = 1;
983:     PetscScalar bdInt = 0.0;
984:     PetscReal   exact = 3.3333333333;

986:     DMGetLabel(dm, "marker", &label);
987:     DMPlexComputeBdIntegral(dm, u, label, 1, &id, bd_integral_2d, &bdInt, NULL);
988:     PetscPrintf(PETSC_COMM_WORLD, "Solution boundary integral: %.4g\n", (double) PetscAbsScalar(bdInt));
990:   }

992:   MatNullSpaceDestroy(&nullSpace);
993:   if (user.jacobianMF) VecDestroy(&userJ.u);
994:   if (A != J) MatDestroy(&A);
995:   MatDestroy(&J);
996:   VecDestroy(&u);
997:   SNESDestroy(&snes);
998:   DMDestroy(&dm);
999:   PetscFree2(user.exactFuncs, user.exactFields);
1000:   PetscFree(user.kgrid);
1001:   PetscFinalize();
1002:   return 0;
1003: }

1005: /*TEST
1006:   # 2D serial P1 test 0-4
1007:   test:
1008:     suffix: 2d_p1_0
1009:     requires: triangle
1010:     args: -run_type test -bc_type dirichlet -dm_plex_interpolate 0 -petscspace_degree 1 -show_initial -dm_plex_print_fem 1

1012:   test:
1013:     suffix: 2d_p1_1
1014:     requires: triangle
1015:     args: -run_type test -bc_type dirichlet -petscspace_degree 1 -show_initial -dm_plex_print_fem 1

1017:   test:
1018:     suffix: 2d_p1_2
1019:     requires: triangle
1020:     args: -run_type test -dm_refine_volume_limit_pre 0.0625 -bc_type dirichlet -petscspace_degree 1 -show_initial -dm_plex_print_fem 1

1022:   test:
1023:     suffix: 2d_p1_neumann_0
1024:     requires: triangle
1025:     args: -dm_coord_space 0 -run_type test -bc_type neumann -dm_plex_boundary_label boundary -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 -dm_view ascii::ascii_info_detail

1027:   test:
1028:     suffix: 2d_p1_neumann_1
1029:     requires: triangle
1030:     args: -run_type test -dm_refine_volume_limit_pre 0.0625 -bc_type neumann -dm_plex_boundary_label boundary -petscspace_degree 1 -show_initial -dm_plex_print_fem 1

1032:   # 2D serial P2 test 5-8
1033:   test:
1034:     suffix: 2d_p2_0
1035:     requires: triangle
1036:     args: -run_type test -bc_type dirichlet -petscspace_degree 2 -show_initial -dm_plex_print_fem 1

1038:   test:
1039:     suffix: 2d_p2_1
1040:     requires: triangle
1041:     args: -run_type test -dm_refine_volume_limit_pre 0.0625 -bc_type dirichlet -petscspace_degree 2 -show_initial -dm_plex_print_fem 1

1043:   test:
1044:     suffix: 2d_p2_neumann_0
1045:     requires: triangle
1046:     args: -dm_coord_space 0 -run_type test -bc_type neumann -dm_plex_boundary_label boundary -petscspace_degree 2 -show_initial -dm_plex_print_fem 1 -dm_view ascii::ascii_info_detail

1048:   test:
1049:     suffix: 2d_p2_neumann_1
1050:     requires: triangle
1051:     args: -dm_coord_space 0 -run_type test -dm_refine_volume_limit_pre 0.0625 -bc_type neumann -dm_plex_boundary_label boundary -petscspace_degree 2 -show_initial -dm_plex_print_fem 1 -dm_view ascii::ascii_info_detail

1053:   test:
1054:     suffix: bd_int_0
1055:     requires: triangle
1056:     args: -run_type test -bc_type dirichlet -petscspace_degree 2 -bd_integral -dm_view -quiet

1058:   test:
1059:     suffix: bd_int_1
1060:     requires: triangle
1061:     args: -run_type test -dm_refine 2 -bc_type dirichlet -petscspace_degree 2 -bd_integral -dm_view -quiet

1063:   # 3D serial P1 test 9-12
1064:   test:
1065:     suffix: 3d_p1_0
1066:     requires: ctetgen
1067:     args: -run_type test -dm_plex_dim 3 -bc_type dirichlet -dm_plex_interpolate 0 -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 -dm_view

1069:   test:
1070:     suffix: 3d_p1_1
1071:     requires: ctetgen
1072:     args: -run_type test -dm_plex_dim 3 -bc_type dirichlet -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 -dm_view

1074:   test:
1075:     suffix: 3d_p1_2
1076:     requires: ctetgen
1077:     args: -run_type test -dm_plex_dim 3 -dm_refine_volume_limit_pre 0.0125 -bc_type dirichlet -petscspace_degree 1 -show_initial -dm_plex_print_fem 1 -dm_view

1079:   test:
1080:     suffix: 3d_p1_neumann_0
1081:     requires: ctetgen
1082:     args: -run_type test -dm_plex_dim 3 -bc_type neumann -dm_plex_boundary_label boundary -petscspace_degree 1 -snes_fd -show_initial -dm_plex_print_fem 1 -dm_view

1084:   # Analytic variable coefficient 13-20
1085:   test:
1086:     suffix: 13
1087:     requires: triangle
1088:     args: -run_type test -variable_coefficient analytic -petscspace_degree 1 -show_initial -dm_plex_print_fem 1
1089:   test:
1090:     suffix: 14
1091:     requires: triangle
1092:     args: -run_type test -dm_refine_volume_limit_pre 0.0625 -variable_coefficient analytic -petscspace_degree 1 -show_initial -dm_plex_print_fem 1
1093:   test:
1094:     suffix: 15
1095:     requires: triangle
1096:     args: -run_type test -variable_coefficient analytic -petscspace_degree 2 -show_initial -dm_plex_print_fem 1
1097:   test:
1098:     suffix: 16
1099:     requires: triangle
1100:     args: -run_type test -dm_refine_volume_limit_pre 0.0625 -variable_coefficient analytic -petscspace_degree 2 -show_initial -dm_plex_print_fem 1
1101:   test:
1102:     suffix: 17
1103:     requires: ctetgen
1104:     args: -run_type test -dm_plex_dim 3 -variable_coefficient analytic -petscspace_degree 1 -show_initial -dm_plex_print_fem 1

1106:   test:
1107:     suffix: 18
1108:     requires: ctetgen
1109:     args: -run_type test -dm_plex_dim 3 -dm_refine_volume_limit_pre 0.0125 -variable_coefficient analytic -petscspace_degree 1 -show_initial -dm_plex_print_fem 1

1111:   test:
1112:     suffix: 19
1113:     requires: ctetgen
1114:     args: -run_type test -dm_plex_dim 3 -variable_coefficient analytic -petscspace_degree 2 -show_initial -dm_plex_print_fem 1

1116:   test:
1117:     suffix: 20
1118:     requires: ctetgen
1119:     args: -run_type test -dm_plex_dim 3 -dm_refine_volume_limit_pre 0.0125 -variable_coefficient analytic -petscspace_degree 2 -show_initial -dm_plex_print_fem 1

1121:   # P1 variable coefficient 21-28
1122:   test:
1123:     suffix: 21
1124:     requires: triangle
1125:     args: -run_type test -variable_coefficient field -petscspace_degree 1 -mat_petscspace_degree 1 -show_initial -dm_plex_print_fem 1

1127:   test:
1128:     suffix: 22
1129:     requires: triangle
1130:     args: -run_type test -dm_refine_volume_limit_pre 0.0625 -variable_coefficient field -petscspace_degree 1 -mat_petscspace_degree 1 -show_initial -dm_plex_print_fem 1

1132:   test:
1133:     suffix: 23
1134:     requires: triangle
1135:     args: -run_type test -variable_coefficient field -petscspace_degree 2 -mat_petscspace_degree 1 -show_initial -dm_plex_print_fem 1

1137:   test:
1138:     suffix: 24
1139:     requires: triangle
1140:     args: -run_type test -dm_refine_volume_limit_pre 0.0625 -variable_coefficient field -petscspace_degree 2 -mat_petscspace_degree 1 -show_initial -dm_plex_print_fem 1

1142:   test:
1143:     suffix: 25
1144:     requires: ctetgen
1145:     args: -run_type test -dm_plex_dim 3 -variable_coefficient field -petscspace_degree 1 -mat_petscspace_degree 1 -show_initial -dm_plex_print_fem 1

1147:   test:
1148:     suffix: 26
1149:     requires: ctetgen
1150:     args: -run_type test -dm_plex_dim 3 -dm_refine_volume_limit_pre 0.0125 -variable_coefficient field -petscspace_degree 1 -mat_petscspace_degree 1 -show_initial -dm_plex_print_fem 1

1152:   test:
1153:     suffix: 27
1154:     requires: ctetgen
1155:     args: -run_type test -dm_plex_dim 3 -variable_coefficient field -petscspace_degree 2 -mat_petscspace_degree 1 -show_initial -dm_plex_print_fem 1

1157:   test:
1158:     suffix: 28
1159:     requires: ctetgen
1160:     args: -run_type test -dm_plex_dim 3 -dm_refine_volume_limit_pre 0.0125 -variable_coefficient field -petscspace_degree 2 -mat_petscspace_degree 1 -show_initial -dm_plex_print_fem 1

1162:   # P0 variable coefficient 29-36
1163:   test:
1164:     suffix: 29
1165:     requires: triangle
1166:     args: -run_type test -variable_coefficient field -petscspace_degree 1 -show_initial -dm_plex_print_fem 1

1168:   test:
1169:     suffix: 30
1170:     requires: triangle
1171:     args: -run_type test -dm_refine_volume_limit_pre 0.0625 -variable_coefficient field -petscspace_degree 1 -show_initial -dm_plex_print_fem 1

1173:   test:
1174:     suffix: 31
1175:     requires: triangle
1176:     args: -run_type test -variable_coefficient field -petscspace_degree 2 -show_initial -dm_plex_print_fem 1

1178:   test:
1179:     requires: triangle
1180:     suffix: 32
1181:     args: -run_type test -dm_refine_volume_limit_pre 0.0625 -variable_coefficient field -petscspace_degree 2 -show_initial -dm_plex_print_fem 1

1183:   test:
1184:     requires: ctetgen
1185:     suffix: 33
1186:     args: -run_type test -dm_plex_dim 3 -variable_coefficient field -petscspace_degree 1 -show_initial -dm_plex_print_fem 1

1188:   test:
1189:     suffix: 34
1190:     requires: ctetgen
1191:     args: -run_type test -dm_plex_dim 3 -dm_refine_volume_limit_pre 0.0125 -variable_coefficient field -petscspace_degree 1 -show_initial -dm_plex_print_fem 1

1193:   test:
1194:     suffix: 35
1195:     requires: ctetgen
1196:     args: -run_type test -dm_plex_dim 3 -variable_coefficient field -petscspace_degree 2 -show_initial -dm_plex_print_fem 1

1198:   test:
1199:     suffix: 36
1200:     requires: ctetgen
1201:     args: -run_type test -dm_plex_dim 3 -dm_refine_volume_limit_pre 0.0125 -variable_coefficient field -petscspace_degree 2 -show_initial -dm_plex_print_fem 1

1203:   # Full solve 39-44
1204:   test:
1205:     suffix: 39
1206:     requires: triangle !single
1207:     args: -run_type full -dm_refine_volume_limit_pre 0.015625 -petscspace_degree 2 -pc_type gamg -pc_gamg_esteig_ksp_type cg -pc_gamg_esteig_ksp_max_it 10 -ksp_rtol 1.0e-10 -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason ::ascii_info_detail
1208:   test:
1209:     suffix: 40
1210:     requires: triangle !single
1211:     args: -run_type full -dm_refine_volume_limit_pre 0.015625 -variable_coefficient nonlinear -petscspace_degree 2 -pc_type svd -ksp_rtol 1.0e-10 -snes_monitor_short -snes_converged_reason ::ascii_info_detail
1212:   test:
1213:     suffix: 41
1214:     requires: triangle !single
1215:     args: -run_type full -dm_refine_volume_limit_pre 0.03125 -variable_coefficient nonlinear -petscspace_degree 1 -snes_type fas -snes_fas_levels 2 -fas_coarse_pc_type svd -fas_coarse_ksp_rtol 1.0e-10 -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -dm_refine_hierarchy 1 -snes_view -fas_levels_1_snes_type newtonls -fas_levels_1_pc_type svd -fas_levels_1_ksp_rtol 1.0e-10 -fas_levels_1_snes_monitor_short
1216:   test:
1217:     suffix: 42
1218:     requires: triangle !single
1219:     args: -run_type full -dm_refine_volume_limit_pre 0.0625 -variable_coefficient nonlinear -petscspace_degree 1 -snes_type fas -snes_fas_levels 3 -fas_coarse_pc_type svd -fas_coarse_ksp_rtol 1.0e-10 -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -dm_refine_hierarchy 2 -snes_view -fas_levels_1_snes_type newtonls -fas_levels_1_pc_type svd -fas_levels_1_ksp_rtol 1.0e-10 -fas_levels_1_snes_monitor_short -fas_levels_2_snes_type newtonls -fas_levels_2_pc_type svd -fas_levels_2_ksp_rtol 1.0e-10 -fas_levels_2_snes_atol 1.0e-11 -fas_levels_2_snes_monitor_short
1220:   test:
1221:     suffix: 43
1222:     requires: triangle !single
1223:     nsize: 2
1224:     args: -run_type full -dm_refine_volume_limit_pre 0.03125 -variable_coefficient nonlinear -petscspace_degree 1 -snes_type fas -snes_fas_levels 2 -fas_coarse_pc_type svd -fas_coarse_ksp_rtol 1.0e-10 -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -dm_refine_hierarchy 1 -snes_view -fas_levels_1_snes_type newtonls -fas_levels_1_pc_type svd -fas_levels_1_ksp_rtol 1.0e-10 -fas_levels_1_snes_monitor_short

1226:   test:
1227:     suffix: 44
1228:     requires: triangle !single
1229:     nsize: 2
1230:     args: -run_type full -dm_refine_volume_limit_pre 0.0625 -variable_coefficient nonlinear -petscspace_degree 1 -snes_type fas -snes_fas_levels 3 -fas_coarse_pc_type svd -fas_coarse_ksp_rtol 1.0e-10 -fas_coarse_snes_monitor_short -snes_monitor_short  -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -dm_refine_hierarchy 2 -dm_plex_print_fem 0 -snes_view -fas_levels_1_snes_type newtonls -fas_levels_1_pc_type svd -fas_levels_1_ksp_rtol 1.0e-10 -fas_levels_1_snes_monitor_short -fas_levels_2_snes_type newtonls -fas_levels_2_pc_type svd -fas_levels_2_ksp_rtol 1.0e-10 -fas_levels_2_snes_atol 1.0e-11 -fas_levels_2_snes_monitor_short

1232:   # These tests use a loose tolerance just to exercise the PtAP operations for MATIS and multiple PCBDDC setup calls inside PCMG
1233:   testset:
1234:     requires: triangle !single
1235:     nsize: 3
1236:     args: -run_type full -petscspace_degree 1 -dm_mat_type is -pc_type mg -pc_mg_levels 2 -mg_coarse_pc_type bddc -pc_mg_galerkin pmat -ksp_rtol 1.0e-2 -snes_converged_reason -dm_refine_hierarchy 2 -snes_max_it 4
1237:     test:
1238:       suffix: gmg_bddc
1239:       filter: sed -e "s/CONVERGED_FNORM_RELATIVE iterations 3/CONVERGED_FNORM_RELATIVE iterations 4/g"
1240:       args: -mg_levels_pc_type jacobi
1241:     test:
1242:       filter: sed -e "s/iterations [0-4]/iterations 4/g"
1243:       suffix: gmg_bddc_lev
1244:       args: -mg_levels_pc_type bddc

1246:   # Restarting
1247:   testset:
1248:     suffix: restart
1249:     requires: hdf5 triangle !complex
1250:     args: -run_type test -bc_type dirichlet -petscspace_degree 1
1251:     test:
1252:       args: -dm_view hdf5:sol.h5 -vec_view hdf5:sol.h5::append
1253:     test:
1254:       args: -dm_plex_filename sol.h5 -dm_plex_name box -restart

1256:   # Periodicity
1257:   test:
1258:     suffix: periodic_0
1259:     requires: triangle
1260:     args: -run_type full -bc_type dirichlet -petscspace_degree 1 -snes_converged_reason ::ascii_info_detail

1262:   test:
1263:     requires: !complex
1264:     suffix: periodic_1
1265:     args: -quiet -run_type test -dm_plex_simplex 0 -dm_plex_box_faces 3,3 -dm_plex_box_bd periodic,periodic -vec_view vtk:test.vtu:vtk_vtu -petscspace_degree 1 -dm_refine 1

1267:   # 2D serial P1 test with field bc
1268:   test:
1269:     suffix: field_bc_2d_p1_0
1270:     requires: triangle
1271:     args: -run_type test -bc_type dirichlet -field_bc -petscspace_degree 1 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1

1273:   test:
1274:     suffix: field_bc_2d_p1_1
1275:     requires: triangle
1276:     args: -run_type test -dm_refine 1 -bc_type dirichlet -field_bc -petscspace_degree 1 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1

1278:   test:
1279:     suffix: field_bc_2d_p1_neumann_0
1280:     requires: triangle
1281:     args: -run_type test -bc_type neumann -dm_plex_boundary_label boundary -field_bc -petscspace_degree 1 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1

1283:   test:
1284:     suffix: field_bc_2d_p1_neumann_1
1285:     requires: triangle
1286:     args: -run_type test -dm_refine 1 -bc_type neumann -dm_plex_boundary_label boundary -field_bc -petscspace_degree 1 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1

1288:   # 3D serial P1 test with field bc
1289:   test:
1290:     suffix: field_bc_3d_p1_0
1291:     requires: ctetgen
1292:     args: -run_type test -dm_plex_dim 3 -bc_type dirichlet -field_bc -petscspace_degree 1 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1

1294:   test:
1295:     suffix: field_bc_3d_p1_1
1296:     requires: ctetgen
1297:     args: -run_type test -dm_plex_dim 3 -dm_refine 1 -bc_type dirichlet -field_bc -petscspace_degree 1 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1

1299:   test:
1300:     suffix: field_bc_3d_p1_neumann_0
1301:     requires: ctetgen
1302:     args: -run_type test -dm_plex_dim 3 -bc_type neumann -dm_plex_boundary_label boundary -field_bc -petscspace_degree 1 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1

1304:   test:
1305:     suffix: field_bc_3d_p1_neumann_1
1306:     requires: ctetgen
1307:     args: -run_type test -dm_plex_dim 3 -dm_refine 1 -bc_type neumann -dm_plex_boundary_label boundary -field_bc -petscspace_degree 1 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1

1309:   # 2D serial P2 test with field bc
1310:   test:
1311:     suffix: field_bc_2d_p2_0
1312:     requires: triangle
1313:     args: -run_type test -bc_type dirichlet -field_bc -petscspace_degree 2 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1

1315:   test:
1316:     suffix: field_bc_2d_p2_1
1317:     requires: triangle
1318:     args: -run_type test -dm_refine 1 -bc_type dirichlet -field_bc -petscspace_degree 2 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1

1320:   test:
1321:     suffix: field_bc_2d_p2_neumann_0
1322:     requires: triangle
1323:     args: -run_type test -bc_type neumann -dm_plex_boundary_label boundary -field_bc -petscspace_degree 2 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1

1325:   test:
1326:     suffix: field_bc_2d_p2_neumann_1
1327:     requires: triangle
1328:     args: -run_type test -dm_refine 1 -bc_type neumann -dm_plex_boundary_label boundary -field_bc -petscspace_degree 2 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1

1330:   # 3D serial P2 test with field bc
1331:   test:
1332:     suffix: field_bc_3d_p2_0
1333:     requires: ctetgen
1334:     args: -run_type test -dm_plex_dim 3 -bc_type dirichlet -field_bc -petscspace_degree 2 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1

1336:   test:
1337:     suffix: field_bc_3d_p2_1
1338:     requires: ctetgen
1339:     args: -run_type test -dm_plex_dim 3 -dm_refine 1 -bc_type dirichlet -field_bc -petscspace_degree 2 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1

1341:   test:
1342:     suffix: field_bc_3d_p2_neumann_0
1343:     requires: ctetgen
1344:     args: -run_type test -dm_plex_dim 3 -bc_type neumann -dm_plex_boundary_label boundary -field_bc -petscspace_degree 2 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1

1346:   test:
1347:     suffix: field_bc_3d_p2_neumann_1
1348:     requires: ctetgen
1349:     args: -run_type test -dm_plex_dim 3 -dm_refine 1 -bc_type neumann -dm_plex_boundary_label boundary -field_bc -petscspace_degree 2 -bc_petscspace_degree 2 -show_initial -dm_plex_print_fem 1

1351:   # Full solve simplex: Convergence
1352:   test:
1353:     suffix: 3d_p1_conv
1354:     requires: ctetgen
1355:     args: -run_type full -dm_plex_dim 3 -dm_refine 1 -bc_type dirichlet -petscspace_degree 1 \
1356:       -snes_convergence_estimate -convest_num_refine 1 -pc_type lu

1358:   # Full solve simplex: PCBDDC
1359:   test:
1360:     suffix: tri_bddc
1361:     requires: triangle !single
1362:     nsize: 5
1363:     args: -run_type full -petscpartitioner_type simple -dm_refine 2 -bc_type dirichlet -petscspace_degree 1 -ksp_type gmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -dm_mat_type is -pc_type bddc -snes_monitor_short -ksp_monitor_short -snes_converged_reason ::ascii_info_detail -ksp_converged_reason -snes_view -show_solution 0

1365:   # Full solve simplex: PCBDDC
1366:   test:
1367:     suffix: tri_parmetis_bddc
1368:     requires: triangle !single parmetis
1369:     nsize: 4
1370:     args: -run_type full -petscpartitioner_type parmetis -dm_refine 2 -bc_type dirichlet -petscspace_degree 1 -ksp_type gmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -dm_mat_type is -pc_type bddc -snes_monitor_short -ksp_monitor_short -snes_converged_reason ::ascii_info_detail -ksp_converged_reason -snes_view -show_solution 0

1372:   testset:
1373:     args: -run_type full -dm_plex_simplex 0 -dm_plex_box_faces 3,3 -petscpartitioner_type simple -dm_refine 2 -bc_type dirichlet -petscspace_degree 2 -dm_mat_type is -pc_type bddc -ksp_type gmres -snes_monitor_short -ksp_monitor_short -snes_view -petscspace_poly_tensor -pc_bddc_corner_selection -ksp_rtol 1.e-9 -pc_bddc_use_edges 0
1374:     nsize: 5
1375:     output_file: output/ex12_quad_bddc.out
1376:     filter: sed -e "s/aijcusparse/aij/g" -e "s/aijviennacl/aij/g" -e "s/factorization: cusparse/factorization: petsc/g"
1377:     test:
1378:       requires: !single
1379:       suffix: quad_bddc
1380:     test:
1381:       requires: !single cuda
1382:       suffix: quad_bddc_cuda
1383:       args: -matis_localmat_type aijcusparse -pc_bddc_dirichlet_pc_factor_mat_solver_type cusparse -pc_bddc_neumann_pc_factor_mat_solver_type cusparse
1384:     test:
1385:       requires: !single viennacl
1386:       suffix: quad_bddc_viennacl
1387:       args: -matis_localmat_type aijviennacl

1389:   # Full solve simplex: ASM
1390:   test:
1391:     suffix: tri_q2q1_asm_lu
1392:     requires: triangle !single
1393:     args: -run_type full -dm_refine 3 -bc_type dirichlet -petscspace_degree 1 -ksp_type gmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -pc_type asm -pc_asm_type restrict -pc_asm_blocks 4 -sub_pc_type lu -snes_monitor_short -ksp_monitor_short -snes_converged_reason ::ascii_info_detail -ksp_converged_reason -snes_view -show_solution 0

1395:   test:
1396:     suffix: tri_q2q1_msm_lu
1397:     requires: triangle !single
1398:     args: -run_type full -dm_refine 3 -bc_type dirichlet -petscspace_degree 1 -ksp_type gmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -pc_type asm -pc_asm_type restrict -pc_asm_local_type multiplicative -pc_asm_blocks 4 -sub_pc_type lu -snes_monitor_short -ksp_monitor_short -snes_converged_reason ::ascii_info_detail -ksp_converged_reason -snes_view -show_solution 0

1400:   test:
1401:     suffix: tri_q2q1_asm_sor
1402:     requires: triangle !single
1403:     args: -run_type full -dm_refine 3 -bc_type dirichlet -petscspace_degree 1 -ksp_type gmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -pc_type asm -pc_asm_type restrict -pc_asm_blocks 4 -sub_pc_type sor -snes_monitor_short -ksp_monitor_short -snes_converged_reason ::ascii_info_detail -ksp_converged_reason -snes_view -show_solution 0

1405:   test:
1406:     suffix: tri_q2q1_msm_sor
1407:     requires: triangle !single
1408:     args: -run_type full -dm_refine 3 -bc_type dirichlet -petscspace_degree 1 -ksp_type gmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -pc_type asm -pc_asm_type restrict -pc_asm_local_type multiplicative -pc_asm_blocks 4 -sub_pc_type sor -snes_monitor_short -ksp_monitor_short -snes_converged_reason ::ascii_info_detail -ksp_converged_reason -snes_view -show_solution 0

1410:   # Full solve simplex: FAS
1411:   test:
1412:     suffix: fas_newton_0
1413:     requires: triangle !single
1414:     args: -run_type full -variable_coefficient nonlinear -petscspace_degree 1 -snes_type fas -snes_fas_levels 2 -fas_coarse_pc_type svd -fas_coarse_ksp_rtol 1.0e-10 -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -dm_refine_hierarchy 1 -snes_view -fas_levels_1_snes_type newtonls -fas_levels_1_pc_type svd -fas_levels_1_ksp_rtol 1.0e-10 -fas_levels_1_snes_monitor_short

1416:   test:
1417:     suffix: fas_newton_1
1418:     requires: triangle !single
1419:     args: -run_type full -dm_refine_hierarchy 3 -petscspace_degree 1 -snes_type fas -snes_fas_levels 3 -fas_coarse_pc_type lu -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_snes_linesearch_type basic -fas_levels_ksp_rtol 1.0e-10 -fas_levels_snes_monitor_short
1420:     filter: sed -e "s/total number of linear solver iterations=14/total number of linear solver iterations=15/g"

1422:   test:
1423:     suffix: fas_ngs_0
1424:     requires: triangle !single
1425:     args: -run_type full -variable_coefficient nonlinear -petscspace_degree 1 -snes_type fas -snes_fas_levels 2 -fas_coarse_pc_type svd -fas_coarse_ksp_rtol 1.0e-10 -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -dm_refine_hierarchy 1 -snes_view -fas_levels_1_snes_type ngs -fas_levels_1_snes_monitor_short

1427:   # These two tests are broken because DMPlexComputeInjectorFEM() only works for regularly refined meshes
1428:   test:
1429:     suffix: fas_newton_coarse_0
1430:     requires: pragmatic triangle
1431:     TODO: broken
1432:     args: -run_type full -variable_coefficient nonlinear -petscspace_degree 1 \
1433:           -dm_refine 2 -dm_coarsen_hierarchy 1 -dm_plex_hash_location -dm_adaptor pragmatic \
1434:           -snes_type fas -snes_fas_levels 2 -snes_converged_reason ::ascii_info_detail -snes_monitor_short -snes_view \
1435:             -fas_coarse_pc_type svd -fas_coarse_ksp_rtol 1.0e-10 -fas_coarse_snes_monitor_short -fas_coarse_snes_linesearch_type basic \
1436:             -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

1438:   test:
1439:     suffix: mg_newton_coarse_0
1440:     requires: triangle pragmatic
1441:     TODO: broken
1442:     args: -run_type full -petscspace_degree 1 \
1443:           -dm_refine 3 -dm_coarsen_hierarchy 3 -dm_plex_hash_location -dm_adaptor pragmatic \
1444:           -snes_atol 1.0e-8 -snes_rtol 0.0 -snes_monitor_short -snes_converged_reason ::ascii_info_detail -snes_view \
1445:             -ksp_type richardson -ksp_atol 1.0e-8 -ksp_rtol 0.0 -ksp_norm_type unpreconditioned -ksp_monitor_true_residual \
1446:               -pc_type mg -pc_mg_levels 4 \
1447:               -mg_levels_ksp_type gmres -mg_levels_pc_type ilu -mg_levels_ksp_max_it 10

1449:   # Full solve tensor
1450:   test:
1451:     suffix: tensor_plex_2d
1452:     args: -run_type test -dm_plex_simplex 0 -bc_type dirichlet -petscspace_degree 1 -dm_refine_hierarchy 2

1454:   test:
1455:     suffix: tensor_p4est_2d
1456:     requires: p4est
1457:     args: -run_type test -dm_plex_simplex 0 -bc_type dirichlet -petscspace_degree 1 -dm_forest_initial_refinement 2 -dm_forest_minimum_refinement 0 -dm_plex_convert_type p4est

1459:   test:
1460:     suffix: tensor_plex_3d
1461:     args: -run_type test -dm_plex_simplex 0 -bc_type dirichlet -petscspace_degree 1 -dm_plex_dim 3 -dm_refine_hierarchy 1 -dm_plex_box_faces 2,2,2

1463:   test:
1464:     suffix: tensor_p4est_3d
1465:     requires: p4est
1466:     args: -run_type test -dm_plex_simplex 0 -bc_type dirichlet -petscspace_degree 1 -dm_forest_initial_refinement 1 -dm_forest_minimum_refinement 0 -dm_plex_dim 3 -dm_plex_convert_type p8est -dm_plex_box_faces 2,2,2

1468:   test:
1469:     suffix: p4est_test_q2_conformal_serial
1470:     requires: p4est
1471:     args: -run_type test -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2

1473:   test:
1474:     suffix: p4est_test_q2_conformal_parallel
1475:     requires: p4est
1476:     nsize: 7
1477:     args: -run_type test -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -petscpartitioner_type simple

1479:   test:
1480:     suffix: p4est_test_q2_conformal_parallel_parmetis
1481:     requires: parmetis p4est
1482:     nsize: 4
1483:     args: -run_type test -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -petscpartitioner_type parmetis

1485:   test:
1486:     suffix: p4est_test_q2_nonconformal_serial
1487:     requires: p4est
1488:     filter: grep -v "CG or CGNE: variant"
1489:     args: -run_type test -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash

1491:   test:
1492:     suffix: p4est_test_q2_nonconformal_parallel
1493:     requires: p4est
1494:     filter: grep -v "CG or CGNE: variant"
1495:     nsize: 7
1496:     args: -run_type test -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type simple

1498:   test:
1499:     suffix: p4est_test_q2_nonconformal_parallel_parmetis
1500:     requires: parmetis p4est
1501:     nsize: 4
1502:     args: -run_type test -petscspace_degree 2 -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type parmetis

1504:   test:
1505:     suffix: p4est_exact_q2_conformal_serial
1506:     requires: p4est !single !complex !__float128
1507:     args: -run_type exact -petscspace_degree 2 -fas_levels_snes_atol 1.e-10 -snes_max_it 1 -snes_type fas -snes_fas_levels 3 -fas_coarse_pc_type none -fas_coarse_ksp_type preonly -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type none -fas_levels_ksp_type preonly -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2

1509:   test:
1510:     suffix: p4est_exact_q2_conformal_parallel
1511:     requires: p4est !single !complex !__float128
1512:     nsize: 4
1513:     args: -run_type exact -petscspace_degree 2 -fas_levels_snes_atol 1.e-10 -snes_max_it 1 -snes_type fas -snes_fas_levels 3 -fas_coarse_pc_type none -fas_coarse_ksp_type preonly -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type none -fas_levels_ksp_type preonly -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2

1515:   test:
1516:     suffix: p4est_exact_q2_conformal_parallel_parmetis
1517:     requires: parmetis p4est !single
1518:     nsize: 4
1519:     args: -run_type exact -petscspace_degree 2 -fas_levels_snes_atol 1.e-10 -snes_max_it 1 -snes_type fas -snes_fas_levels 3 -fas_coarse_pc_type none -fas_coarse_ksp_type preonly -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type none -fas_levels_ksp_type preonly -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -petscpartitioner_type parmetis

1521:   test:
1522:     suffix: p4est_exact_q2_nonconformal_serial
1523:     requires: p4est
1524:     args: -run_type exact -petscspace_degree 2 -fas_levels_snes_atol 1.e-10 -snes_max_it 1 -snes_type fas -snes_fas_levels 3 -fas_coarse_pc_type none -fas_coarse_ksp_type preonly -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type none -fas_levels_ksp_type preonly -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash

1526:   test:
1527:     suffix: p4est_exact_q2_nonconformal_parallel
1528:     requires: p4est
1529:     nsize: 7
1530:     args: -run_type exact -petscspace_degree 2 -fas_levels_snes_atol 1.e-10 -snes_max_it 1 -snes_type fas -snes_fas_levels 3 -fas_coarse_pc_type none -fas_coarse_ksp_type preonly -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type none -fas_levels_ksp_type preonly -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type simple

1532:   test:
1533:     suffix: p4est_exact_q2_nonconformal_parallel_parmetis
1534:     requires: parmetis p4est
1535:     nsize: 4
1536:     args: -run_type exact -petscspace_degree 2 -fas_levels_snes_atol 1.e-10 -snes_max_it 1 -snes_type fas -snes_fas_levels 3 -fas_coarse_pc_type none -fas_coarse_ksp_type preonly -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type none -fas_levels_ksp_type preonly -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type parmetis

1538:   test:
1539:     suffix: p4est_full_q2_nonconformal_serial
1540:     requires: p4est !single
1541:     filter: grep -v "variant HERMITIAN"
1542:     args: -run_type full -petscspace_degree 2 -snes_max_it 20 -snes_type fas -snes_fas_levels 3 -fas_coarse_pc_type jacobi -fas_coarse_ksp_type cg -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type jacobi -fas_levels_ksp_type cg -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash

1544:   test:
1545:     suffix: p4est_full_q2_nonconformal_parallel
1546:     requires: p4est !single
1547:     filter: grep -v "variant HERMITIAN"
1548:     nsize: 7
1549:     args: -run_type full -petscspace_degree 2 -snes_max_it 20 -snes_type fas -snes_fas_levels 3 -fas_coarse_pc_type jacobi -fas_coarse_ksp_type cg -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type jacobi -fas_levels_ksp_type cg -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type simple

1551:   test:
1552:     suffix: p4est_full_q2_nonconformal_parallel_bddcfas
1553:     requires: p4est !single
1554:     filter: grep -v "variant HERMITIAN"
1555:     nsize: 7
1556:     args: -run_type full -petscspace_degree 2 -snes_max_it 20 -snes_type fas -snes_fas_levels 3 -dm_mat_type is -fas_coarse_pc_type bddc -fas_coarse_ksp_type cg -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type bddc -fas_levels_ksp_type cg -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type simple

1558:   test:
1559:     suffix: p4est_full_q2_nonconformal_parallel_bddc
1560:     requires: p4est !single
1561:     filter: grep -v "variant HERMITIAN"
1562:     nsize: 7
1563:     args: -run_type full -petscspace_degree 2 -snes_max_it 20 -snes_type newtonls -dm_mat_type is -pc_type bddc -ksp_type cg -snes_monitor_short -snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type simple

1565:   test:
1566:     TODO: broken
1567:     suffix: p4est_fas_q2_conformal_serial
1568:     requires: p4est !complex !__float128
1569:     args: -run_type full -variable_coefficient nonlinear -petscspace_degree 2 -snes_max_it 20 -snes_type fas -snes_fas_levels 3 -pc_type jacobi -ksp_type gmres -fas_coarse_pc_type svd -fas_coarse_ksp_type gmres -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type svd -fas_levels_ksp_type gmres -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_refine_hierarchy 3

1571:   test:
1572:     TODO: broken
1573:     suffix: p4est_fas_q2_nonconformal_serial
1574:     requires: p4est
1575:     args: -run_type full -variable_coefficient nonlinear -petscspace_degree 2 -snes_max_it 20 -snes_type fas -snes_fas_levels 3 -pc_type jacobi -ksp_type gmres -fas_coarse_pc_type jacobi -fas_coarse_ksp_type gmres -fas_coarse_ksp_monitor_true_residual -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type jacobi -fas_levels_ksp_type gmres -fas_levels_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash

1577:   test:
1578:     suffix: fas_newton_0_p4est
1579:     requires: p4est !single !__float128
1580:     args: -run_type full -variable_coefficient nonlinear -petscspace_degree 1 -snes_type fas -snes_fas_levels 2 -fas_coarse_pc_type svd -fas_coarse_ksp_rtol 1.0e-10 -fas_coarse_snes_monitor_short -snes_monitor_short -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_1_snes_type newtonls -fas_levels_1_pc_type svd -fas_levels_1_ksp_rtol 1.0e-10 -fas_levels_1_snes_monitor_short -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash

1582:   # Full solve simplicial AMR
1583:   test:
1584:     suffix: tri_p1_adapt_init_pragmatic
1585:     requires: pragmatic
1586:     args: -run_type exact -dm_refine 5 -bc_type dirichlet -petscspace_degree 1 -variable_coefficient ball -snes_converged_reason ::ascii_info_detail -pc_type lu -snes_adapt_initial 1 -adaptor_target_num 4000 -dm_plex_metric_h_max 0.5 -dm_adaptor pragmatic

1588:   test:
1589:     suffix: tri_p2_adapt_init_pragmatic
1590:     requires: pragmatic
1591:     args: -run_type exact -dm_refine 5 -bc_type dirichlet -petscspace_degree 2 -variable_coefficient ball -snes_converged_reason ::ascii_info_detail -pc_type lu -snes_adapt_initial 1 -adaptor_target_num 4000 -dm_plex_metric_h_max 0.5 -dm_adaptor pragmatic

1593:   test:
1594:     suffix: tri_p1_adapt_init_mmg
1595:     requires: mmg
1596:     args: -run_type exact -dm_refine 5 -bc_type dirichlet -petscspace_degree 1 -variable_coefficient ball -snes_converged_reason ::ascii_info_detail -pc_type lu -snes_adapt_initial 1 -adaptor_target_num 4000 -dm_plex_metric_h_max 0.5 -dm_adaptor mmg

1598:   test:
1599:     suffix: tri_p2_adapt_init_mmg
1600:     requires: mmg
1601:     args: -run_type exact -dm_refine 5 -bc_type dirichlet -petscspace_degree 2 -variable_coefficient ball -snes_converged_reason ::ascii_info_detail -pc_type lu -snes_adapt_initial 1 -adaptor_target_num 4000 -dm_plex_metric_h_max 0.5 -dm_adaptor mmg

1603:   test:
1604:     suffix: tri_p1_adapt_seq_pragmatic
1605:     requires: pragmatic
1606:     args: -run_type exact -dm_refine 5 -bc_type dirichlet -petscspace_degree 1 -variable_coefficient ball -snes_converged_reason ::ascii_info_detail -pc_type lu -snes_adapt_sequence 2 -adaptor_target_num 4000 -dm_plex_metric_h_max 0.5 -dm_adaptor pragmatic

1608:   test:
1609:     suffix: tri_p2_adapt_seq_pragmatic
1610:     requires: pragmatic
1611:     args: -run_type exact -dm_refine 5 -bc_type dirichlet -petscspace_degree 2 -variable_coefficient ball -snes_converged_reason ::ascii_info_detail -pc_type lu -snes_adapt_sequence 2 -adaptor_target_num 4000 -dm_plex_metric_h_max 0.5 -dm_adaptor pragmatic

1613:   test:
1614:     suffix: tri_p1_adapt_seq_mmg
1615:     requires: mmg
1616:     args: -run_type exact -dm_refine 5 -bc_type dirichlet -petscspace_degree 1 -variable_coefficient ball -snes_converged_reason ::ascii_info_detail -pc_type lu -snes_adapt_sequence 2 -adaptor_target_num 4000 -dm_plex_metric_h_max 0.5 -dm_adaptor mmg

1618:   test:
1619:     suffix: tri_p2_adapt_seq_mmg
1620:     requires: mmg
1621:     args: -run_type exact -dm_refine 5 -bc_type dirichlet -petscspace_degree 2 -variable_coefficient ball -snes_converged_reason ::ascii_info_detail -pc_type lu -snes_adapt_sequence 2 -adaptor_target_num 4000 -dm_plex_metric_h_max 0.5 -dm_adaptor mmg

1623:   test:
1624:     suffix: tri_p1_adapt_analytic_pragmatic
1625:     requires: pragmatic
1626:     args: -run_type exact -dm_refine 3 -bc_type dirichlet -petscspace_degree 1 -variable_coefficient cross -snes_adapt_initial 4 -adaptor_target_num 500 -adaptor_monitor -dm_plex_metric_h_min 0.0001 -dm_plex_metric_h_max 0.05 -dm_adaptor pragmatic

1628:   test:
1629:     suffix: tri_p2_adapt_analytic_pragmatic
1630:     requires: pragmatic
1631:     args: -run_type exact -dm_refine 3 -bc_type dirichlet -petscspace_degree 2 -variable_coefficient cross -snes_adapt_initial 4 -adaptor_target_num 500 -adaptor_monitor -dm_plex_metric_h_min 0.0001 -dm_plex_metric_h_max 0.05 -dm_adaptor pragmatic

1633:   test:
1634:     suffix: tri_p1_adapt_analytic_mmg
1635:     requires: mmg
1636:     args: -run_type exact -dm_refine 3 -bc_type dirichlet -petscspace_degree 1 -variable_coefficient cross -snes_adapt_initial 4 -adaptor_target_num 500 -adaptor_monitor -dm_plex_metric_h_max 0.5 -dm_adaptor mmg

1638:   test:
1639:     suffix: tri_p2_adapt_analytic_mmg
1640:     requires: mmg
1641:     args: -run_type exact -dm_refine 3 -bc_type dirichlet -petscspace_degree 2 -variable_coefficient cross -snes_adapt_initial 4 -adaptor_target_num 500 -adaptor_monitor -dm_plex_metric_h_max 0.5 -dm_adaptor mmg

1643:   test:
1644:     suffix: tri_p1_adapt_uniform_pragmatic
1645:     requires: pragmatic tetgen
1646:     nsize: 2
1647:     args: -run_type full -dm_plex_box_faces 8,8,8 -bc_type dirichlet -petscspace_degree 1 -variable_coefficient none -snes_converged_reason ::ascii_info_detail -ksp_type cg -pc_type sor -snes_adapt_sequence 3 -adaptor_target_num 400 -dm_plex_metric_h_max 0.5 -dm_plex_dim 3 -dm_adaptor pragmatic
1648:     timeoutfactor: 2

1650:   test:
1651:     suffix: tri_p2_adapt_uniform_pragmatic
1652:     requires: pragmatic tetgen
1653:     nsize: 2
1654:     args: -run_type full -dm_plex_box_faces 8,8,8 -bc_type dirichlet -petscspace_degree 2 -variable_coefficient none -snes_converged_reason ::ascii_info_detail -ksp_type cg -pc_type sor -snes_adapt_sequence 1 -adaptor_target_num 400 -dm_plex_metric_h_max 0.5 -dm_plex_dim 3 -dm_adaptor pragmatic
1655:     timeoutfactor: 1

1657:   test:
1658:     suffix: tri_p1_adapt_uniform_mmg
1659:     requires: mmg tetgen
1660:     args: -run_type full -dm_plex_box_faces 4,4,4 -bc_type dirichlet -petscspace_degree 1 -variable_coefficient none -snes_converged_reason ::ascii_info_detail -ksp_type cg -pc_type sor -snes_adapt_sequence 3 -adaptor_target_num 400 -dm_plex_metric_h_max 0.5 -dm_plex_dim 3 -dm_adaptor mmg
1661:     timeoutfactor: 2

1663:   test:
1664:     suffix: tri_p2_adapt_uniform_mmg
1665:     requires: mmg tetgen
1666:     args: -run_type full -dm_plex_box_faces 4,4,4 -bc_type dirichlet -petscspace_degree 2 -variable_coefficient none -snes_converged_reason ::ascii_info_detail -ksp_type cg -pc_type sor -snes_adapt_sequence 1 -adaptor_target_num 400 -dm_plex_metric_h_max 0.5 -dm_plex_dim 3 -dm_adaptor mmg
1667:     timeoutfactor: 1

1669:   test:
1670:     suffix: tri_p1_adapt_uniform_parmmg
1671:     requires: parmmg tetgen
1672:     nsize: 2
1673:     args: -run_type full -dm_plex_box_faces 8,8,8 -bc_type dirichlet -petscspace_degree 1 -variable_coefficient none -snes_converged_reason ::ascii_info_detail -ksp_type cg -pc_type sor -snes_adapt_sequence 3 -adaptor_target_num 400 -dm_plex_metric_h_max 0.5 -dm_plex_dim 3 -dm_adaptor parmmg
1674:     timeoutfactor: 2

1676:   test:
1677:     suffix: tri_p2_adapt_uniform_parmmg
1678:     requires: parmmg tetgen
1679:     nsize: 2
1680:     args: -run_type full -dm_plex_box_faces 8,8,8 -bc_type dirichlet -petscspace_degree 2 -variable_coefficient none -snes_converged_reason ::ascii_info_detail -ksp_type cg -pc_type sor -snes_adapt_sequence 1 -adaptor_target_num 400 -dm_plex_metric_h_max 0.5 -dm_plex_dim 3 -dm_adaptor parmmg
1681:     timeoutfactor: 1

1683:   # Full solve tensor AMR
1684:   test:
1685:     suffix: quad_q1_adapt_0
1686:     requires: p4est
1687:     args: -run_type exact -dm_plex_simplex 0 -dm_plex_convert_type p4est -bc_type dirichlet -petscspace_degree 1 -variable_coefficient ball -snes_converged_reason ::ascii_info_detail -pc_type lu -dm_forest_initial_refinement 4 -snes_adapt_initial 1 -dm_view
1688:     filter: grep -v DM_

1690:   test:
1691:     suffix: amr_0
1692:     nsize: 5
1693:     args: -run_type test -petscpartitioner_type simple -dm_plex_simplex 0 -bc_type dirichlet -petscspace_degree 1 -dm_refine 1

1695:   test:
1696:     suffix: amr_1
1697:     requires: p4est !complex
1698:     args: -run_type test -dm_plex_simplex 0 -bc_type dirichlet -petscspace_degree 1 -dm_plex_convert_type p4est -dm_p4est_refine_pattern center -dm_forest_maximum_refinement 5 -dm_view vtk:amr.vtu:vtk_vtu -vec_view vtk:amr.vtu:vtk_vtu:append

1700:   test:
1701:     suffix: p4est_solve_bddc
1702:     requires: p4est !complex
1703:     args: -run_type full -variable_coefficient nonlinear -nonzero_initial_guess 1 -petscspace_degree 2 -snes_max_it 20 -snes_type newtonls -dm_mat_type is -pc_type bddc -ksp_type cg -snes_monitor_short -ksp_monitor -snes_linesearch_type bt -snes_converged_reason -snes_view -dm_plex_simplex 0 -petscspace_poly_tensor -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type simple -pc_bddc_detect_disconnected
1704:     nsize: 4

1706:   test:
1707:     suffix: p4est_solve_fas
1708:     requires: p4est
1709:     args: -run_type full -variable_coefficient nonlinear -nonzero_initial_guess 1 -petscspace_degree 2 -snes_max_it 10 -snes_type fas -snes_linesearch_type bt -snes_fas_levels 3 -fas_coarse_snes_type newtonls -fas_coarse_snes_linesearch_type basic -fas_coarse_ksp_type cg -fas_coarse_pc_type jacobi -fas_coarse_snes_monitor_short -fas_levels_snes_max_it 4 -fas_levels_snes_type newtonls -fas_levels_snes_linesearch_type bt -fas_levels_ksp_type cg -fas_levels_pc_type jacobi -fas_levels_snes_monitor_short -fas_levels_cycle_snes_linesearch_type bt -snes_monitor_short -snes_converged_reason -snes_view -dm_plex_simplex 0 -petscspace_poly_tensor -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash
1710:     nsize: 4
1711:     TODO: identical machine two runs produce slightly different solver trackers

1713:   test:
1714:     suffix: p4est_convergence_test_1
1715:     requires: p4est
1716:     args:  -quiet -run_type test -petscspace_degree 1 -dm_plex_simplex 0 -petscspace_poly_tensor -dm_plex_convert_type p4est -dm_forest_minimum_refinement 2 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash
1717:     nsize: 4

1719:   test:
1720:     suffix: p4est_convergence_test_2
1721:     requires: p4est
1722:     args: -quiet -run_type test -petscspace_degree 1 -dm_plex_simplex 0 -petscspace_poly_tensor -dm_plex_convert_type p4est -dm_forest_minimum_refinement 3 -dm_forest_initial_refinement 3 -dm_forest_maximum_refinement 5 -dm_p4est_refine_pattern hash

1724:   test:
1725:     suffix: p4est_convergence_test_3
1726:     requires: p4est
1727:     args: -quiet -run_type test -petscspace_degree 1 -dm_plex_simplex 0 -petscspace_poly_tensor -dm_plex_convert_type p4est -dm_forest_minimum_refinement 4 -dm_forest_initial_refinement 4 -dm_forest_maximum_refinement 6 -dm_p4est_refine_pattern hash

1729:   test:
1730:     suffix: p4est_convergence_test_4
1731:     requires: p4est
1732:     args: -quiet -run_type test -petscspace_degree 1 -dm_plex_simplex 0 -petscspace_poly_tensor -dm_plex_convert_type p4est -dm_forest_minimum_refinement 5 -dm_forest_initial_refinement 5 -dm_forest_maximum_refinement 7 -dm_p4est_refine_pattern hash
1733:     timeoutfactor: 5

1735:   # Serial tests with GLVis visualization
1736:   test:
1737:     suffix: glvis_2d_tet_p1
1738:     args: -quiet -run_type test -bc_type dirichlet -petscspace_degree 1 -vec_view glvis: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -dm_plex_boundary_label marker -dm_plex_gmsh_periodic 0 -dm_coord_space 0
1739:   test:
1740:     suffix: glvis_2d_tet_p2
1741:     args: -quiet -run_type test -bc_type dirichlet -petscspace_degree 2 -vec_view glvis: -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -dm_plex_boundary_label marker -dm_plex_gmsh_periodic 0 -dm_coord_space 0
1742:   test:
1743:     suffix: glvis_2d_hex_p1
1744:     args: -quiet -run_type test -bc_type dirichlet -petscspace_degree 1 -vec_view glvis: -dm_plex_simplex 0 -dm_refine 1 -dm_coord_space 0
1745:   test:
1746:     suffix: glvis_2d_hex_p2
1747:     args: -quiet -run_type test -bc_type dirichlet -petscspace_degree 2 -vec_view glvis: -dm_plex_simplex 0 -dm_refine 1 -dm_coord_space 0
1748:   test:
1749:     suffix: glvis_2d_hex_p2_p4est
1750:     requires: p4est
1751:     args: -quiet -run_type test -bc_type dirichlet -petscspace_degree 2 -vec_view glvis: -dm_plex_simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 1 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -viewer_glvis_dm_plex_enable_ncmesh
1752:   test:
1753:     suffix: glvis_2d_tet_p0
1754:     args: -run_type exact -guess_vec_view glvis: -nonzero_initial_guess 1 -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -dm_plex_boundary_label marker -petscspace_degree 0 -dm_coord_space 0
1755:   test:
1756:     suffix: glvis_2d_hex_p0
1757:     args: -run_type exact -guess_vec_view glvis: -nonzero_initial_guess 1 -dm_plex_box_faces 5,7 -dm_plex_simplex 0 -petscspace_degree 0 -dm_coord_space 0

1759:   # PCHPDDM tests
1760:   testset:
1761:     nsize: 4
1762:     requires: hpddm slepc !single defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
1763:     args: -run_type test -run_test_check_ksp -quiet -petscspace_degree 1 -petscpartitioner_type simple -bc_type none -dm_plex_simplex 0 -pc_type hpddm -pc_hpddm_levels_1_sub_pc_type lu -pc_hpddm_levels_1_eps_nev 2 -pc_hpddm_coarse_p 1 -pc_hpddm_coarse_pc_type svd -ksp_rtol 1.e-10 -pc_hpddm_levels_1_st_pc_factor_shift_type INBLOCKS -ksp_converged_reason
1764:     test:
1765:       suffix: quad_singular_hpddm
1766:       args: -dm_plex_box_faces 6,7
1767:     test:
1768:       requires: p4est
1769:       suffix: p4est_singular_2d_hpddm
1770:       args: -dm_plex_convert_type p4est -dm_forest_minimum_refinement 1 -dm_forest_initial_refinement 3 -dm_forest_maximum_refinement 3
1771:     test:
1772:       requires: p4est
1773:       suffix: p4est_nc_singular_2d_hpddm
1774:       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
1775:   testset:
1776:     nsize: 4
1777:     requires: hpddm slepc triangle !single defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
1778:     args: -run_type full -petscpartitioner_type simple -dm_refine 2 -bc_type dirichlet -petscspace_degree 2 -ksp_type gmres -ksp_gmres_restart 100 -pc_type hpddm -snes_monitor_short -ksp_monitor_short -snes_converged_reason ::ascii_info_detail -ksp_converged_reason -snes_view -show_solution 0 -pc_type hpddm -pc_hpddm_levels_1_sub_pc_type lu -pc_hpddm_levels_1_eps_nev 4 -pc_hpddm_coarse_p 2 -pc_hpddm_coarse_pc_type redundant -ksp_rtol 1.e-1
1779:     test:
1780:       args: -pc_hpddm_coarse_mat_type baij -options_left no
1781:       suffix: tri_hpddm_reuse_baij
1782:     test:
1783:       requires: !complex
1784:       suffix: tri_hpddm_reuse
1785:   testset:
1786:     nsize: 4
1787:     requires: hpddm slepc !single defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
1788:     args: -run_type full -petscpartitioner_type simple -dm_plex_box_faces 7,5 -dm_refine 2 -dm_plex_simplex 0 -bc_type dirichlet -petscspace_degree 2 -ksp_type gmres -ksp_gmres_restart 100 -pc_type hpddm -snes_monitor_short -ksp_monitor_short -snes_converged_reason ::ascii_info_detail -ksp_converged_reason -snes_view -show_solution 0 -pc_type hpddm -pc_hpddm_levels_1_sub_pc_type lu -pc_hpddm_levels_1_eps_nev 4 -pc_hpddm_coarse_p 2 -pc_hpddm_coarse_pc_type redundant -ksp_rtol 1.e-1
1789:     test:
1790:       args: -pc_hpddm_coarse_mat_type baij -options_left no
1791:       suffix: quad_hpddm_reuse_baij
1792:     test:
1793:       requires: !complex
1794:       suffix: quad_hpddm_reuse
1795:   testset:
1796:     nsize: 4
1797:     requires: hpddm slepc !single defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
1798:     args: -run_type full -petscpartitioner_type simple -dm_plex_box_faces 7,5 -dm_refine 2 -dm_plex_simplex 0 -bc_type dirichlet -petscspace_degree 1 -ksp_type gmres -ksp_gmres_restart 100 -pc_type hpddm -snes_monitor_short -ksp_monitor_short -snes_converged_reason ::ascii_info_detail -ksp_converged_reason -snes_view -show_solution 0 -pc_type hpddm -pc_hpddm_levels_1_sub_pc_type lu -pc_hpddm_levels_1_eps_threshold 0.1 -pc_hpddm_coarse_p 2 -pc_hpddm_coarse_pc_type redundant -ksp_rtol 1.e-1
1799:     test:
1800:       args: -pc_hpddm_coarse_mat_type baij -options_left no
1801:       suffix: quad_hpddm_reuse_threshold_baij
1802:     test:
1803:       requires: !complex
1804:       suffix: quad_hpddm_reuse_threshold
1805:   testset:
1806:     nsize: 4
1807:     requires: hpddm slepc parmetis !single defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
1808:     filter: sed -e "s/linear solver iterations=17/linear solver iterations=16/g"
1809:     args: -run_type full -petscpartitioner_type parmetis -dm_refine 3 -bc_type dirichlet -petscspace_degree 1 -ksp_type gmres -ksp_gmres_restart 100 -pc_type hpddm -snes_monitor_short -snes_converged_reason ::ascii_info_detail -snes_view -show_solution 0 -pc_type hpddm -pc_hpddm_levels_1_sub_pc_type icc -pc_hpddm_levels_1_eps_nev 20 -pc_hpddm_coarse_p 2 -pc_hpddm_coarse_pc_type redundant -ksp_rtol 1.e-10 -dm_plex_filename ${PETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -dm_plex_boundary_label marker -pc_hpddm_levels_1_sub_pc_factor_levels 3 -variable_coefficient ball -dm_plex_gmsh_periodic 0
1810:     test:
1811:       args: -pc_hpddm_coarse_mat_type baij -options_left no
1812:       filter: grep -v "      total: nonzeros=" | grep -v "      rows=" | sed -e "s/total number of linear solver iterations=1[5-7]/total number of linear solver iterations=16/g"
1813:       suffix: tri_parmetis_hpddm_baij
1814:     test:
1815:       filter: grep -v "      total: nonzeros=" | grep -v "      rows=" | sed -e "s/total number of linear solver iterations=1[5-7]/total number of linear solver iterations=16/g"
1816:       requires: !complex
1817:       suffix: tri_parmetis_hpddm

1819:   # 2D serial P1 tests for adaptive MG
1820:   test:
1821:     suffix: 2d_p1_adaptmg_0
1822:     requires: triangle bamg
1823:     args: -dm_refine_hierarchy 3 -dm_plex_box_faces 4,4 -bc_type dirichlet -petscspace_degree 1 \
1824:           -variable_coefficient checkerboard_0 -mat_petscspace_degree 0 -div 16 -k 3 \
1825:           -snes_max_it 1 -ksp_converged_reason \
1826:           -ksp_rtol 1e-8 -pc_type mg
1827:   # -ksp_monitor_true_residual -ksp_converged_reason -mg_levels_ksp_monitor_true_residual -pc_mg_mesp_monitor -dm_adapt_interp_view_fine draw -dm_adapt_interp_view_coarse draw -draw_pause 1
1828:   test:
1829:     suffix: 2d_p1_adaptmg_1
1830:     requires: triangle bamg
1831:     args: -dm_refine_hierarchy 3 -dm_plex_box_faces 4,4 -bc_type dirichlet -petscspace_degree 1 \
1832:           -variable_coefficient checkerboard_0 -mat_petscspace_degree 0 -div 16 -k 3 \
1833:           -snes_max_it 1 -ksp_converged_reason \
1834:           -ksp_rtol 1e-8 -pc_type mg -pc_mg_galerkin -pc_mg_adapt_interp -pc_mg_adapt_interp_coarse_space eigenvector -pc_mg_adapt_interp_n 1 \
1835:             -pc_mg_mesp_ksp_type richardson -pc_mg_mesp_ksp_richardson_self_scale -pc_mg_mesp_ksp_max_it 100 -pc_mg_mesp_pc_type none

1837: TEST*/