Actual source code: ex12.c

petsc-3.9.4 2018-09-11
Report Typos and Errors
  1: static char help[] = "Poisson Problem in 2d and 3d with simplicial finite elements.\n\
  2: We solve the Poisson problem in a rectangular\n\
  3: domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\
  4: This example supports discretized auxiliary fields (conductivity) as well as\n\
  5: multilevel nonlinear solvers.\n\n\n";

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

 10:   -dm_adapt_view hdf5:$PWD/adapt.h5 -sol_adapt_view hdf5:$PWD/adapt.h5::append -dm_adapt_pre_view hdf5:$PWD/orig.h5 -sol_adapt_pre_view hdf5:$PWD/orig.h5::append

 12: Information on refinement:

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

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

 23: typedef enum {NEUMANN, DIRICHLET, NONE} BCType;
 24: typedef enum {RUN_FULL, RUN_EXACT, RUN_TEST, RUN_PERF} RunType;
 25: typedef enum {COEFF_NONE, COEFF_ANALYTIC, COEFF_FIELD, COEFF_NONLINEAR, COEFF_CIRCLE, COEFF_CROSS} CoeffType;

 27: typedef struct {
 28:   PetscInt       debug;             /* The debugging level */
 29:   RunType        runType;           /* Whether to run tests, or solve the full problem */
 30:   PetscBool      jacobianMF;        /* Whether to calculate the Jacobian action on the fly */
 31:   PetscLogEvent  createMeshEvent;
 32:   PetscBool      showInitial, showSolution, restart, check, quiet, nonzInit;
 33:   /* Domain and mesh definition */
 34:   PetscInt       dim;               /* The topological mesh dimension */
 35:   DMBoundaryType periodicity[3];    /* The domain periodicity */
 36:   PetscInt       cells[3];          /* The initial domain division */
 37:   char           filename[2048];    /* The optional mesh file */
 38:   PetscBool      interpolate;       /* Generate intermediate mesh elements */
 39:   PetscReal      refinementLimit;   /* The largest allowable cell volume */
 40:   PetscBool      viewHierarchy;     /* Whether to view the hierarchy */
 41:   PetscBool      simplex;           /* Simplicial mesh */
 42:   /* Problem definition */
 43:   BCType         bcType;
 44:   CoeffType      variableCoefficient;
 45:   PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx);
 46:   PetscBool      fieldBC;
 47:   void           (**exactFields)(PetscInt, PetscInt, PetscInt,
 48:                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
 49:                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
 50:                                  PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);
 51:   /* Solver */
 52:   PC            pcmg;              /* This is needed for error monitoring */
 53: } AppCtx;

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

 61: static PetscErrorCode ecks(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
 62: {
 63:   u[0] = x[0];
 64:   return 0;
 65: }

 67: /*
 68:   In 2D for Dirichlet conditions, we use exact solution:

 70:     u = x^2 + y^2
 71:     f = 4

 73:   so that

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

 77:   For Neumann conditions, we have

 79:     -\nabla u \cdot -\hat y |_{y=0} =  (2y)|_{y=0} =  0 (bottom)
 80:     -\nabla u \cdot  \hat y |_{y=1} = -(2y)|_{y=1} = -2 (top)
 81:     -\nabla u \cdot -\hat x |_{x=0} =  (2x)|_{x=0} =  0 (left)
 82:     -\nabla u \cdot  \hat x |_{x=1} = -(2x)|_{x=1} = -2 (right)

 84:   Which we can express as

 86:     \nabla u \cdot  \hat n|_\Gamma = {2 x, 2 y} \cdot \hat n = 2 (x + y)
 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 circle_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 = PetscSinScalar(alpha*xy) * (alpha*PetscAbsScalar(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_circle_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:   const PetscReal alpha   = 500.;
136:   const PetscReal radius2 = PetscSqr(0.15);
137:   const PetscReal r2      = PetscSqr(x[0] - 0.5) + PetscSqr(x[1] - 0.5);
138:   const PetscReal xi      = alpha*(radius2 - r2);

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

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

151:   f0[0] = PetscSinScalar(alpha*xy) * (alpha*PetscAbsScalar(xy) < 2*PETSC_PI ? 1 : 0.01);
152: }

154: static void f0_bd_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
155:                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
156:                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
157:                     PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
158: {
159:   PetscInt d;
160:   for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += -n[d]*2.0*x[d];
161: }

163: static void f1_bd_zero(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 f1[])
167: {
168:   PetscInt comp;
169:   for (comp = 0; comp < dim; ++comp) f1[comp] = 0.0;
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: void f0_analytic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
259:                    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
260:                    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
261:                    PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
262: {
263:   f0[0] = 6.0*(x[0] + x[1]);
264: }

266: /* gradU[comp*dim+d] = {u_x, u_y} or {u_x, u_y, u_z} */
267: void f1_analytic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
268:                    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
269:                    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
270:                    PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
271: {
272:   PetscInt d;
273:   for (d = 0; d < dim; ++d) f1[d] = (x[0] + x[1])*u_x[d];
274: }

276: void f1_field_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
277:                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
278:                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
279:                 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
280: {
281:   PetscInt d;
282:   for (d = 0; d < dim; ++d) f1[d] = a[0]*u_x[d];
283: }

285: /* < \nabla v, \nabla u + {\nabla u}^T >
286:    This just gives \nabla u, give the perdiagonal for the transpose */
287: void g3_analytic_uu(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, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
291: {
292:   PetscInt d;
293:   for (d = 0; d < dim; ++d) g3[d*dim+d] = x[0] + x[1];
294: }

296: void g3_field_uu(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, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
300: {
301:   PetscInt d;
302:   for (d = 0; d < dim; ++d) g3[d*dim+d] = a[0];
303: }

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

308:     u  = x^2 + y^2
309:     f  = 16 (x^2 + y^2)
310:     nu = 1/2 |grad u|^2

312:   so that

314:     -\div \nu \grad u + f = -16 (x^2 + y^2) + 16 (x^2 + y^2) = 0
315: */
316: void f0_analytic_nonlinear_u(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, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
320: {
321:   f0[0] = 16.0*(x[0]*x[0] + x[1]*x[1]);
322: }

324: /* gradU[comp*dim+d] = {u_x, u_y} or {u_x, u_y, u_z} */
325: void f1_analytic_nonlinear_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
326:                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
327:                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
328:                              PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
329: {
330:   PetscScalar nu = 0.0;
331:   PetscInt    d;
332:   for (d = 0; d < dim; ++d) nu += u_x[d]*u_x[d];
333:   for (d = 0; d < dim; ++d) f1[d] = 0.5*nu*u_x[d];
334: }

336: /*
337:   grad (u + eps w) - grad u = eps grad w

339:   1/2 |grad (u + eps w)|^2 grad (u + eps w) - 1/2 |grad u|^2 grad u
340: = 1/2 (|grad u|^2 + 2 eps <grad u,grad w>) (grad u + eps grad w) - 1/2 |grad u|^2 grad u
341: = 1/2 (eps |grad u|^2 grad w + 2 eps <grad u,grad w> grad u)
342: = eps (1/2 |grad u|^2 grad w + grad u <grad u,grad w>)
343: */
344: void g3_analytic_nonlinear_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
345:                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
346:                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
347:                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
348: {
349:   PetscScalar nu = 0.0;
350:   PetscInt    d, e;
351:   for (d = 0; d < dim; ++d) nu += u_x[d]*u_x[d];
352:   for (d = 0; d < dim; ++d) {
353:     g3[d*dim+d] = 0.5*nu;
354:     for (e = 0; e < dim; ++e) {
355:       g3[d*dim+e] += u_x[d]*u_x[e];
356:     }
357:   }
358: }

360: /*
361:   In 3D for Dirichlet conditions we use exact solution:

363:     u = 2/3 (x^2 + y^2 + z^2)
364:     f = 4

366:   so that

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

370:   For Neumann conditions, we have

372:     -\nabla u \cdot -\hat z |_{z=0} =  (2z)|_{z=0} =  0 (bottom)
373:     -\nabla u \cdot  \hat z |_{z=1} = -(2z)|_{z=1} = -2 (top)
374:     -\nabla u \cdot -\hat y |_{y=0} =  (2y)|_{y=0} =  0 (front)
375:     -\nabla u \cdot  \hat y |_{y=1} = -(2y)|_{y=1} = -2 (back)
376:     -\nabla u \cdot -\hat x |_{x=0} =  (2x)|_{x=0} =  0 (left)
377:     -\nabla u \cdot  \hat x |_{x=1} = -(2x)|_{x=1} = -2 (right)

379:   Which we can express as

381:     \nabla u \cdot  \hat n|_\Gamma = {2 x, 2 y, 2z} \cdot \hat n = 2 (x + y + z)
382: */
383: static PetscErrorCode quadratic_u_3d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
384: {
385:   *u = 2.0*(x[0]*x[0] + x[1]*x[1] + x[2]*x[2])/3.0;
386:   return 0;
387: }

389: static void quadratic_u_field_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux,
390:                                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
391:                                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
392:                                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar uexact[])
393: {
394:   uexact[0] = a[0];
395: }

397: static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
398: {
399:   const char    *bcTypes[3]  = {"neumann", "dirichlet", "none"};
400:   const char    *runTypes[4] = {"full", "exact", "test", "perf"};
401:   const char    *coeffTypes[6] = {"none", "analytic", "field", "nonlinear", "circle", "cross"};
402:   PetscInt       bd, bc, run, coeff, n;
403:   PetscBool      flg;

407:   options->debug               = 0;
408:   options->runType             = RUN_FULL;
409:   options->dim                 = 2;
410:   options->periodicity[0]      = DM_BOUNDARY_NONE;
411:   options->periodicity[1]      = DM_BOUNDARY_NONE;
412:   options->periodicity[2]      = DM_BOUNDARY_NONE;
413:   options->cells[0]            = 2;
414:   options->cells[1]            = 2;
415:   options->cells[2]            = 2;
416:   options->filename[0]         = '\0';
417:   options->interpolate         = PETSC_FALSE;
418:   options->refinementLimit     = 0.0;
419:   options->bcType              = DIRICHLET;
420:   options->variableCoefficient = COEFF_NONE;
421:   options->fieldBC             = PETSC_FALSE;
422:   options->jacobianMF          = PETSC_FALSE;
423:   options->showInitial         = PETSC_FALSE;
424:   options->showSolution        = PETSC_FALSE;
425:   options->restart             = PETSC_FALSE;
426:   options->check               = PETSC_FALSE;
427:   options->viewHierarchy       = PETSC_FALSE;
428:   options->simplex             = PETSC_TRUE;
429:   options->quiet               = PETSC_FALSE;
430:   options->nonzInit            = PETSC_FALSE;

432:   PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX");
433:   PetscOptionsInt("-debug", "The debugging level", "ex12.c", options->debug, &options->debug, NULL);
434:   run  = options->runType;
435:   PetscOptionsEList("-run_type", "The run type", "ex12.c", runTypes, 3, runTypes[options->runType], &run, NULL);

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

439:   PetscOptionsInt("-dim", "The topological mesh dimension", "ex12.c", options->dim, &options->dim, NULL);
440:   bd = options->periodicity[0];
441:   PetscOptionsEList("-x_periodicity", "The x-boundary periodicity", "ex12.c", DMBoundaryTypes, 5, DMBoundaryTypes[options->periodicity[0]], &bd, NULL);
442:   options->periodicity[0] = (DMBoundaryType) bd;
443:   bd = options->periodicity[1];
444:   PetscOptionsEList("-y_periodicity", "The y-boundary periodicity", "ex12.c", DMBoundaryTypes, 5, DMBoundaryTypes[options->periodicity[1]], &bd, NULL);
445:   options->periodicity[1] = (DMBoundaryType) bd;
446:   bd = options->periodicity[2];
447:   PetscOptionsEList("-z_periodicity", "The z-boundary periodicity", "ex12.c", DMBoundaryTypes, 5, DMBoundaryTypes[options->periodicity[2]], &bd, NULL);
448:   options->periodicity[2] = (DMBoundaryType) bd;
449:   n = 3;
450:   PetscOptionsIntArray("-cells", "The initial mesh division", "ex12.c", options->cells, &n, NULL);
451:   PetscOptionsString("-f", "Mesh filename to read", "ex12.c", options->filename, options->filename, sizeof(options->filename), &flg);
452:   PetscOptionsBool("-interpolate", "Generate intermediate mesh elements", "ex12.c", options->interpolate, &options->interpolate, NULL);
453:   PetscOptionsReal("-refinement_limit", "The largest allowable cell volume", "ex12.c", options->refinementLimit, &options->refinementLimit, NULL);
454:   bc   = options->bcType;
455:   PetscOptionsEList("-bc_type","Type of boundary condition","ex12.c",bcTypes,3,bcTypes[options->bcType],&bc,NULL);
456:   options->bcType = (BCType) bc;
457:   coeff = options->variableCoefficient;
458:   PetscOptionsEList("-variable_coefficient","Type of variable coefficent","ex12.c",coeffTypes,6,coeffTypes[options->variableCoefficient],&coeff,NULL);
459:   options->variableCoefficient = (CoeffType) coeff;

461:   PetscOptionsBool("-field_bc", "Use a field representation for the BC", "ex12.c", options->fieldBC, &options->fieldBC, NULL);
462:   PetscOptionsBool("-jacobian_mf", "Calculate the action of the Jacobian on the fly", "ex12.c", options->jacobianMF, &options->jacobianMF, NULL);
463:   PetscOptionsBool("-show_initial", "Output the initial guess for verification", "ex12.c", options->showInitial, &options->showInitial, NULL);
464:   PetscOptionsBool("-show_solution", "Output the solution for verification", "ex12.c", options->showSolution, &options->showSolution, NULL);
465:   PetscOptionsBool("-restart", "Read in the mesh and solution from a file", "ex12.c", options->restart, &options->restart, NULL);
466:   PetscOptionsBool("-check", "Compare with default integration routines", "ex12.c", options->check, &options->check, NULL);
467:   PetscOptionsBool("-dm_view_hierarchy", "View the coarsened hierarchy", "ex12.c", options->viewHierarchy, &options->viewHierarchy, NULL);
468:   PetscOptionsBool("-simplex", "Simplicial (true) or tensor (false) mesh", "ex12.c", options->simplex, &options->simplex, NULL);
469:   PetscOptionsBool("-quiet", "Don't print any vecs", "ex12.c", options->quiet, &options->quiet, NULL);
470:   PetscOptionsBool("-nonzero_initial_guess", "nonzero intial guess", "ex12.c", options->nonzInit, &options->nonzInit, NULL);
471:   PetscOptionsEnd();
472:   PetscLogEventRegister("CreateMesh", DM_CLASSID, &options->createMeshEvent);
473:   return(0);
474: }

476: static PetscErrorCode CreateBCLabel(DM dm, const char name[])
477: {
478:   DMLabel        label;

482:   DMCreateLabel(dm, name);
483:   DMGetLabel(dm, name, &label);
484:   DMPlexMarkBoundaryFaces(dm, 1, label);
485:   DMPlexLabelComplete(dm, label);
486:   return(0);
487: }

489: static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
490: {
491:   PetscInt       dim             = user->dim;
492:   const char    *filename        = user->filename;
493:   PetscBool      interpolate     = user->interpolate;
494:   PetscReal      refinementLimit = user->refinementLimit;
495:   size_t         len;

499:   PetscLogEventBegin(user->createMeshEvent,0,0,0,0);
500:   PetscStrlen(filename, &len);
501:   if (!len) {
502:     PetscInt d;

504:     if (user->periodicity[0] || user->periodicity[1] || user->periodicity[2]) for (d = 0; d < dim; ++d) user->cells[d] = PetscMax(user->cells[d], 3);
505:     DMPlexCreateBoxMesh(comm, dim, user->simplex, user->cells, NULL, NULL, user->periodicity, interpolate, dm);
506:     PetscObjectSetName((PetscObject) *dm, "Mesh");
507:   } else {
508:     DMPlexCreateFromFile(comm, filename, interpolate, dm);
509:     DMPlexSetRefinementUniform(*dm, PETSC_FALSE);
510:   }
511:   {
512:     PetscPartitioner part;
513:     DM               refinedMesh     = NULL;
514:     DM               distributedMesh = NULL;

516:     /* Refine mesh using a volume constraint */
517:     if (refinementLimit > 0.0) {
518:       DMPlexSetRefinementLimit(*dm, refinementLimit);
519:       DMRefine(*dm, comm, &refinedMesh);
520:       if (refinedMesh) {
521:         const char *name;

523:         PetscObjectGetName((PetscObject) *dm,         &name);
524:         PetscObjectSetName((PetscObject) refinedMesh,  name);
525:         DMDestroy(dm);
526:         *dm  = refinedMesh;
527:       }
528:     }
529:     /* Distribute mesh over processes */
530:     DMPlexGetPartitioner(*dm,&part);
531:     PetscPartitionerSetFromOptions(part);
532:     DMPlexDistribute(*dm, 0, NULL, &distributedMesh);
533:     if (distributedMesh) {
534:       DMDestroy(dm);
535:       *dm  = distributedMesh;
536:     }
537:   }
538:   if (user->bcType == NEUMANN) {
539:     DMLabel   label;

541:     DMCreateLabel(*dm, "boundary");
542:     DMGetLabel(*dm, "boundary", &label);
543:     DMPlexMarkBoundaryFaces(*dm, 1, label);
544:   } else if (user->bcType == DIRICHLET) {
545:     PetscBool hasLabel;

547:     DMHasLabel(*dm,"marker",&hasLabel);
548:     if (!hasLabel) {CreateBCLabel(*dm, "marker");}
549:   }
550:   {
551:     char      convType[256];
552:     PetscBool flg;

554:     PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX");
555:     PetscOptionsFList("-dm_plex_convert_type","Convert DMPlex to another format","ex12",DMList,DMPLEX,convType,256,&flg);
556:     PetscOptionsEnd();
557:     if (flg) {
558:       DM dmConv;

560:       DMConvert(*dm,convType,&dmConv);
561:       if (dmConv) {
562:         DMDestroy(dm);
563:         *dm  = dmConv;
564:       }
565:     }
566:   }
567:   DMLocalizeCoordinates(*dm); /* needed for periodic */
568:   DMSetFromOptions(*dm);
569:   DMViewFromOptions(*dm, NULL, "-dm_view");
570:   if (user->viewHierarchy) {
571:     DM       cdm = *dm;
572:     PetscInt i   = 0;
573:     char     buf[256];

575:     while (cdm) {
576:       DMSetUp(cdm);
577:       DMGetCoarseDM(cdm, &cdm);
578:       ++i;
579:     }
580:     cdm = *dm;
581:     while (cdm) {
582:       PetscViewer       viewer;
583:       PetscBool   isHDF5, isVTK;

585:       --i;
586:       PetscViewerCreate(comm,&viewer);
587:       PetscViewerSetType(viewer,PETSCVIEWERHDF5);
588:       PetscViewerSetOptionsPrefix(viewer,"hierarchy_");
589:       PetscViewerSetFromOptions(viewer);
590:       PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&isHDF5);
591:       PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isVTK);
592:       if (isHDF5) {
593:         PetscSNPrintf(buf, 256, "ex12-%d.h5", i);
594:       } else if (isVTK) {
595:         PetscSNPrintf(buf, 256, "ex12-%d.vtu", i);
596:         PetscViewerPushFormat(viewer,PETSC_VIEWER_VTK_VTU);
597:       } else {
598:         PetscSNPrintf(buf, 256, "ex12-%d", i);
599:       }
600:       PetscViewerFileSetMode(viewer,FILE_MODE_WRITE);
601:       PetscViewerFileSetName(viewer,buf);
602:       DMView(cdm, viewer);
603:       PetscViewerDestroy(&viewer);
604:       DMGetCoarseDM(cdm, &cdm);
605:     }
606:   }
607:   PetscLogEventEnd(user->createMeshEvent,0,0,0,0);
608:   return(0);
609: }

611: static PetscErrorCode SetupProblem(PetscDS prob, AppCtx *user)
612: {
613:   const PetscInt id = 1;

617:   switch (user->variableCoefficient) {
618:   case COEFF_NONE:
619:     if (user->periodicity[0]) {
620:       if (user->periodicity[1]) {
621:         PetscDSSetResidual(prob, 0, f0_xytrig_u, f1_u);
622:         PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu);
623:       } else {
624:         PetscDSSetResidual(prob, 0, f0_xtrig_u,  f1_u);
625:         PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu);
626:       }
627:     } else {
628:       PetscDSSetResidual(prob, 0, f0_u, f1_u);
629:       PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu);
630:     }
631:     break;
632:   case COEFF_ANALYTIC:
633:     PetscDSSetResidual(prob, 0, f0_analytic_u, f1_analytic_u);
634:     PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_analytic_uu);
635:     break;
636:   case COEFF_FIELD:
637:     PetscDSSetResidual(prob, 0, f0_analytic_u, f1_field_u);
638:     PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_field_uu);
639:     break;
640:   case COEFF_NONLINEAR:
641:     PetscDSSetResidual(prob, 0, f0_analytic_nonlinear_u, f1_analytic_nonlinear_u);
642:     PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_analytic_nonlinear_uu);
643:     break;
644:   case COEFF_CIRCLE:
645:     PetscDSSetResidual(prob, 0, f0_circle_u, f1_u);
646:     PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu);
647:     break;
648:   case COEFF_CROSS:
649:     PetscDSSetResidual(prob, 0, f0_cross_u, f1_u);
650:     PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu);
651:     break;
652:   default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid variable coefficient type %d", user->variableCoefficient);
653:   }
654:   switch (user->dim) {
655:   case 2:
656:     switch (user->variableCoefficient) {
657:     case COEFF_CIRCLE:
658:       user->exactFuncs[0]  = circle_u_2d;break;
659:     case COEFF_CROSS:
660:       user->exactFuncs[0]  = cross_u_2d;break;
661:     default:
662:       if (user->periodicity[0]) {
663:         if (user->periodicity[1]) {
664:           user->exactFuncs[0] = xytrig_u_2d;
665:         } else {
666:           user->exactFuncs[0] = xtrig_u_2d;
667:         }
668:       } else {
669:         user->exactFuncs[0]  = quadratic_u_2d;
670:         user->exactFields[0] = quadratic_u_field_2d;
671:       }
672:     }
673:     if (user->bcType == NEUMANN) {PetscDSSetBdResidual(prob, 0, f0_bd_u, f1_bd_zero);}
674:     break;
675:   case 3:
676:     user->exactFuncs[0]  = quadratic_u_3d;
677:     user->exactFields[0] = quadratic_u_field_3d;
678:     if (user->bcType == NEUMANN) {PetscDSSetBdResidual(prob, 0, f0_bd_u, f1_bd_zero);}
679:     break;
680:   default:
681:     SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %d", user->dim);
682:   }
683:   PetscDSAddBoundary(prob, user->bcType == DIRICHLET ? (user->fieldBC ? DM_BC_ESSENTIAL_FIELD : DM_BC_ESSENTIAL) : DM_BC_NATURAL,
684:                             "wall", user->bcType == DIRICHLET ? "marker" : "boundary", 0, 0, NULL,
685:                             user->fieldBC ? (void (*)()) user->exactFields[0] : (void (*)()) user->exactFuncs[0], 1, &id, user);
686:   PetscDSSetExactSolution(prob, 0, user->exactFuncs[0]);
687:   PetscDSSetFromOptions(prob);
688:   return(0);
689: }

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

698:   DMCreateLocalVector(dmAux, &nu);
699:   DMProjectFunctionLocal(dmAux, 0.0, matFuncs, NULL, INSERT_ALL_VALUES, nu);
700:   PetscObjectCompose((PetscObject) dm, "A", (PetscObject) nu);
701:   VecDestroy(&nu);
702:   return(0);
703: }

705: static PetscErrorCode SetupBC(DM dm, DM dmAux, AppCtx *user)
706: {
707:   PetscErrorCode (*bcFuncs[1])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx);
708:   Vec            uexact;
709:   PetscInt       dim;

713:   DMGetDimension(dm, &dim);
714:   if (dim == 2) bcFuncs[0] = quadratic_u_2d;
715:   else          bcFuncs[0] = quadratic_u_3d;
716:   DMCreateLocalVector(dmAux, &uexact);
717:   DMProjectFunctionLocal(dmAux, 0.0, bcFuncs, NULL, INSERT_ALL_VALUES, uexact);
718:   PetscObjectCompose((PetscObject) dm, "A", (PetscObject) uexact);
719:   VecDestroy(&uexact);
720:   return(0);
721: }

723: static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
724: {
725:   DM             cdm   = dm;
726:   const PetscInt dim   = user->dim;
727:   PetscFE        feAux = NULL;
728:   PetscFE        feCh  = NULL;
729:   PetscFE        fe;
730:   PetscDS        prob, probAux = NULL, probCh = NULL;
731:   PetscBool      simplex = user->simplex;

735:   /* Create finite element */
736:   PetscFECreateDefault(dm, dim, 1, simplex, NULL, -1, &fe);
737:   PetscObjectSetName((PetscObject) fe, "potential");
738:   if (user->variableCoefficient == COEFF_FIELD) {
739:     PetscQuadrature q;

741:     PetscFECreateDefault(dm, dim, 1, simplex, "mat_", -1, &feAux);
742:     PetscFEGetQuadrature(fe, &q);
743:     PetscFESetQuadrature(feAux, q);
744:     PetscDSCreate(PetscObjectComm((PetscObject)dm),&probAux);
745:     PetscDSSetDiscretization(probAux, 0, (PetscObject) feAux);
746:   } else if (user->fieldBC) {
747:     PetscQuadrature q;

749:     PetscFECreateDefault(dm, dim, 1, simplex, "bc_", -1, &feAux);
750:     PetscFEGetQuadrature(fe, &q);
751:     PetscFESetQuadrature(feAux, q);
752:     PetscDSCreate(PetscObjectComm((PetscObject)dm),&probAux);
753:     PetscDSSetDiscretization(probAux, 0, (PetscObject) feAux);
754:   }
755:   if (user->check) {
756:     PetscFECreateDefault(dm, dim, 1, simplex, "ch_", -1, &feCh);
757:     PetscDSCreate(PetscObjectComm((PetscObject)dm),&probCh);
758:     PetscDSSetDiscretization(probCh, 0, (PetscObject) feCh);
759:   }
760:   /* Set discretization and boundary conditions for each mesh */
761:   DMGetDS(dm, &prob);
762:   PetscDSSetDiscretization(prob, 0, (PetscObject) fe);
763:   SetupProblem(prob, user);
764:   while (cdm) {
765:     DM coordDM;

767:     DMSetDS(cdm,prob);
768:     DMGetCoordinateDM(cdm,&coordDM);
769:     if (feAux) {
770:       DM      dmAux;

772:       DMClone(cdm, &dmAux);
773:       DMSetCoordinateDM(dmAux, coordDM);
774:       DMSetDS(dmAux, probAux);
775:       PetscObjectCompose((PetscObject) dm, "dmAux", (PetscObject) dmAux);
776:       if (user->fieldBC) {SetupBC(cdm, dmAux, user);}
777:       else               {SetupMaterial(cdm, dmAux, user);}
778:       DMDestroy(&dmAux);
779:     }
780:     if (feCh) {
781:       DM      dmCh;

783:       DMClone(cdm, &dmCh);
784:       DMSetCoordinateDM(dmCh, coordDM);
785:       DMSetDS(dmCh, probCh);
786:       PetscObjectCompose((PetscObject) dm, "dmCh", (PetscObject) dmCh);
787:       DMDestroy(&dmCh);
788:     }
789:     if (user->bcType == DIRICHLET) {
790:       PetscBool hasLabel;

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

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

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

810:   Collective on KSP

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

818:   Level: intermediate

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

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

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

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

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

890:   Collective on SNES

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

898:   Level: intermediate

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

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

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

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

959:   PetscMalloc2(1, &user.exactFuncs, 1, &user.exactFields);
960:   SetupDiscretization(dm, &user);

962:   DMCreateGlobalVector(dm, &u);
963:   PetscObjectSetName((PetscObject) u, "potential");

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

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

979:     userJ.dm   = dm;
980:     userJ.J    = J;
981:     userJ.user = &user;

983:     DMCreateLocalVector(dm, &userJ.u);
984:     if (user.fieldBC) {DMProjectFieldLocal(dm, 0.0, userJ.u, user.exactFields, INSERT_BC_VALUES, userJ.u);}
985:     else              {DMProjectFunctionLocal(dm, 0.0, user.exactFuncs, NULL, INSERT_BC_VALUES, userJ.u);}
986:     MatShellSetContext(A, &userJ);
987:   } else {
988:     A = J;
989:   }
990:   if (user.bcType == NEUMANN) {
991:     MatNullSpaceCreate(PetscObjectComm((PetscObject) dm), PETSC_TRUE, 0, NULL, &nullSpace);
992:     MatSetNullSpace(A, nullSpace);
993:   }

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

998:   SNESSetFromOptions(snes);

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

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

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

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

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

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

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

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

1120:   if (user.bcType == NEUMANN) {MatNullSpaceDestroy(&nullSpace);}
1121:   if (user.jacobianMF) {VecDestroy(&userJ.u);}
1122:   if (A != J) {MatDestroy(&A);}
1123:   MatDestroy(&J);
1124:   VecDestroy(&u);
1125:   SNESDestroy(&snes);
1126:   DMDestroy(&dm);
1127:   PetscFree2(user.exactFuncs, user.exactFields);
1128:   PetscFinalize();
1129:   return ierr;
1130: }

1132: /*TEST
1133:   build:
1134:     requires: triangle
1135:   # 2D serial P1 test 0-4
1136:   test:
1137:     suffix: 0
1138:     requires: triangle
1139:     args: -run_type test -refinement_limit 0.0    -bc_type dirichlet -interpolate 0 -petscspace_order 1 -show_initial -dm_plex_print_fem 1

1141:   test:
1142:     suffix: 1
1143:     requires: triangle
1144:     args: -run_type test -refinement_limit 0.0    -bc_type dirichlet -interpolate 1 -petscspace_order 1 -show_initial -dm_plex_print_fem 1

1146:   test:
1147:     suffix: 2
1148:     requires: triangle
1149:     args: -run_type test -refinement_limit 0.0625 -bc_type dirichlet -interpolate 1 -petscspace_order 1 -show_initial -dm_plex_print_fem 1

1151:   test:
1152:     suffix: 3
1153:     requires: triangle
1154:     args: -run_type test -refinement_limit 0.0    -bc_type neumann   -interpolate 1 -petscspace_order 1 -show_initial -dm_plex_print_fem 1 -dm_view ascii::ascii_info_detail

1156:   test:
1157:     suffix: 4
1158:     requires: triangle
1159:     args: -run_type test -refinement_limit 0.0625 -bc_type neumann   -interpolate 1 -petscspace_order 1 -show_initial -dm_plex_print_fem 1

1161:   # 2D serial P2 test 5-8
1162:   test:
1163:     suffix: 5
1164:     requires: triangle
1165:     args: -run_type test -refinement_limit 0.0    -bc_type dirichlet -interpolate 1 -petscspace_order 2 -show_initial -dm_plex_print_fem 1

1167:   test:
1168:     suffix: 6
1169:     requires: triangle
1170:     args: -run_type test -refinement_limit 0.0625 -bc_type dirichlet -interpolate 1 -petscspace_order 2 -show_initial -dm_plex_print_fem 1

1172:   test:
1173:     suffix: 7
1174:     requires: triangle
1175:     args: -run_type test -refinement_limit 0.0    -bc_type neumann   -interpolate 1 -petscspace_order 2 -show_initial -dm_plex_print_fem 1 -dm_view ascii::ascii_info_detail

1177:   test:
1178:     suffix: 8
1179:     requires: triangle
1180:     args: -run_type test -refinement_limit 0.0625 -bc_type neumann   -interpolate 1 -petscspace_order 2 -show_initial -dm_plex_print_fem 1 -dm_view ascii::ascii_info_detail

1182:   # 3D serial P1 test 9-12
1183:   test:
1184:     suffix: 9
1185:     requires: ctetgen
1186:     args: -run_type test -dim 3 -refinement_limit 0.0    -bc_type dirichlet -interpolate 0 -petscspace_order 1 -show_initial -dm_plex_print_fem 1 -dm_view -cells 1,1,1

1188:   test:
1189:     suffix: 10
1190:     requires: ctetgen
1191:     args: -run_type test -dim 3 -refinement_limit 0.0    -bc_type dirichlet -interpolate 1 -petscspace_order 1 -show_initial -dm_plex_print_fem 1 -dm_view -cells 1,1,1

1193:   test:
1194:     suffix: 11
1195:     requires: ctetgen
1196:     args: -run_type test -dim 3 -refinement_limit 0.0125 -bc_type dirichlet -interpolate 1 -petscspace_order 1 -show_initial -dm_plex_print_fem 1 -dm_view -cells 1,1,1

1198:   test:
1199:     suffix: 12
1200:     requires: ctetgen
1201:     args: -run_type test -dim 3 -refinement_limit 0.0    -bc_type neumann   -interpolate 1 -petscspace_order 1 -snes_fd -show_initial -dm_plex_print_fem 1 -dm_view -cells 1,1,1

1203:   # Analytic variable coefficient 13-20
1204:   test:
1205:     suffix: 13
1206:     requires: triangle
1207:     args: -run_type test -refinement_limit 0.0    -variable_coefficient analytic -interpolate 1 -petscspace_order 1 -show_initial -dm_plex_print_fem 1
1208:   test:
1209:     suffix: 14
1210:     requires: triangle
1211:     args: -run_type test -refinement_limit 0.0625 -variable_coefficient analytic -interpolate 1 -petscspace_order 1 -show_initial -dm_plex_print_fem 1
1212:   test:
1213:     suffix: 15
1214:     requires: triangle
1215:     args: -run_type test -refinement_limit 0.0    -variable_coefficient analytic -interpolate 1 -petscspace_order 2 -show_initial -dm_plex_print_fem 1
1216:   test:
1217:     suffix: 16
1218:     requires: triangle
1219:     args: -run_type test -refinement_limit 0.0625 -variable_coefficient analytic -interpolate 1 -petscspace_order 2 -show_initial -dm_plex_print_fem 1
1220:   test:
1221:     suffix: 17
1222:     requires: hdf5 ctetgen
1223:     args: -run_type test -dim 3 -refinement_limit 0.0    -variable_coefficient analytic -interpolate 1 -petscspace_order 1 -show_initial -dm_plex_print_fem 1 -cells 1,1,1

1225:   test:
1226:     suffix: 18
1227:     requires: ctetgen
1228:     args: -run_type test -dim 3 -refinement_limit 0.0125 -variable_coefficient analytic -interpolate 1 -petscspace_order 1 -show_initial -dm_plex_print_fem 1 -cells 1,1,1

1230:   test:
1231:     suffix: 19
1232:     requires: ctetgen
1233:     args: -run_type test -dim 3 -refinement_limit 0.0    -variable_coefficient analytic -interpolate 1 -petscspace_order 2 -show_initial -dm_plex_print_fem 1 -cells 1,1,1

1235:   test:
1236:     suffix: 20
1237:     requires: ctetgen
1238:     args: -run_type test -dim 3 -refinement_limit 0.0125 -variable_coefficient analytic -interpolate 1 -petscspace_order 2 -show_initial -dm_plex_print_fem 1 -cells 1,1,1

1240:   # P1 variable coefficient 21-28
1241:   test:
1242:     suffix: 21
1243:     requires: triangle
1244:     args: -run_type test -refinement_limit 0.0    -variable_coefficient field    -interpolate 1 -petscspace_order 1 -mat_petscspace_order 1 -show_initial -dm_plex_print_fem 1

1246:   test:
1247:     suffix: 22
1248:     requires: triangle
1249:     args: -run_type test -refinement_limit 0.0625 -variable_coefficient field    -interpolate 1 -petscspace_order 1 -mat_petscspace_order 1 -show_initial -dm_plex_print_fem 1

1251:   test:
1252:     suffix: 23
1253:     requires: triangle
1254:     args: -run_type test -refinement_limit 0.0    -variable_coefficient field    -interpolate 1 -petscspace_order 2 -mat_petscspace_order 1 -show_initial -dm_plex_print_fem 1

1256:   test:
1257:     suffix: 24
1258:     requires: triangle
1259:     args: -run_type test -refinement_limit 0.0625 -variable_coefficient field    -interpolate 1 -petscspace_order 2 -mat_petscspace_order 1 -show_initial -dm_plex_print_fem 1

1261:   test:
1262:     suffix: 25
1263:     requires: ctetgen
1264:     args: -run_type test -dim 3 -refinement_limit 0.0    -variable_coefficient field    -interpolate 1 -petscspace_order 1 -mat_petscspace_order 1 -show_initial -dm_plex_print_fem 1 -cells 1,1,1

1266:   test:
1267:     suffix: 26
1268:     requires: ctetgen
1269:     args: -run_type test -dim 3 -refinement_limit 0.0125 -variable_coefficient field    -interpolate 1 -petscspace_order 1 -mat_petscspace_order 1 -show_initial -dm_plex_print_fem 1 -cells 1,1,1

1271:   test:
1272:     suffix: 27
1273:     requires: ctetgen
1274:     args: -run_type test -dim 3 -refinement_limit 0.0    -variable_coefficient field    -interpolate 1 -petscspace_order 2 -mat_petscspace_order 1 -show_initial -dm_plex_print_fem 1 -cells 1,1,1

1276:   test:
1277:     suffix: 28
1278:     requires: ctetgen
1279:     args: -run_type test -dim 3 -refinement_limit 0.0125 -variable_coefficient field    -interpolate 1 -petscspace_order 2 -mat_petscspace_order 1 -show_initial -dm_plex_print_fem 1 -cells 1,1,1

1281:   # P0 variable coefficient 29-36
1282:   test:
1283:     suffix: 29
1284:     requires: triangle
1285:     args: -run_type test -refinement_limit 0.0    -variable_coefficient field    -interpolate 1 -petscspace_order 1 -show_initial -dm_plex_print_fem 1

1287:   test:
1288:     suffix: 30
1289:     requires: triangle
1290:     args: -run_type test -refinement_limit 0.0625 -variable_coefficient field    -interpolate 1 -petscspace_order 1 -show_initial -dm_plex_print_fem 1

1292:   test:
1293:     suffix: 31
1294:     requires: triangle
1295:     args: -run_type test -refinement_limit 0.0    -variable_coefficient field    -interpolate 1 -petscspace_order 2 -show_initial -dm_plex_print_fem 1

1297:   test:
1298:     requires: triangle
1299:     suffix: 32
1300:     args: -run_type test -refinement_limit 0.0625 -variable_coefficient field    -interpolate 1 -petscspace_order 2 -show_initial -dm_plex_print_fem 1

1302:   test:
1303:     requires: ctetgen
1304:     suffix: 33
1305:     args: -run_type test -dim 3 -refinement_limit 0.0    -variable_coefficient field    -interpolate 1 -petscspace_order 1 -show_initial -dm_plex_print_fem 1 -cells 1,1,1

1307:   test:
1308:     suffix: 34
1309:     requires: ctetgen
1310:     args: -run_type test -dim 3 -refinement_limit 0.0125 -variable_coefficient field    -interpolate 1 -petscspace_order 1 -show_initial -dm_plex_print_fem 1 -cells 1,1,1

1312:   test:
1313:     suffix: 35
1314:     requires: ctetgen
1315:     args: -run_type test -dim 3 -refinement_limit 0.0    -variable_coefficient field    -interpolate 1 -petscspace_order 2 -show_initial -dm_plex_print_fem 1 -cells 1,1,1

1317:   test:
1318:     suffix: 36
1319:     requires: ctetgen
1320:     args: -run_type test -dim 3 -refinement_limit 0.0125 -variable_coefficient field    -interpolate 1 -petscspace_order 2 -show_initial -dm_plex_print_fem 1 -cells 1,1,1

1322:   # Full solve 39-44
1323:   test:
1324:     suffix: 39
1325:     requires: triangle !single
1326:     args: -run_type full -refinement_limit 0.015625 -interpolate 1 -petscspace_order 2 -pc_type gamg -ksp_rtol 1.0e-10 -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason ::ascii_info_detail
1327:   test:
1328:     suffix: 40
1329:     requires: triangle !single
1330:     args: -run_type full -refinement_limit 0.015625 -variable_coefficient nonlinear -interpolate 1 -petscspace_order 2 -pc_type svd -ksp_rtol 1.0e-10 -snes_monitor_short -snes_converged_reason ::ascii_info_detail
1331:   test:
1332:     suffix: 41
1333:     requires: triangle !single
1334:     args: -run_type full -refinement_limit 0.03125 -variable_coefficient nonlinear -interpolate 1 -petscspace_order 1 -snes_type fas -snes_fas_levels 2 -pc_type svd -ksp_rtol 1.0e-10 -fas_coarse_pc_type svd -fas_coarse_ksp_rtol 1.0e-10 -fas_coarse_snes_monitor_short -snes_monitor_short -snes_linesearch_type basic -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
1335:   test:
1336:     suffix: 42
1337:     requires: triangle !single
1338:     args: -run_type full -refinement_limit 0.0625 -variable_coefficient nonlinear -interpolate 1 -petscspace_order 1 -snes_type fas -snes_fas_levels 3 -pc_type svd -ksp_rtol 1.0e-10 -fas_coarse_pc_type svd -fas_coarse_ksp_rtol 1.0e-10 -fas_coarse_snes_monitor_short -snes_monitor_short -snes_linesearch_type basic -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
1339:   test:
1340:     suffix: 43
1341:     requires: triangle !single
1342:     nsize: 2
1343:     args: -run_type full -refinement_limit 0.03125 -variable_coefficient nonlinear -interpolate 1 -petscspace_order 1 -snes_type fas -snes_fas_levels 2 -pc_type svd -ksp_rtol 1.0e-10 -fas_coarse_pc_type svd -fas_coarse_ksp_rtol 1.0e-10 -fas_coarse_snes_monitor_short -snes_monitor_short -snes_linesearch_type basic -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

1345:   test:
1346:     suffix: 44
1347:     requires: triangle !single
1348:     nsize: 2
1349:     args: -run_type full -refinement_limit 0.0625 -variable_coefficient nonlinear -interpolate 1 -petscspace_order 1 -snes_type fas -snes_fas_levels 3 -pc_type svd -ksp_rtol 1.0e-10 -fas_coarse_pc_type svd -fas_coarse_ksp_rtol 1.0e-10 -fas_coarse_snes_monitor_short -snes_monitor_short -snes_linesearch_type basic -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

1351:   # Restarting
1352:   testset:
1353:     suffix: restart
1354:     requires: hdf5 triangle !complex
1355:     args: -run_type test -refinement_limit 0.0    -bc_type dirichlet -interpolate 1 -petscspace_order 1
1356:     test:
1357:       args: -dm_view hdf5:sol.h5 -vec_view hdf5:sol.h5::append
1358:     test:
1359:       args: -f sol.h5 -restart

1361:   # Periodicity
1362:   test:
1363:     suffix: periodic_0
1364:     requires: triangle
1365:     args: -run_type full -refinement_limit 0.0    -bc_type dirichlet -interpolate 1 -petscspace_order 1 -snes_converged_reason ::ascii_info_detail

1367:   # 2D serial P1 test with field bc
1368:   test:
1369:     suffix: field_bc_p1_0
1370:     requires: triangle
1371:     args: -run_type test              -interpolate 1 -bc_type dirichlet -field_bc -petscspace_order 1 -bc_petscspace_order 2 -show_initial -dm_plex_print_fem 1

1373:   test:
1374:     suffix: field_bc_p1_1
1375:     requires: triangle
1376:     args: -run_type test -dm_refine 1 -interpolate 1 -bc_type dirichlet -field_bc -petscspace_order 1 -bc_petscspace_order 2 -show_initial -dm_plex_print_fem 1

1378:   test:
1379:     suffix: field_bc_p1_2
1380:     requires: triangle
1381:     args: -run_type test              -interpolate 1 -bc_type neumann   -field_bc -petscspace_order 1 -bc_petscspace_order 2 -show_initial -dm_plex_print_fem 1

1383:   test:
1384:     suffix: field_bc_p1_3
1385:     requires: triangle
1386:     args: -run_type test -dm_refine 1 -interpolate 1 -bc_type neumann   -field_bc -petscspace_order 1 -bc_petscspace_order 2 -show_initial -dm_plex_print_fem 1

1388:   # 3D serial P1 test with field bc
1389:   test:
1390:     suffix: field_bc_p1_4
1391:     requires: ctetgen
1392:     args: -run_type test -dim 3              -interpolate 1 -bc_type dirichlet -field_bc -petscspace_order 1 -bc_petscspace_order 2 -show_initial -dm_plex_print_fem 1 -cells 1,1,1

1394:   test:
1395:     suffix: field_bc_p1_5
1396:     requires: ctetgen
1397:     args: -run_type test -dim 3 -dm_refine 1 -interpolate 1 -bc_type dirichlet -field_bc -petscspace_order 1 -bc_petscspace_order 2 -show_initial -dm_plex_print_fem 1 -cells 1,1,1

1399:   test:
1400:     suffix: field_bc_p1_6
1401:     requires: ctetgen
1402:     args: -run_type test -dim 3              -interpolate 1 -bc_type neumann   -field_bc -petscspace_order 1 -bc_petscspace_order 2 -show_initial -dm_plex_print_fem 1 -cells 1,1,1

1404:   test:
1405:     suffix: field_bc_p1_7
1406:     requires: ctetgen
1407:     args: -run_type test -dim 3 -dm_refine 1 -interpolate 1 -bc_type neumann   -field_bc -petscspace_order 1 -bc_petscspace_order 2 -show_initial -dm_plex_print_fem 1 -cells 1,1,1

1409:   # 2D serial P2 test with field bc
1410:   test:
1411:     suffix: field_bc_p2_0
1412:     requires: triangle
1413:     args: -run_type test              -interpolate 1 -bc_type dirichlet -field_bc -petscspace_order 2 -bc_petscspace_order 2 -show_initial -dm_plex_print_fem 1

1415:   test:
1416:     suffix: field_bc_p2_1
1417:     requires: triangle
1418:     args: -run_type test -dm_refine 1 -interpolate 1 -bc_type dirichlet -field_bc -petscspace_order 2 -bc_petscspace_order 2 -show_initial -dm_plex_print_fem 1

1420:   test:
1421:     suffix: field_bc_p2_2
1422:     requires: triangle
1423:     args: -run_type test              -interpolate 1 -bc_type neumann   -field_bc -petscspace_order 2 -bc_petscspace_order 2 -show_initial -dm_plex_print_fem 1

1425:   test:
1426:     suffix: field_bc_p2_3
1427:     requires: triangle
1428:     args: -run_type test -dm_refine 1 -interpolate 1 -bc_type neumann   -field_bc -petscspace_order 2 -bc_petscspace_order 2 -show_initial -dm_plex_print_fem 1

1430:   # 3D serial P2 test with field bc
1431:   test:
1432:     suffix: field_bc_p2_4
1433:     requires: ctetgen
1434:     args: -run_type test -dim 3              -interpolate 1 -bc_type dirichlet -field_bc -petscspace_order 2 -bc_petscspace_order 2 -show_initial -dm_plex_print_fem 1 -cells 1,1,1

1436:   test:
1437:     suffix: field_bc_p2_5
1438:     requires: ctetgen
1439:     args: -run_type test -dim 3 -dm_refine 1 -interpolate 1 -bc_type dirichlet -field_bc -petscspace_order 2 -bc_petscspace_order 2 -show_initial -dm_plex_print_fem 1 -cells 1,1,1

1441:   test:
1442:     suffix: field_bc_p2_6
1443:     requires: ctetgen
1444:     args: -run_type test -dim 3              -interpolate 1 -bc_type neumann   -field_bc -petscspace_order 2 -bc_petscspace_order 2 -show_initial -dm_plex_print_fem 1 -cells 1,1,1

1446:   test:
1447:     suffix: field_bc_p2_7
1448:     requires: ctetgen
1449:     args: -run_type test -dim 3 -dm_refine 1 -interpolate 1 -bc_type neumann   -field_bc -petscspace_order 2 -bc_petscspace_order 2 -show_initial -dm_plex_print_fem 1 -cells 1,1,1

1451:   # Full solve simplex: Convergence
1452:   test:
1453:     suffix: tet_conv_p1_r0
1454:     requires: ctetgen
1455:     args: -run_type full -dim 3 -dm_refine 0 -bc_type dirichlet -interpolate 1 -petscspace_order 1 -dm_view -snes_converged_reason ::ascii_info_detail -pc_type lu -cells 1,1,1
1456:   test:
1457:     suffix: tet_conv_p1_r2
1458:     requires: ctetgen
1459:     args: -run_type full -dim 3 -dm_refine 2 -bc_type dirichlet -interpolate 1 -petscspace_order 1 -dm_view -snes_converged_reason ::ascii_info_detail -pc_type lu -cells 1,1,1
1460:   test:
1461:     suffix: tet_conv_p1_r3
1462:     requires: ctetgen
1463:     args: -run_type full -dim 3 -dm_refine 3 -bc_type dirichlet -interpolate 1 -petscspace_order 1 -dm_view -snes_converged_reason ::ascii_info_detail -pc_type lu -cells 1,1,1
1464:   test:
1465:     suffix: tet_conv_p2_r0
1466:     requires: ctetgen
1467:     args: -run_type full -dim 3 -dm_refine 0 -bc_type dirichlet -interpolate 1 -petscspace_order 2 -dm_view -snes_converged_reason ::ascii_info_detail -pc_type lu -cells 1,1,1
1468:   test:
1469:     suffix: tet_conv_p2_r2
1470:     requires: ctetgen
1471:     args: -run_type full -dim 3 -dm_refine 2 -bc_type dirichlet -interpolate 1 -petscspace_order 2 -dm_view -snes_converged_reason ::ascii_info_detail -pc_type lu -cells 1,1,1

1473:   # Full solve simplex: BDDC
1474:   test:
1475:     suffix: tri_bddc
1476:     requires: triangle !single
1477:     nsize: 5
1478:     args: -run_type full -petscpartitioner_type simple -dm_refine 3 -bc_type dirichlet -interpolate 1 -petscspace_order 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

1480:   # Full solve simplex: BDDC
1481:   test:
1482:     suffix: tri_bddc_parmetis
1483:     requires: hdf5 triangle !single
1484:     nsize: 4
1485:     args: -run_type full -petscpartitioner_type parmetis -dm_refine 3 -bc_type dirichlet -interpolate 1 -petscspace_order 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

1487:   test:
1488:     suffix: quad_bddc
1489:     args: -run_type full -petscpartitioner_type simple -dm_refine 2 -bc_type dirichlet -interpolate 1 -petscspace_order 2 -dm_mat_type is -pc_type bddc -ksp_type gmres -snes_monitor_short -ksp_monitor_short -snes_view -simplex 0 -petscspace_poly_tensor -pc_bddc_corner_selection -cells 3,3 -ksp_rtol 1.e-9 -pc_bddc_use_edges 0
1490:     nsize: 9

1492:   # Full solve simplex: ASM
1493:   test:
1494:     suffix: tri_q2q1_asm_lu
1495:     requires: triangle !single
1496:     args: -run_type full -dm_refine 3 -bc_type dirichlet -interpolate 1 -petscspace_order 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

1498:   test:
1499:     suffix: tri_q2q1_msm_lu
1500:     requires: triangle !single
1501:     args: -run_type full -dm_refine 3 -bc_type dirichlet -interpolate 1 -petscspace_order 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

1503:   test:
1504:     suffix: tri_q2q1_asm_sor
1505:     requires: triangle !single
1506:     args: -run_type full -dm_refine 3 -bc_type dirichlet -interpolate 1 -petscspace_order 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

1508:   test:
1509:     suffix: tri_q2q1_msm_sor
1510:     requires: triangle !single
1511:     args: -run_type full -dm_refine 3 -bc_type dirichlet -interpolate 1 -petscspace_order 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

1513:   # Full solve simplex: FAS
1514:   test:
1515:     suffix: fas_newton_0
1516:     requires: triangle !single
1517:     args: -run_type full -variable_coefficient nonlinear -interpolate 1 -petscspace_order 1 -snes_type fas -snes_fas_levels 2 -pc_type svd -ksp_rtol 1.0e-10 -fas_coarse_pc_type svd -fas_coarse_ksp_rtol 1.0e-10 -fas_coarse_snes_monitor_short -snes_monitor_short -snes_linesearch_type basic -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

1519:   test:
1520:     suffix: fas_newton_1
1521:     requires: triangle !single
1522:     args: -run_type full -dm_refine_hierarchy 3 -interpolate 1 -petscspace_order 1 -snes_type fas -snes_fas_levels 3 -ksp_rtol 1.0e-10 -fas_coarse_pc_type lu -fas_coarse_snes_monitor_short -snes_monitor_short -snes_linesearch_type basic -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

1524:   test:
1525:     suffix: fas_ngs_0
1526:     requires: triangle !single
1527:     args: -run_type full -variable_coefficient nonlinear -interpolate 1 -petscspace_order 1 -snes_type fas -snes_fas_levels 2 -pc_type svd -ksp_rtol 1.0e-10 -fas_coarse_pc_type svd -fas_coarse_ksp_rtol 1.0e-10 -fas_coarse_snes_monitor_short -snes_monitor_short -snes_linesearch_type basic -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

1529:   test:
1530:     suffix: fas_newton_coarse_0
1531:     requires: pragmatic triangle
1532:     TODO: broken
1533:     args: -run_type full -dm_refine 2 -dm_plex_hash_location -variable_coefficient nonlinear -interpolate 1 -petscspace_order 1 -snes_type fas -snes_fas_levels 2 -pc_type svd -ksp_rtol 1.0e-10 -fas_coarse_pc_type svd -fas_coarse_ksp_rtol 1.0e-10 -fas_coarse_snes_monitor_short -snes_monitor_short -snes_linesearch_type basic -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -dm_coarsen_hierarchy 1 -snes_view -fas_levels_1_snes_type newtonls -fas_levels_1_pc_type svd -fas_levels_1_ksp_rtol 1.0e-10 -fas_levels_1_snes_monitor_short

1535:   test:
1536:     suffix: mg_newton_coarse_0
1537:     requires: triangle pragmatic
1538:     args: -run_type full -dm_refine 3 -interpolate 1 -petscspace_order 1 -snes_monitor_short -ksp_monitor_true_residual -snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -dm_coarsen_hierarchy 3 -dm_plex_hash_location -snes_view -dm_view -ksp_type richardson -pc_type mg  -pc_mg_levels 4 -snes_atol 1.0e-8 -ksp_atol 1.0e-8 -snes_rtol 0.0 -ksp_rtol 0.0 -ksp_norm_type unpreconditioned -mg_levels_ksp_type gmres -mg_levels_pc_type ilu -mg_levels_ksp_max_it 10

1540:   test:
1541:     suffix: mg_newton_coarse_1
1542:     requires: triangle pragmatic
1543:     args: -run_type full -dm_refine 5 -interpolate 1 -petscspace_order 1 -dm_coarsen_hierarchy 5 -dm_plex_hash_location -dm_plex_separate_marker -dm_plex_coarsen_bd_label marker -dm_plex_remesh_bd -ksp_type richardson -ksp_rtol 1.0e-12 -pc_type mg -pc_mg_levels 3 -mg_levels_ksp_max_it 2 -snes_converged_reason ::ascii_info_detail -snes_monitor -ksp_monitor_true_residual -mg_levels_ksp_monitor_true_residual -dm_view -ksp_view

1545:   test:
1546:     suffix: mg_newton_coarse_2
1547:     requires: triangle pragmatic
1548:     args: -run_type full -dm_refine 5 -interpolate 1 -petscspace_order 1 -dm_coarsen_hierarchy 5 -dm_plex_hash_location -dm_plex_separate_marker -dm_plex_remesh_bd -ksp_type richardson -ksp_rtol 1.0e-12 -pc_type mg -pc_mg_levels 3 -mg_levels_ksp_max_it 2 -snes_converged_reason ::ascii_info_detail -snes_monitor -ksp_monitor_true_residual -mg_levels_ksp_monitor_true_residual -dm_view -ksp_view

1550:   # Full solve tensor
1551:   test:
1552:     suffix: tensor_plex_2d
1553:     args: -run_type test -refinement_limit 0.0 -simplex 0 -interpolate -bc_type dirichlet -petscspace_order 1 -dm_refine_hierarchy 2 -cells 2,2

1555:   test:
1556:     suffix: tensor_p4est_2d
1557:     requires: p4est
1558:     args: -run_type test -refinement_limit 0.0 -simplex 0 -interpolate -bc_type dirichlet -petscspace_order 1 -dm_forest_initial_refinement 2 -dm_forest_minimum_refinement 0 -dm_plex_convert_type p4est -cells 2,2

1560:   test:
1561:     suffix: tensor_plex_3d
1562:     args: -run_type test -refinement_limit 0.0 -simplex 0 -interpolate -bc_type dirichlet -petscspace_order 1 -dim 3 -dm_refine_hierarchy 1 -cells 2,2,2

1564:   test:
1565:     suffix: tensor_p4est_3d
1566:     requires: p4est
1567:     args: -run_type test -refinement_limit 0.0 -simplex 0 -interpolate -bc_type dirichlet -petscspace_order 1 -dm_forest_initial_refinement 1 -dm_forest_minimum_refinement 0 -dim 3 -dm_plex_convert_type p8est -cells 2,2,2

1569:   test:
1570:     suffix: p4est_test_q2_conformal_serial
1571:     requires: p4est
1572:     args: -run_type test -interpolate 1 -petscspace_order 2 -simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -cells 2,2

1574:   test:
1575:     suffix: p4est_test_q2_conformal_parallel
1576:     requires: p4est
1577:     nsize: 7
1578:     args: -run_type test -interpolate 1 -petscspace_order 2 -simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -petscpartitioner_type simple -cells 2,2

1580:   test:
1581:     suffix: p4est_test_q2_conformal_parallel_parmetis
1582:     requires: hdf5 p4est
1583:     nsize: 4
1584:     args: -run_type test -interpolate 1 -petscspace_order 2 -simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -petscpartitioner_type parmetis -cells 2,2

1586:   test:
1587:     suffix: p4est_test_q2_nonconformal_serial
1588:     requires: p4est
1589:     filter: grep -v "CG or CGNE: variant"
1590:     args: -run_type test -interpolate 1 -petscspace_order 2 -simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -cells 2,2

1592:   test:
1593:     suffix: p4est_test_q2_nonconformal_parallel
1594:     requires: p4est
1595:     filter: grep -v "CG or CGNE: variant"
1596:     nsize: 7
1597:     args: -run_type test -interpolate 1 -petscspace_order 2 -simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type simple -cells 2,2

1599:   test:
1600:     suffix: p4est_test_q2_nonconformal_parallel_parmetis
1601:     requires: hdf5 p4est
1602:     nsize: 4
1603:     args: -run_type test -interpolate 1 -petscspace_order 2 -simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type parmetis -cells 2,2

1605:   test:
1606:     suffix: p4est_exact_q2_conformal_serial
1607:     requires: p4est !single
1608:     args: -run_type exact -interpolate 1 -petscspace_order 2 -snes_max_it 1 -snes_type fas -snes_fas_levels 3 -pc_type none -ksp_type preonly -fas_coarse_pc_type none -fas_coarse_ksp_type preonly -fas_coarse_snes_monitor_short -snes_monitor_short -snes_linesearch_type basic -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type none -fas_levels_ksp_type preonly -fas_levels_snes_monitor_short -simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -cells 2,2

1610:   test:
1611:     suffix: p4est_exact_q2_conformal_parallel
1612:     requires: p4est !single
1613:     nsize: 7
1614:     args: -run_type exact -interpolate 1 -petscspace_order 2 -snes_max_it 1 -snes_type fas -snes_fas_levels 3 -pc_type none -ksp_type preonly -fas_coarse_pc_type none -fas_coarse_ksp_type preonly -fas_coarse_snes_monitor_short -snes_monitor_short -snes_linesearch_type basic -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type none -fas_levels_ksp_type preonly -fas_levels_snes_monitor_short -simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -cells 2,2

1616:   test:
1617:     suffix: p4est_exact_q2_conformal_parallel_parmetis
1618:     requires: hdf5 p4est !single
1619:     nsize: 4
1620:     args: -run_type exact -interpolate 1 -petscspace_order 2 -snes_max_it 1 -snes_type fas -snes_fas_levels 3 -pc_type none -ksp_type preonly -fas_coarse_pc_type none -fas_coarse_ksp_type preonly -fas_coarse_snes_monitor_short -snes_monitor_short -snes_linesearch_type basic -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type none -fas_levels_ksp_type preonly -fas_levels_snes_monitor_short -simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -petscpartitioner_type parmetis  -cells 2,2

1622:   test:
1623:     suffix: p4est_exact_q2_nonconformal_serial
1624:     requires: p4est
1625:     args: -run_type exact -interpolate 1 -petscspace_order 2 -snes_max_it 1 -snes_type fas -snes_fas_levels 3 -pc_type none -ksp_type preonly -fas_coarse_pc_type none -fas_coarse_ksp_type preonly -fas_coarse_snes_monitor_short -snes_monitor_short -snes_linesearch_type basic -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type none -fas_levels_ksp_type preonly -fas_levels_snes_monitor_short -simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -cells 2,2

1627:   test:
1628:     suffix: p4est_exact_q2_nonconformal_parallel
1629:     requires: p4est
1630:     nsize: 7
1631:     args: -run_type exact -interpolate 1 -petscspace_order 2 -snes_max_it 1 -snes_type fas -snes_fas_levels 3 -pc_type none -ksp_type preonly -fas_coarse_pc_type none -fas_coarse_ksp_type preonly -fas_coarse_snes_monitor_short -snes_monitor_short -snes_linesearch_type basic -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type none -fas_levels_ksp_type preonly -fas_levels_snes_monitor_short -simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type simple -cells 2,2

1633:   test:
1634:     suffix: p4est_exact_q2_nonconformal_parallel_parmetis
1635:     requires: hdf5 p4est
1636:     nsize: 4
1637:     args: -run_type exact -interpolate 1 -petscspace_order 2 -snes_max_it 1 -snes_type fas -snes_fas_levels 3 -pc_type none -ksp_type preonly -fas_coarse_pc_type none -fas_coarse_ksp_type preonly -fas_coarse_snes_monitor_short -snes_monitor_short -snes_linesearch_type basic -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type none -fas_levels_ksp_type preonly -fas_levels_snes_monitor_short -simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type parmetis -cells 2,2

1639:   test:
1640:     suffix: p4est_full_q2_nonconformal_serial
1641:     requires: p4est !single
1642:     filter: grep -v "variant HERMITIAN"
1643:     args: -run_type full -interpolate 1 -petscspace_order 2 -snes_max_it 20 -snes_type fas -snes_fas_levels 3 -pc_type jacobi -ksp_type cg -fas_coarse_pc_type jacobi -fas_coarse_ksp_type cg -fas_coarse_snes_monitor_short -snes_monitor_short -snes_linesearch_type basic -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type jacobi -fas_levels_ksp_type cg -fas_levels_snes_monitor_short -simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -cells 2,2

1645:   test:
1646:     suffix: p4est_full_q2_nonconformal_parallel
1647:     requires: p4est !single
1648:     filter: grep -v "variant HERMITIAN"
1649:     nsize: 7
1650:     args: -run_type full -interpolate 1 -petscspace_order 2 -snes_max_it 20 -snes_type fas -snes_fas_levels 3 -pc_type jacobi -ksp_type cg -fas_coarse_pc_type jacobi -fas_coarse_ksp_type cg -fas_coarse_snes_monitor_short -snes_monitor_short -snes_linesearch_type basic -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type jacobi -fas_levels_ksp_type cg -fas_levels_snes_monitor_short -simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type simple -cells 2,2

1652:   test:
1653:     suffix: p4est_full_q2_nonconformal_parallel_bddcfas
1654:     requires: p4est
1655:     filter: grep -v "variant HERMITIAN"
1656:     nsize: 7
1657:     args: -run_type full -interpolate 1 -petscspace_order 2 -snes_max_it 20 -snes_type fas -snes_fas_levels 3 -dm_mat_type is -pc_type bddc -ksp_type cg -fas_coarse_pc_type bddc -fas_coarse_ksp_type cg -fas_coarse_snes_monitor_short -snes_monitor_short -snes_linesearch_type basic -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type bddc -fas_levels_ksp_type cg -fas_levels_snes_monitor_short -simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type simple -cells 2,2

1659:   test:
1660:     suffix: p4est_full_q2_nonconformal_parallel_bddc
1661:     requires: p4est
1662:     filter: grep -v "variant HERMITIAN"
1663:     nsize: 7
1664:     args: -run_type full -interpolate 1 -petscspace_order 2 -snes_max_it 20 -snes_type newtonls -dm_mat_type is -pc_type bddc -ksp_type cg -snes_monitor_short -snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type simple -cells 2,2

1666:   test:
1667:     suffix: p4est_fas_q2_conformal_serial
1668:     requires: p4est
1669:     args: -run_type full -variable_coefficient nonlinear -interpolate 1 -petscspace_order 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 -snes_linesearch_type basic -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type svd -fas_levels_ksp_type gmres -fas_levels_snes_monitor_short -simplex 0 -dm_refine_hierarchy 3 -cells 2,2
1670:     TODO: broken

1672:   test:
1673:     suffix: p4est_fas_q2_nonconformal_serial
1674:     requires: p4est broken
1675:     args: -run_type full -variable_coefficient nonlinear -interpolate 1 -petscspace_order 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 -snes_linesearch_type basic -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_snes_type newtonls -fas_levels_pc_type jacobi -fas_levels_ksp_type gmres -fas_levels_snes_monitor_short -simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -cells 2,2

1677:   test:
1678:     suffix: fas_newton_0_p4est
1679:     requires: p4est !single
1680:     args: -run_type full -variable_coefficient nonlinear -interpolate 1 -petscspace_order 1 -snes_type fas -snes_fas_levels 2 -pc_type svd -ksp_rtol 1.0e-10 -fas_coarse_pc_type svd -fas_coarse_ksp_rtol 1.0e-10 -fas_coarse_snes_monitor_short -snes_monitor_short -snes_linesearch_type basic -fas_coarse_snes_linesearch_type basic -snes_converged_reason ::ascii_info_detail -snes_view -fas_levels_1_snes_type newtonls -fas_levels_1_pc_type svd -fas_levels_1_ksp_rtol 1.0e-10 -fas_levels_1_snes_monitor_short -simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -cells 2,2

1682:   # Full solve simplicial AMR
1683:   test:
1684:     suffix: tri_p1_adapt_0
1685:     requires: pragmatic
1686:     args: -run_type exact -dim 2 -dm_refine 5 -bc_type dirichlet -interpolate 1 -petscspace_order 1 -variable_coefficient circle -snes_converged_reason ::ascii_info_detail -pc_type lu -adaptor_refinement_factor 1.0 -dm_view -dm_adapt_view -snes_adapt_initial

1688:   test:
1689:     suffix: tri_p1_adapt_1
1690:     requires: pragmatic
1691:     args: -run_type exact -dim 2 -dm_refine 5 -bc_type dirichlet -interpolate 1 -petscspace_order 1 -variable_coefficient circle -snes_converged_reason ::ascii_info_detail -pc_type lu -adaptor_refinement_factor 1.0 -dm_view -dm_adapt_iter_view -dm_adapt_view -snes_adapt_sequence 2

1693:   test:
1694:     suffix: tri_p1_adapt_analytic_0
1695:     requires: pragmatic
1696:     args: -run_type exact -dim 2 -dm_refine 3 -bc_type dirichlet -interpolate 1 -petscspace_order 1 -variable_coefficient cross -snes_adapt_initial 4 -adaptor_target_num 500 -adaptor_monitor -dm_view -dm_adapt_iter_view

1698:   # Full solve tensor AMR
1699:   test:
1700:     suffix: quad_q1_adapt_0
1701:     requires: p4est
1702:     args: -run_type exact -dim 2 -simplex 0 -dm_plex_convert_type p4est -bc_type dirichlet -interpolate 1 -petscspace_order 1 -variable_coefficient circle -snes_converged_reason ::ascii_info_detail -pc_type lu -dm_forest_initial_refinement 4   -snes_adapt_initial -dm_view
1703:     filter: grep -v DM_

1705:   test:
1706:     suffix: amr_0
1707:     nsize: 5
1708:     args: -run_type test -petscpartitioner_type simple -refinement_limit 0.0 -simplex 0 -interpolate -bc_type dirichlet -petscspace_order 1 -dm_refine 1 -cells 2,2

1710:   test:
1711:     suffix: amr_1
1712:     requires: p4est !complex
1713:     args: -run_type test -refinement_limit 0.0 -simplex 0 -interpolate -bc_type dirichlet -petscspace_order 1 -dm_plex_convert_type p4est -dm_p4est_refine_pattern center -dm_forest_maximum_refinement 5 -dm_view vtk:amr.vtu:vtk_vtu -vec_view vtk:amr.vtu:vtk_vtu:append -cells 2,2

1715:   test:
1716:     suffix: p4est_solve_bddc
1717:     requires: p4est
1718:     args: -run_type full -variable_coefficient nonlinear -nonzero_initial_guess 1 -interpolate 1 -petscspace_order 2 -snes_max_it 20 -snes_type newtonls -dm_mat_type is -pc_type bddc -ksp_type cg -snes_monitor_short -ksp_monitor -snes_linesearch_type bt -snes_converged_reason -snes_view -simplex 0 -petscspace_poly_tensor -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -petscpartitioner_type simple -pc_bddc_detect_disconnected
1719:     nsize: 4

1721:   test:
1722:     suffix: p4est_solve_fas
1723:     requires: p4est
1724:     args: -run_type full -variable_coefficient nonlinear -nonzero_initial_guess 1 -interpolate 1 -petscspace_order 2 -snes_max_it 10 -snes_type fas -snes_linesearch_type bt -snes_fas_levels 3 -fas_coarse_snes_type newtonls -fas_coarse_snes_linesearch_type basic -fas_coarse_ksp_type cg -fas_coarse_pc_type jacobi -fas_coarse_snes_monitor_short -fas_levels_snes_max_it 4 -fas_levels_snes_type newtonls -fas_levels_snes_linesearch_type bt -fas_levels_ksp_type cg -fas_levels_pc_type jacobi -fas_levels_snes_monitor_short -fas_levels_cycle_snes_linesearch_type bt -snes_monitor_short -snes_converged_reason -snes_view -simplex 0 -petscspace_poly_tensor -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash
1725:     nsize: 4
1726:     TODO: identical machine two runs produce slightly different solver trackers

1728:   test:
1729:     suffix: p4est_convergence_test_1
1730:     requires: p4est
1731:     args:  -quiet -run_type test -interpolate 1 -petscspace_order 1 -simplex 0 -petscspace_poly_tensor -dm_plex_convert_type p4est -dm_forest_minimum_refinement 2 -dm_forest_initial_refinement 2 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash
1732:     nsize: 4

1734:   test:
1735:     suffix: p4est_convergence_test_2
1736:     requires: p4est
1737:     args: -quiet -run_type test -interpolate 1 -petscspace_order 1 -simplex 0 -petscspace_poly_tensor -dm_plex_convert_type p4est -dm_forest_minimum_refinement 3 -dm_forest_initial_refinement 3 -dm_forest_maximum_refinement 5 -dm_p4est_refine_pattern hash

1739:   test:
1740:     suffix: p4est_convergence_test_3
1741:     requires: p4est
1742:     args: -quiet -run_type test -interpolate 1 -petscspace_order 1 -simplex 0 -petscspace_poly_tensor -dm_plex_convert_type p4est -dm_forest_minimum_refinement 4 -dm_forest_initial_refinement 4 -dm_forest_maximum_refinement 6 -dm_p4est_refine_pattern hash

1744:   test:
1745:     suffix: p4est_convergence_test_4
1746:     requires: p4est
1747:     args: -quiet -run_type test -interpolate 1 -petscspace_order 1 -simplex 0 -petscspace_poly_tensor -dm_plex_convert_type p4est -dm_forest_minimum_refinement 5 -dm_forest_initial_refinement 5 -dm_forest_maximum_refinement 7 -dm_p4est_refine_pattern hash
1748:     timeoutfactor: 5

1750:   # Serial tests with GLVis visualization
1751:   test:
1752:     suffix: glvis_2d_tet_p1
1753:     args: -quiet -run_type test -interpolate 1 -bc_type dirichlet -petscspace_order 1 -vec_view glvis: -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh
1754:   test:
1755:     suffix: glvis_2d_tet_p2
1756:     args: -quiet -run_type test -interpolate 1 -bc_type dirichlet -petscspace_order 2 -vec_view glvis: -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh
1757:   test:
1758:     suffix: glvis_2d_hex_p1
1759:     args: -quiet -run_type test -interpolate 1 -bc_type dirichlet -petscspace_order 1 -vec_view glvis: -simplex 0 -dm_refine 1
1760:   test:
1761:     suffix: glvis_2d_hex_p2
1762:     args: -quiet -run_type test -interpolate 1 -bc_type dirichlet -petscspace_order 2 -vec_view glvis: -simplex 0 -dm_refine 1
1763:   test:
1764:     suffix: glvis_2d_hex_p2_p4est
1765:     requires: p4est
1766:     args: -quiet -run_type test -interpolate 1 -bc_type dirichlet -petscspace_order 2 -vec_view glvis: -simplex 0 -dm_plex_convert_type p4est -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 1 -dm_forest_maximum_refinement 4 -dm_p4est_refine_pattern hash -cells 2,2 -viewer_glvis_dm_plex_enable_ncmesh

1768: TEST*/