Actual source code: ex62.c

petsc-master 2017-10-15
Report Typos and Errors
  1: static char help[] = "Stokes Problem in 2d and 3d with simplicial finite elements.\n\
  2: We solve the Stokes problem in a rectangular\n\
  3: domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\n\n";

  5: /*
  6: The isoviscous Stokes problem, which we discretize using the finite
  7: element method on an unstructured mesh. The weak form equations are

  9:   < \nabla v, \nabla u + {\nabla u}^T > - < \nabla\cdot v, p > + < v, f > = 0
 10:   < q, \nabla\cdot u >                                                    = 0

 12: We start with homogeneous Dirichlet conditions. We will expand this as the set
 13: of test problems is developed.

 15: Discretization:

 17: We use PetscFE to generate a tabulation of the finite element basis functions
 18: at quadrature points. We can currently generate an arbitrary order Lagrange
 19: element.

 21: Field Data:

 23:   DMPLEX data is organized by point, and the closure operation just stacks up the
 24: data from each sieve point in the closure. Thus, for a P_2-P_1 Stokes element, we
 25: have

 27:   cl{e} = {f e_0 e_1 e_2 v_0 v_1 v_2}
 28:   x     = [u_{e_0} v_{e_0} u_{e_1} v_{e_1} u_{e_2} v_{e_2} u_{v_0} v_{v_0} p_{v_0} u_{v_1} v_{v_1} p_{v_1} u_{v_2} v_{v_2} p_{v_2}]

 30: The problem here is that we would like to loop over each field separately for
 31: integration. Therefore, the closure visitor in DMPlexVecGetClosure() reorders
 32: the data so that each field is contiguous

 34:   x'    = [u_{e_0} v_{e_0} u_{e_1} v_{e_1} u_{e_2} v_{e_2} u_{v_0} v_{v_0} u_{v_1} v_{v_1} u_{v_2} v_{v_2} p_{v_0} p_{v_1} p_{v_2}]

 36: Likewise, DMPlexVecSetClosure() takes data partitioned by field, and correctly
 37: puts it into the Sieve ordering.

 39: Next Steps:

 41: - Refine and show convergence of correct order automatically (use femTest.py)
 42: - Fix InitialGuess for arbitrary disc (means making dual application work again)
 43: - Redo slides from GUCASTutorial for this new example

 45: For tensor product meshes, see SNES ex67, ex72
 46: */

 48:  #include <petscdmplex.h>
 49:  #include <petscsnes.h>
 50:  #include <petscds.h>

 52: PetscInt spatialDim = 0;

 54: typedef enum {NEUMANN, DIRICHLET} BCType;
 55: typedef enum {RUN_FULL, RUN_TEST} RunType;

 57: typedef struct {
 58:   PetscInt      debug;             /* The debugging level */
 59:   RunType       runType;           /* Whether to run tests, or solve the full problem */
 60:   PetscLogEvent createMeshEvent;
 61:   PetscBool     showInitial, showSolution, showError;
 62:   /* Domain and mesh definition */
 63:   PetscInt      dim;               /* The topological mesh dimension */
 64:   PetscBool     interpolate;       /* Generate intermediate mesh elements */
 65:   PetscBool     simplex;           /* Use simplices or tensor product cells */
 66:   PetscReal     refinementLimit;   /* The largest allowable cell volume */
 67:   PetscBool     testPartition;     /* Use a fixed partitioning for testing */
 68:   /* Problem definition */
 69:   BCType        bcType;
 70:   PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
 71: } AppCtx;

 73: PetscErrorCode zero_scalar(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
 74: {
 75:   u[0] = 0.0;
 76:   return 0;
 77: }
 78: PetscErrorCode zero_vector(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
 79: {
 80:   PetscInt d;
 81:   for (d = 0; d < spatialDim; ++d) u[d] = 0.0;
 82:   return 0;
 83: }

 85: /*
 86:   In 2D we use exact solution:

 88:     u = x^2 + y^2
 89:     v = 2 x^2 - 2xy
 90:     p = x + y - 1
 91:     f_x = f_y = 3

 93:   so that

 95:     -\Delta u + \nabla p + f = <-4, -4> + <1, 1> + <3, 3> = 0
 96:     \nabla \cdot u           = 2x - 2x                    = 0
 97: */
 98: PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
 99: {
100:   u[0] = x[0]*x[0] + x[1]*x[1];
101:   u[1] = 2.0*x[0]*x[0] - 2.0*x[0]*x[1];
102:   return 0;
103: }

105: PetscErrorCode linear_p_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *p, void *ctx)
106: {
107:   *p = x[0] + x[1] - 1.0;
108:   return 0;
109: }
110: PetscErrorCode constant_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *p, void *ctx)
111: {
112:   *p = 1.0;
113:   return 0;
114: }

116: void f0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
117:           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
118:           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
119:           PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
120: {
121:   PetscInt c;
122:   for (c = 0; c < dim; ++c) f0[c] = 3.0;
123: }

125: /* gradU[comp*dim+d] = {u_x, u_y, v_x, v_y} or {u_x, u_y, u_z, v_x, v_y, v_z, w_x, w_y, w_z}
126:    u[Ncomp]          = {p} */
127: void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
128:           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
129:           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
130:           PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
131: {
132:   const PetscInt Ncomp = dim;
133:   PetscInt       comp, d;

135:   for (comp = 0; comp < Ncomp; ++comp) {
136:     for (d = 0; d < dim; ++d) {
137:       /* f1[comp*dim+d] = 0.5*(gradU[comp*dim+d] + gradU[d*dim+comp]); */
138:       f1[comp*dim+d] = u_x[comp*dim+d];
139:     }
140:     f1[comp*dim+comp] -= u[Ncomp];
141:   }
142: }

144: /* gradU[comp*dim+d] = {u_x, u_y, v_x, v_y} or {u_x, u_y, u_z, v_x, v_y, v_z, w_x, w_y, w_z} */
145: void f0_p(PetscInt dim, PetscInt Nf, PetscInt NfAux,
146:           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
147:           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
148:           PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
149: {
150:   PetscInt d;
151:   for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d*dim+d];
152: }

154: void f1_p(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[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
158: {
159:   PetscInt d;
160:   for (d = 0; d < dim; ++d) f1[d] = 0.0;
161: }

163: /* < q, \nabla\cdot u >
164:    NcompI = 1, NcompJ = dim */
165: void g1_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
166:            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
167:            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
168:            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
169: {
170:   PetscInt d;
171:   for (d = 0; d < dim; ++d) g1[d*dim+d] = 1.0; /* \frac{\partial\phi^{u_d}}{\partial x_d} */
172: }

174: /* -< \nabla\cdot v, p >
175:     NcompI = dim, NcompJ = 1 */
176: void g2_up(PetscInt dim, PetscInt Nf, PetscInt NfAux,
177:            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
178:            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
179:            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
180: {
181:   PetscInt d;
182:   for (d = 0; d < dim; ++d) g2[d*dim+d] = -1.0; /* \frac{\partial\psi^{u_d}}{\partial x_d} */
183: }

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

195:   for (compI = 0; compI < Ncomp; ++compI) {
196:     for (d = 0; d < dim; ++d) {
197:       g3[((compI*Ncomp+compI)*dim+d)*dim+d] = 1.0;
198:     }
199:   }
200: }

202: /*
203:   In 3D we use exact solution:

205:     u = x^2 + y^2
206:     v = y^2 + z^2
207:     w = x^2 + y^2 - 2(x+y)z
208:     p = x + y + z - 3/2
209:     f_x = f_y = f_z = 3

211:   so that

213:     -\Delta u + \nabla p + f = <-4, -4, -4> + <1, 1, 1> + <3, 3, 3> = 0
214:     \nabla \cdot u           = 2x + 2y - 2(x + y)                   = 0
215: */
216: PetscErrorCode quadratic_u_3d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
217: {
218:   u[0] = x[0]*x[0] + x[1]*x[1];
219:   u[1] = x[1]*x[1] + x[2]*x[2];
220:   u[2] = x[0]*x[0] + x[1]*x[1] - 2.0*(x[0] + x[1])*x[2];
221:   return 0;
222: }

224: PetscErrorCode linear_p_3d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *p, void *ctx)
225: {
226:   *p = x[0] + x[1] + x[2] - 1.5;
227:   return 0;
228: }

230: PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
231: {
232:   const char    *bcTypes[2]  = {"neumann", "dirichlet"};
233:   const char    *runTypes[2] = {"full", "test"};
234:   PetscInt       bc, run;

238:   options->debug           = 0;
239:   options->runType         = RUN_FULL;
240:   options->dim             = 2;
241:   options->interpolate     = PETSC_FALSE;
242:   options->simplex         = PETSC_TRUE;
243:   options->refinementLimit = 0.0;
244:   options->testPartition   = PETSC_FALSE;
245:   options->bcType          = DIRICHLET;
246:   options->showInitial     = PETSC_FALSE;
247:   options->showSolution    = PETSC_TRUE;
248:   options->showError       = PETSC_FALSE;

250:   PetscOptionsBegin(comm, "", "Stokes Problem Options", "DMPLEX");
251:   PetscOptionsInt("-debug", "The debugging level", "ex62.c", options->debug, &options->debug, NULL);
252:   run  = options->runType;
253:   PetscOptionsEList("-run_type", "The run type", "ex62.c", runTypes, 2, runTypes[options->runType], &run, NULL);

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

257:   PetscOptionsInt("-dim", "The topological mesh dimension", "ex62.c", options->dim, &options->dim, NULL);
258:   spatialDim = options->dim;
259:   PetscOptionsBool("-interpolate", "Generate intermediate mesh elements", "ex62.c", options->interpolate, &options->interpolate, NULL);
260:   PetscOptionsBool("-simplex", "Use simplices or tensor product cells", "ex62.c", options->simplex, &options->simplex, NULL);
261:   PetscOptionsReal("-refinement_limit", "The largest allowable cell volume", "ex62.c", options->refinementLimit, &options->refinementLimit, NULL);
262:   PetscOptionsBool("-test_partition", "Use a fixed partition for testing", "ex62.c", options->testPartition, &options->testPartition, NULL);
263:   bc   = options->bcType;
264:   PetscOptionsEList("-bc_type","Type of boundary condition","ex62.c",bcTypes,2,bcTypes[options->bcType],&bc,NULL);

266:   options->bcType = (BCType) bc;

268:   PetscOptionsBool("-show_initial", "Output the initial guess for verification", "ex62.c", options->showInitial, &options->showInitial, NULL);
269:   PetscOptionsBool("-show_solution", "Output the solution for verification", "ex62.c", options->showSolution, &options->showSolution, NULL);
270:   PetscOptionsBool("-show_error", "Output the error for verification", "ex62.c", options->showError, &options->showError, NULL);
271:   PetscOptionsEnd();

273:   PetscLogEventRegister("CreateMesh", DM_CLASSID, &options->createMeshEvent);
274:   return(0);
275: }

277: PetscErrorCode DMVecViewLocal(DM dm, Vec v)
278: {
279:   Vec            lv;

283:   DMGetLocalVector(dm, &lv);
284:   DMGlobalToLocalBegin(dm, v, INSERT_VALUES, lv);
285:   DMGlobalToLocalEnd(dm, v, INSERT_VALUES, lv);
286:   DMPrintLocalVec(dm, "Local function", 0.0, lv);
287:   DMRestoreLocalVector(dm, &lv);
288:   return(0);
289: }

291: PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
292: {
293:   PetscInt       dim             = user->dim;
294:   PetscBool      interpolate     = user->interpolate;
295:   PetscReal      refinementLimit = user->refinementLimit;
296:   const PetscInt cells[3]        = {3, 3, 3};

300:   PetscLogEventBegin(user->createMeshEvent,0,0,0,0);
301:   DMPlexCreateBoxMesh(comm, dim, user->simplex, user->simplex ? NULL : cells, NULL, NULL, NULL, interpolate, dm);
302:   {
303:     DM refinedMesh     = NULL;
304:     DM distributedMesh = NULL;

306:     /* Refine mesh using a volume constraint */
307:     DMPlexSetRefinementLimit(*dm, refinementLimit);
308:     if (user->simplex) {DMRefine(*dm, comm, &refinedMesh);}
309:     if (refinedMesh) {
310:       DMDestroy(dm);
311:       *dm  = refinedMesh;
312:     }
313:     /* Setup test partitioning */
314:     if (user->testPartition) {
315:       PetscInt         triSizes_n2[2]       = {4, 4};
316:       PetscInt         triPoints_n2[8]      = {3, 5, 6, 7, 0, 1, 2, 4};
317:       PetscInt         triSizes_n3[3]       = {2, 3, 3};
318:       PetscInt         triPoints_n3[8]      = {3, 5, 1, 6, 7, 0, 2, 4};
319:       PetscInt         triSizes_n5[5]       = {1, 2, 2, 1, 2};
320:       PetscInt         triPoints_n5[8]      = {3, 5, 6, 4, 7, 0, 1, 2};
321:       PetscInt         triSizes_ref_n2[2]   = {8, 8};
322:       PetscInt         triPoints_ref_n2[16] = {1, 5, 6, 7, 10, 11, 14, 15, 0, 2, 3, 4, 8, 9, 12, 13};
323:       PetscInt         triSizes_ref_n3[3]   = {5, 6, 5};
324:       PetscInt         triPoints_ref_n3[16] = {1, 7, 10, 14, 15, 2, 6, 8, 11, 12, 13, 0, 3, 4, 5, 9};
325:       PetscInt         triSizes_ref_n5[5]   = {3, 4, 3, 3, 3};
326:       PetscInt         triPoints_ref_n5[16] = {1, 7, 10, 2, 11, 13, 14, 5, 6, 15, 0, 8, 9, 3, 4, 12};
327:       const PetscInt  *sizes = NULL;
328:       const PetscInt  *points = NULL;
329:       PetscPartitioner part;
330:       PetscInt         cEnd;
331:       PetscMPIInt      rank, size;

333:       MPI_Comm_rank(comm, &rank);
334:       MPI_Comm_size(comm, &size);
335:       DMPlexGetHeightStratum(*dm, 0, NULL, &cEnd);
336:       if (!rank) {
337:         if (dim == 2 && user->simplex && size == 2 && cEnd == 8) {
338:            sizes = triSizes_n2; points = triPoints_n2;
339:         } else if (dim == 2 && user->simplex && size == 3 && cEnd == 8) {
340:           sizes = triSizes_n3; points = triPoints_n3;
341:         } else if (dim == 2 && user->simplex && size == 5 && cEnd == 8) {
342:           sizes = triSizes_n5; points = triPoints_n5;
343:         } else if (dim == 2 && user->simplex && size == 2 && cEnd == 16) {
344:            sizes = triSizes_ref_n2; points = triPoints_ref_n2;
345:         } else if (dim == 2 && user->simplex && size == 3 && cEnd == 16) {
346:           sizes = triSizes_ref_n3; points = triPoints_ref_n3;
347:         } else if (dim == 2 && user->simplex && size == 5 && cEnd == 16) {
348:           sizes = triSizes_ref_n5; points = triPoints_ref_n5;
349:         } else SETERRQ(comm, PETSC_ERR_ARG_WRONG, "No stored partition matching run parameters");
350:       }
351:       DMPlexGetPartitioner(*dm, &part);
352:       PetscPartitionerSetType(part, PETSCPARTITIONERSHELL);
353:       PetscPartitionerShellSetPartition(part, size, sizes, points);
354:     } else {
355:       PetscPartitioner part;

357:       DMPlexGetPartitioner(*dm, &part);
358:       PetscPartitionerSetFromOptions(part);
359:     }
360:     /* Distribute mesh over processes */
361:     DMPlexDistribute(*dm, 0, NULL, &distributedMesh);
362:     if (distributedMesh) {
363:       DMDestroy(dm);
364:       *dm  = distributedMesh;
365:     }
366:   }
367:   DMSetFromOptions(*dm);
368:   DMViewFromOptions(*dm, NULL, "-dm_view");
369:   PetscLogEventEnd(user->createMeshEvent,0,0,0,0);
370:   return(0);
371: }

373: PetscErrorCode SetupProblem(PetscDS prob, AppCtx *user)
374: {
375:   const PetscInt id = 1;

379:   PetscDSSetResidual(prob, 0, f0_u, f1_u);
380:   PetscDSSetResidual(prob, 1, f0_p, f1_p);
381:   PetscDSSetJacobian(prob, 0, 0, NULL, NULL,  NULL,  g3_uu);
382:   PetscDSSetJacobian(prob, 0, 1, NULL, NULL,  g2_up, NULL);
383:   PetscDSSetJacobian(prob, 1, 0, NULL, g1_pu, NULL,  NULL);
384:   switch (user->dim) {
385:   case 2:
386:     user->exactFuncs[0] = quadratic_u_2d;
387:     user->exactFuncs[1] = linear_p_2d;
388:     break;
389:   case 3:
390:     user->exactFuncs[0] = quadratic_u_3d;
391:     user->exactFuncs[1] = linear_p_3d;
392:     break;
393:   default:
394:     SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %d", user->dim);
395:   }
396:   PetscDSAddBoundary(prob, user->bcType == DIRICHLET ? DM_BC_ESSENTIAL : DM_BC_NATURAL, "wall", user->bcType == NEUMANN ? "boundary" : "marker", 0, 0, NULL, (void (*)(void)) user->exactFuncs[0], 1, &id, user);
397:   return(0);
398: }

400: PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
401: {
402:   DM              cdm   = dm;
403:   const PetscInt  dim   = user->dim;
404:   PetscFE         fe[2];
405:   PetscQuadrature q;
406:   PetscDS         prob;
407:   PetscErrorCode  ierr;

410:   /* Create finite element */
411:   PetscFECreateDefault(dm, dim, dim, user->simplex, "vel_", PETSC_DEFAULT, &fe[0]);
412:   PetscObjectSetName((PetscObject) fe[0], "velocity");
413:   PetscFEGetQuadrature(fe[0], &q);
414:   PetscFECreateDefault(dm, dim, 1, user->simplex, "pres_", PETSC_DEFAULT, &fe[1]);
415:   PetscFESetQuadrature(fe[1], q);
416:   PetscObjectSetName((PetscObject) fe[1], "pressure");
417:   /* Set discretization and boundary conditions for each mesh */
418:   DMGetDS(dm, &prob);
419:   PetscDSSetDiscretization(prob, 0, (PetscObject) fe[0]);
420:   PetscDSSetDiscretization(prob, 1, (PetscObject) fe[1]);
421:   SetupProblem(prob, user);
422:   while (cdm) {
423:     DMSetDS(cdm, prob);
424:     DMGetCoarseDM(cdm, &cdm);
425:   }
426:   PetscFEDestroy(&fe[0]);
427:   PetscFEDestroy(&fe[1]);
428:   return(0);
429: }

431: PetscErrorCode CreatePressureNullSpace(DM dm, AppCtx *user, Vec *v, MatNullSpace *nullSpace)
432: {
433:   Vec              vec;
434:   PetscErrorCode (*funcs[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void* ctx) = {zero_vector, constant_p};
435:   PetscErrorCode   ierr;

438:   DMGetGlobalVector(dm, &vec);
439:   DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec);
440:   VecNormalize(vec, NULL);
441:   if (user->debug) {
442:     PetscPrintf(PetscObjectComm((PetscObject)dm), "Pressure Null Space\n");
443:     VecView(vec, PETSC_VIEWER_STDOUT_WORLD);
444:   }
445:   MatNullSpaceCreate(PetscObjectComm((PetscObject)dm), PETSC_FALSE, 1, &vec, nullSpace);
446:   if (v) {
447:     DMCreateGlobalVector(dm, v);
448:     VecCopy(vec, *v);
449:   }
450:   DMRestoreGlobalVector(dm, &vec);
451:   /* New style for field null spaces */
452:   {
453:     PetscObject  pressure;
454:     MatNullSpace nullSpacePres;

456:     DMGetField(dm, 1, &pressure);
457:     MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nullSpacePres);
458:     PetscObjectCompose(pressure, "nullspace", (PetscObject) nullSpacePres);
459:     MatNullSpaceDestroy(&nullSpacePres);
460:   }
461:   return(0);
462: }

464: int main(int argc, char **argv)
465: {
466:   SNES           snes;                 /* nonlinear solver */
467:   DM             dm;                   /* problem definition */
468:   Vec            u,r;                  /* solution, residual vectors */
469:   Mat            A,J;                  /* Jacobian matrix */
470:   MatNullSpace   nullSpace;            /* May be necessary for pressure */
471:   AppCtx         user;                 /* user-defined work context */
472:   PetscInt       its;                  /* iterations for convergence */
473:   PetscReal      error         = 0.0;  /* L_2 error in the solution */
474:   PetscReal      ferrors[2];

477:   PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
478:   ProcessOptions(PETSC_COMM_WORLD, &user);
479:   SNESCreate(PETSC_COMM_WORLD, &snes);
480:   CreateMesh(PETSC_COMM_WORLD, &user, &dm);
481:   SNESSetDM(snes, dm);
482:   DMSetApplicationContext(dm, &user);

484:   PetscMalloc(2 * sizeof(void (*)(const PetscReal[], PetscScalar *, void *)), &user.exactFuncs);
485:   SetupDiscretization(dm, &user);
486:   DMPlexCreateClosureIndex(dm, NULL);

488:   DMCreateGlobalVector(dm, &u);
489:   VecDuplicate(u, &r);

491:   DMCreateMatrix(dm, &J);
492:   A = J;
493:   CreatePressureNullSpace(dm, &user, NULL, &nullSpace);
494:   MatSetNullSpace(A, nullSpace);

496:   DMPlexSetSNESLocalFEM(dm,&user,&user,&user);
497:   SNESSetJacobian(snes, A, J, NULL, NULL);

499:   SNESSetFromOptions(snes);

501:   DMProjectFunction(dm, 0.0, user.exactFuncs, NULL, INSERT_ALL_VALUES, u);
502:   if (user.showInitial) {DMVecViewLocal(dm, u);}
503:   if (user.runType == RUN_FULL) {
504:     PetscErrorCode (*initialGuess[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void* ctx) = {zero_vector, zero_scalar};

506:     DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u);
507:     if (user.debug) {
508:       PetscPrintf(PETSC_COMM_WORLD, "Initial guess\n");
509:       VecView(u, PETSC_VIEWER_STDOUT_WORLD);
510:     }
511:     SNESSolve(snes, NULL, u);
512:     SNESGetIterationNumber(snes, &its);
513:     PetscPrintf(PETSC_COMM_WORLD, "Number of SNES iterations = %D\n", its);
514:     DMComputeL2Diff(dm, 0.0, user.exactFuncs, NULL, u, &error);
515:     DMComputeL2FieldDiff(dm, 0.0, user.exactFuncs, NULL, u, ferrors);
516:     PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: %g [%g, %g]\n", error, ferrors[0], ferrors[1]);
517:     if (user.showError) {
518:       Vec r;
519:       DMGetGlobalVector(dm, &r);
520:       DMProjectFunction(dm, 0.0, user.exactFuncs, NULL, INSERT_ALL_VALUES, r);
521:       VecAXPY(r, -1.0, u);
522:       PetscPrintf(PETSC_COMM_WORLD, "Solution Error\n");
523:       VecView(r, PETSC_VIEWER_STDOUT_WORLD);
524:       DMRestoreGlobalVector(dm, &r);
525:     }
526:     if (user.showSolution) {
527:       PetscPrintf(PETSC_COMM_WORLD, "Solution\n");
528:       VecChop(u, 3.0e-9);
529:       VecView(u, PETSC_VIEWER_STDOUT_WORLD);
530:     }
531:   } else {
532:     PetscReal res = 0.0;

534:     /* Check discretization error */
535:     PetscPrintf(PETSC_COMM_WORLD, "Initial guess\n");
536:     VecView(u, PETSC_VIEWER_STDOUT_WORLD);
537:     DMComputeL2Diff(dm, 0.0, user.exactFuncs, NULL, u, &error);
538:     if (error >= 1.0e-11) {PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: %g\n", error);}
539:     else                  {PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: < 1.0e-11\n");}
540:     /* Check residual */
541:     SNESComputeFunction(snes, u, r);
542:     PetscPrintf(PETSC_COMM_WORLD, "Initial Residual\n");
543:     VecChop(r, 1.0e-10);
544:     VecView(r, PETSC_VIEWER_STDOUT_WORLD);
545:     VecNorm(r, NORM_2, &res);
546:     PetscPrintf(PETSC_COMM_WORLD, "L_2 Residual: %g\n", res);
547:     /* Check Jacobian */
548:     {
549:       Vec          b;
550:       PetscBool    isNull;

552:       SNESComputeJacobian(snes, u, A, A);
553:       MatNullSpaceTest(nullSpace, J, &isNull);
554:       VecDuplicate(u, &b);
555:       VecSet(r, 0.0);
556:       SNESComputeFunction(snes, r, b);
557:       MatMult(A, u, r);
558:       VecAXPY(r, 1.0, b);
559:       VecDestroy(&b);
560:       PetscPrintf(PETSC_COMM_WORLD, "Au - b = Au + F(0)\n");
561:       VecChop(r, 1.0e-10);
562:       VecView(r, PETSC_VIEWER_STDOUT_WORLD);
563:       VecNorm(r, NORM_2, &res);
564:       PetscPrintf(PETSC_COMM_WORLD, "Linear L_2 Residual: %g\n", res);
565:     }
566:   }
567:   VecViewFromOptions(u, NULL, "-sol_vec_view");

569:   MatNullSpaceDestroy(&nullSpace);
570:   if (A != J) {MatDestroy(&A);}
571:   MatDestroy(&J);
572:   VecDestroy(&u);
573:   VecDestroy(&r);
574:   SNESDestroy(&snes);
575:   DMDestroy(&dm);
576:   PetscFree(user.exactFuncs);
577:   PetscFinalize();
578:   return ierr;
579: }

581: /*TEST

583:   # 2D serial P1 tests 0-3
584:   test:
585:     suffix: 0
586:     requires: triangle
587:     args: -run_type test -refinement_limit 0.0    -bc_type dirichlet -interpolate 0 -vel_petscspace_order 1 -pres_petscspace_order 1 -show_initial -dm_plex_print_fem 1
588:   test:
589:     suffix: 1
590:     requires: triangle
591:     args: -run_type test -refinement_limit 0.0    -bc_type dirichlet -interpolate 1 -vel_petscspace_order 1 -pres_petscspace_order 1 -show_initial -dm_plex_print_fem 1
592:   test:
593:     suffix: 2
594:     requires: triangle
595:     args: -run_type test -refinement_limit 0.0625 -bc_type dirichlet -interpolate 0 -vel_petscspace_order 1 -pres_petscspace_order 1 -show_initial -dm_plex_print_fem 1
596:   test:
597:     suffix: 3
598:     requires: triangle
599:     args: -run_type test -refinement_limit 0.0625 -bc_type dirichlet -interpolate 1 -vel_petscspace_order 1 -pres_petscspace_order 1 -show_initial -dm_plex_print_fem 1
600:   # 2D serial P2 tests 4-5
601:   test:
602:     suffix: 4
603:     requires: triangle
604:     args: -run_type test -refinement_limit 0.0    -bc_type dirichlet -interpolate 1 -vel_petscspace_order 2 -pres_petscspace_order 1 -show_initial -dm_plex_print_fem 1
605:   test:
606:     suffix: 5
607:     requires: triangle
608:     args: -run_type test -refinement_limit 0.0625 -bc_type dirichlet -interpolate 1 -vel_petscspace_order 2 -pres_petscspace_order 1 -show_initial -dm_plex_print_fem 1
609:   # Parallel tests 6-17
610:   test:
611:     suffix: 6
612:     requires: triangle
613:     nsize: 2
614:     args: -run_type test -refinement_limit 0.0    -test_partition -bc_type dirichlet -interpolate 0 -vel_petscspace_order 1 -pres_petscspace_order 1 -dm_plex_print_fem 1
615:   test:
616:     suffix: 7
617:     requires: triangle
618:     nsize: 3
619:     args: -run_type test -refinement_limit 0.0    -test_partition -bc_type dirichlet -interpolate 0 -vel_petscspace_order 1 -pres_petscspace_order 1 -dm_plex_print_fem 1
620:   test:
621:     suffix: 8
622:     requires: triangle
623:     nsize: 5
624:     args: -run_type test -refinement_limit 0.0    -test_partition -bc_type dirichlet -interpolate 0 -vel_petscspace_order 1 -pres_petscspace_order 1 -dm_plex_print_fem 1
625:   test:
626:     suffix: 9
627:     requires: triangle
628:     nsize: 2
629:     args: -run_type test -refinement_limit 0.0    -test_partition -bc_type dirichlet -interpolate 1 -vel_petscspace_order 1 -pres_petscspace_order 1 -dm_plex_print_fem 1
630:   test:
631:     suffix: 10
632:     requires: triangle
633:     nsize: 3
634:     args: -run_type test -refinement_limit 0.0    -test_partition -bc_type dirichlet -interpolate 1 -vel_petscspace_order 1 -pres_petscspace_order 1 -dm_plex_print_fem 1
635:   test:
636:     suffix: 11
637:     requires: triangle
638:     nsize: 5
639:     args: -run_type test -refinement_limit 0.0    -test_partition -bc_type dirichlet -interpolate 1 -vel_petscspace_order 1 -pres_petscspace_order 1 -dm_plex_print_fem 1
640:   test:
641:     suffix: 12
642:     requires: triangle
643:     nsize: 2
644:     args: -run_type test -refinement_limit 0.0625 -test_partition -bc_type dirichlet -interpolate 0 -vel_petscspace_order 1 -pres_petscspace_order 1 -dm_plex_print_fem 1
645:   test:
646:     suffix: 13
647:     requires: triangle
648:     nsize: 3
649:     args: -run_type test -refinement_limit 0.0625 -test_partition -bc_type dirichlet -interpolate 0 -vel_petscspace_order 1 -pres_petscspace_order 1 -dm_plex_print_fem 1
650:   test:
651:     suffix: 14
652:     requires: triangle
653:     nsize: 5
654:     args: -run_type test -refinement_limit 0.0625 -test_partition -bc_type dirichlet -interpolate 0 -vel_petscspace_order 1 -pres_petscspace_order 1 -dm_plex_print_fem 1
655:   test:
656:     suffix: 15
657:     requires: triangle
658:     nsize: 2
659:     args: -run_type test -refinement_limit 0.0625 -test_partition -bc_type dirichlet -interpolate 1 -vel_petscspace_order 1 -pres_petscspace_order 1 -dm_plex_print_fem 1
660:   test:
661:     suffix: 16
662:     requires: triangle
663:     nsize: 3
664:     args: -run_type test -refinement_limit 0.0625 -test_partition -bc_type dirichlet -interpolate 1 -vel_petscspace_order 1 -pres_petscspace_order 1 -dm_plex_print_fem 1
665:   test:
666:     suffix: 17
667:     requires: triangle
668:     nsize: 5
669:     args: -run_type test -refinement_limit 0.0625 -test_partition -bc_type dirichlet -interpolate 1 -vel_petscspace_order 1 -pres_petscspace_order 1 -dm_plex_print_fem 1
670:   # 3D serial P1 tests 43-46
671:   test:
672:     suffix: 43
673:     requires: ctetgen
674:     args: -run_type test -dim 3 -refinement_limit 0.0    -bc_type dirichlet -interpolate 0 -vel_petscspace_order 1 -pres_petscspace_order 1 -show_initial -dm_plex_print_fem 1
675:   test:
676:     suffix: 44
677:     requires: ctetgen
678:     args: -run_type test -dim 3 -refinement_limit 0.0    -bc_type dirichlet -interpolate 1 -vel_petscspace_order 1 -pres_petscspace_order 1 -show_initial -dm_plex_print_fem 1
679:   test:
680:     suffix: 45
681:     requires: ctetgen
682:     args: -run_type test -dim 3 -refinement_limit 0.0125 -bc_type dirichlet -interpolate 0 -vel_petscspace_order 1 -pres_petscspace_order 1 -show_initial -dm_plex_print_fem 1
683:   test:
684:     suffix: 46
685:     requires: ctetgen
686:     args: -run_type test -dim 3 -refinement_limit 0.0125 -bc_type dirichlet -interpolate 1 -vel_petscspace_order 1 -pres_petscspace_order 1 -show_initial -dm_plex_print_fem 1
687:   # Full solutions 18-29
688:   test:
689:     suffix: 18
690:     requires: triangle
691:     filter:  sed -e "s/total number of linear solver iterations=11/total number of linear solver iterations=12/g"
692:     args: -run_type full -refinement_limit 0.0625 -bc_type dirichlet -interpolate 0 -vel_petscspace_order 1 -pres_petscspace_order 1 -pc_type jacobi -ksp_rtol 1.0e-9 -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view
693:   test:
694:     suffix: 19
695:     requires: triangle
696:     nsize: 2
697:     args: -run_type full -petscpartitioner_type simple -refinement_limit 0.0625 -bc_type dirichlet -interpolate 0 -vel_petscspace_order 1 -pres_petscspace_order 1 -pc_type jacobi -ksp_rtol 1.0e-9 -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view
698:   test:
699:     suffix: 20
700:     requires: triangle
701:     filter:  sed -e "s/total number of linear solver iterations=11/total number of linear solver iterations=12/g"
702:     nsize: 3
703:     args: -run_type full -petscpartitioner_type simple -refinement_limit 0.0625 -bc_type dirichlet -interpolate 0 -vel_petscspace_order 1 -pres_petscspace_order 1 -pc_type jacobi -ksp_rtol 1.0e-9 -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view
704:   test:
705:     suffix: 21
706:     requires: triangle
707:     filter:  sed -e "s/total number of linear solver iterations=11/total number of linear solver iterations=12/g"
708:     nsize: 5
709:     args: -run_type full -petscpartitioner_type simple -refinement_limit 0.0625 -bc_type dirichlet -interpolate 0 -vel_petscspace_order 1 -pres_petscspace_order 1 -pc_type jacobi -ksp_rtol 1.0e-9 -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view
710:   test:
711:     suffix: 22
712:     requires: triangle
713:     filter:  sed -e "s/total number of linear solver iterations=11/total number of linear solver iterations=12/g"
714:     args: -run_type full -refinement_limit 0.0625 -bc_type dirichlet -interpolate 1 -vel_petscspace_order 1 -pres_petscspace_order 1 -pc_type jacobi -ksp_rtol 1.0e-9 -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view
715:   test:
716:     suffix: 23
717:     requires: triangle
718:     nsize: 2
719:     args: -run_type full -petscpartitioner_type simple -refinement_limit 0.0625 -bc_type dirichlet -interpolate 1 -vel_petscspace_order 1 -pres_petscspace_order 1 -pc_type jacobi -ksp_rtol 1.0e-9 -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view
720:   test:
721:     suffix: 24
722:     requires: triangle
723:     nsize: 3
724:     filter:  sed -e "s/total number of linear solver iterations=11/total number of linear solver iterations=12/g"
725:     args: -run_type full -petscpartitioner_type simple -refinement_limit 0.0625 -bc_type dirichlet -interpolate 1 -vel_petscspace_order 1 -pres_petscspace_order 1 -pc_type jacobi -ksp_rtol 1.0e-9 -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view
726:   test:
727:     suffix: 25
728:     requires: triangle
729:     filter:  sed -e "s/total number of linear solver iterations=11/total number of linear solver iterations=12/g"
730:     nsize: 5
731:     args: -run_type full -petscpartitioner_type simple -refinement_limit 0.0625 -bc_type dirichlet -interpolate 1 -vel_petscspace_order 1 -pres_petscspace_order 1 -pc_type jacobi -ksp_rtol 1.0e-9 -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view
732:   test:
733:     suffix: 26
734:     requires: triangle
735:     args: -run_type full -refinement_limit 0.0625 -bc_type dirichlet -interpolate 1 -vel_petscspace_order 2 -pres_petscspace_order 1 -pc_type jacobi -ksp_rtol 1.0e-9 -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view
736:   test:
737:     suffix: 27
738:     requires: triangle
739:     nsize: 2
740:     args: -run_type full -petscpartitioner_type simple -refinement_limit 0.0625 -bc_type dirichlet -interpolate 1 -vel_petscspace_order 2 -pres_petscspace_order 1 -pc_type jacobi -ksp_rtol 1.0e-9 -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view
741:   test:
742:     suffix: 28
743:     requires: triangle
744:     nsize: 3
745:     args: -run_type full -petscpartitioner_type simple -refinement_limit 0.0625 -bc_type dirichlet -interpolate 1 -vel_petscspace_order 2 -pres_petscspace_order 1 -pc_type jacobi -ksp_rtol 1.0e-9 -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view
746:   test:
747:     suffix: 29
748:     requires: triangle
749:     nsize: 5
750:     args: -run_type full -petscpartitioner_type simple -refinement_limit 0.0625 -bc_type dirichlet -interpolate 1 -vel_petscspace_order 2 -pres_petscspace_order 1 -pc_type jacobi -ksp_rtol 1.0e-9 -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view
751:   # Full solutions with quads
752:   #   FULL Schur with LU/Jacobi
753:   test:
754:     suffix: quad_q2q1_full
755:     requires: !single
756:     args: -run_type full -simplex 0 -refinement_limit 0.00625 -bc_type dirichlet -interpolate 1 -vel_petscspace_order 2 -pres_petscspace_order 1 -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_velocity_ksp_type gmres -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_pc_type jacobi -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view -show_solution 0
757:   test:
758:     suffix: quad_q2p1_full
759:     requires: !single
760:     args: -run_type full -simplex 0 -refinement_limit 0.00625 -bc_type dirichlet -interpolate 1 -vel_petscspace_order 2 -pres_petscspace_order 1 -pres_petscspace_poly_tensor 0 -pres_petscdualspace_lagrange_continuity 0 -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_velocity_ksp_type gmres -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_pc_type jacobi -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view -show_solution 0
761:   # Stokes preconditioners 30-36
762:   #   Jacobi
763:   test:
764:     suffix: 30
765:     requires: triangle
766:     filter:  sed -e "s/total number of linear solver iterations=756/total number of linear solver iterations=757/g"
767:     args: -run_type full -refinement_limit 0.00625 -bc_type dirichlet -interpolate 1 -vel_petscspace_order 2 -pres_petscspace_order 1 -ksp_gmres_restart 100 -pc_type jacobi -ksp_rtol 1.0e-9 -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view -show_solution 0
768:   #  Block diagonal \begin{pmatrix} A & 0 \\ 0 & I \end{pmatrix}
769:   test:
770:     suffix: 31
771:     requires: triangle
772:     args: -run_type full -refinement_limit 0.00625 -bc_type dirichlet -interpolate 1 -vel_petscspace_order 2 -pres_petscspace_order 1 -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-4 -pc_type fieldsplit -pc_fieldsplit_type additive -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_pc_type jacobi -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view -show_solution 0
773:   #  Block triangular \begin{pmatrix} A & B \\ 0 & I \end{pmatrix}
774:   test:
775:     suffix: 32
776:     requires: triangle
777:     args: -run_type full -refinement_limit 0.00625 -bc_type dirichlet -interpolate 1 -vel_petscspace_order 2 -pres_petscspace_order 1 -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -pc_type fieldsplit -pc_fieldsplit_type multiplicative -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_pc_type jacobi -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view -show_solution 0
778:   #  Diagonal Schur complement \begin{pmatrix} A & 0 \\ 0 & S \end{pmatrix}
779:   test:
780:     suffix: 33
781:     requires: triangle
782:     args: -run_type full -refinement_limit 0.00625 -bc_type dirichlet -interpolate 1 -vel_petscspace_order 2 -pres_petscspace_order 1 -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type diag -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_velocity_ksp_type gmres -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_pc_type jacobi -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view -show_solution 0
783:   #  Upper triangular Schur complement \begin{pmatrix} A & B \\ 0 & S \end{pmatrix}
784:   test:
785:     suffix: 34
786:     requires: triangle
787:     args: -run_type full -refinement_limit 0.00625 -bc_type dirichlet -interpolate 1 -vel_petscspace_order 2 -pres_petscspace_order 1 -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type upper -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_velocity_ksp_type gmres -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_pc_type jacobi -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view -show_solution 0
788:   #  Lower triangular Schur complement \begin{pmatrix} A & B \\ 0 & S \end{pmatrix}
789:   test:
790:     suffix: 35
791:     requires: triangle
792:     args: -run_type full -refinement_limit 0.00625 -bc_type dirichlet -interpolate 1 -vel_petscspace_order 2 -pres_petscspace_order 1 -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type lower -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_velocity_ksp_type gmres -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_pc_type jacobi -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view -show_solution 0
793:   #  Full Schur complement \begin{pmatrix} I & 0 \\ B^T A^{-1} & I \end{pmatrix} \begin{pmatrix} A & 0 \\ 0 & S \end{pmatrix} \begin{pmatrix} I & A^{-1} B \\ 0 & I \end{pmatrix}
794:   test:
795:     suffix: 36
796:     requires: triangle
797:     args: -run_type full -refinement_limit 0.00625 -bc_type dirichlet -interpolate 1 -vel_petscspace_order 2 -pres_petscspace_order 1 -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_velocity_ksp_type gmres -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_pc_type jacobi -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view -show_solution 0
798:   #  SIMPLE \begin{pmatrix} I & 0 \\ B^T A^{-1} & I \end{pmatrix} \begin{pmatrix} A & 0 \\ 0 & B^T diag(A)^{-1} B \end{pmatrix} \begin{pmatrix} I & diag(A)^{-1} B \\ 0 & I \end{pmatrix}
799:   test:
800:     suffix: pc_simple
801:     requires: triangle
802:     args: -run_type full -refinement_limit 0.00625 -bc_type dirichlet -interpolate 1 -vel_petscspace_order 2 -pres_petscspace_order 1 -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_velocity_ksp_type gmres -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_pc_type jacobi -fieldsplit_pressure_inner_ksp_type preonly -fieldsplit_pressure_inner_pc_type jacobi -fieldsplit_pressure_upper_ksp_type preonly -fieldsplit_pressure_upper_pc_type jacobi -snes_error_if_not_converged -ksp_error_if_not_converged -snes_view -show_solution 0
803:   #  SIMPLEC \begin{pmatrix} I & 0 \\ B^T A^{-1} & I \end{pmatrix} \begin{pmatrix} A & 0 \\ 0 & B^T rowsum(A)^{-1} B \end{pmatrix} \begin{pmatrix} I & rowsum(A)^{-1} B \\ 0 & I \end{pmatrix}
804:   test:
805:     suffix: pc_simplec
806:     requires: triangle
807:     args: -run_type full -refinement_limit 0.00625 -bc_type dirichlet -interpolate 1 -vel_petscspace_order 2 -pres_petscspace_order 1 -ksp_type fgmres -ksp_max_it 5 -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_velocity_ksp_type gmres -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_pc_type jacobi -fieldsplit_pressure_inner_ksp_type preonly -fieldsplit_pressure_inner_pc_type jacobi -fieldsplit_pressure_inner_pc_jacobi_type rowsum -fieldsplit_pressure_upper_ksp_type preonly -fieldsplit_pressure_upper_pc_type jacobi -fieldsplit_pressure_upper_pc_jacobi_type rowsum -snes_converged_reason -ksp_converged_reason -snes_view -show_solution 0
808:   # FETI-DP solvers
809:   test:
810:     suffix: fetidp_2d_tri
811:     requires: triangle suitesparse
812:     nsize: 5
813:     args: -run_type full -dm_refine 2 -bc_type dirichlet -interpolate 1 -vel_petscspace_order 2 -pres_petscspace_order 1 -snes_view -snes_error_if_not_converged -show_solution 0 -dm_mat_type is -ksp_type fetidp -ksp_rtol 1.0e-8 -ksp_fetidp_saddlepoint -fetidp_ksp_type cg -fetidp_fieldsplit_p_ksp_max_it 1 -fetidp_fieldsplit_p_ksp_type richardson -fetidp_fieldsplit_p_ksp_richardson_scale 200 -fetidp_fieldsplit_p_pc_type none -ksp_fetidp_saddlepoint_flip 1 -fetidp_bddc_pc_bddc_vertex_size 2 -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_solver_package umfpack -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_solver_package umfpack -petscpartitioner_type simple -fetidp_fieldsplit_lag_ksp_type preonly
814:   test:
815:     suffix: fetidp_3d_tet
816:     requires: ctetgen mumps suitesparse
817:     nsize: 5
818:     args: -run_type full -dm_refine 2 -bc_type dirichlet -interpolate 1 -vel_petscspace_order 2 -pres_petscspace_order 1 -snes_view -snes_error_if_not_converged -show_solution 0 -dm_mat_type is -ksp_type fetidp -ksp_rtol 1.0e-8 -ksp_fetidp_saddlepoint -fetidp_ksp_type cg -fetidp_fieldsplit_p_ksp_max_it 1 -fetidp_fieldsplit_p_ksp_type richardson -fetidp_fieldsplit_p_ksp_richardson_scale 1000 -fetidp_fieldsplit_p_pc_type none -ksp_fetidp_saddlepoint_flip 1 -fetidp_bddc_pc_bddc_vertex_size 3 -fetidp_bddc_pc_bddc_use_deluxe_scaling -fetidp_bddc_pc_bddc_benign_trick -fetidp_bddc_pc_bddc_deluxe_singlemat -dim 3 -fetidp_pc_discrete_harmonic -fetidp_harmonic_pc_factor_mat_solver_package cholmod -fetidp_harmonic_pc_type cholesky -fetidp_bddelta_pc_factor_mat_solver_package umfpack -fetidp_fieldsplit_lag_ksp_type preonly
819:     filter: sed -s "s/linear solver iterations=10[0-9]/linear solver iterations=100/g"
820:   test:
821:     suffix: fetidp_2d_quad
822:     requires: suitesparse
823:     filter: grep -v "CG or CGNE: variant" | sed -e "s/BDDC: Graph max count: 9223372036854775807/BDDC: Graph max count: 2147483647/g"
824:     nsize: 5
825:     args: -run_type full -dm_refine 2 -bc_type dirichlet -interpolate 1 -vel_petscspace_order 2 -pres_petscspace_order 1 -snes_view -snes_error_if_not_converged -show_solution 0 -dm_mat_type is -ksp_type fetidp -ksp_rtol 1.0e-8 -ksp_fetidp_saddlepoint -fetidp_ksp_type cg -fetidp_fieldsplit_p_ksp_max_it 1 -fetidp_fieldsplit_p_ksp_type richardson -fetidp_fieldsplit_p_ksp_richardson_scale 200 -fetidp_fieldsplit_p_pc_type none -ksp_fetidp_saddlepoint_flip 1 -fetidp_bddc_pc_bddc_vertex_size 2 -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_solver_package umfpack -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_solver_package umfpack -simplex 0 -petscpartitioner_type simple -fetidp_fieldsplit_lag_ksp_type preonly
826:   test:
827:     suffix: fetidp_3d_hex
828:     requires: suitesparse
829:     nsize: 5
830:     args: -run_type full -dm_refine 1 -bc_type dirichlet -interpolate 1 -vel_petscspace_order 2 -pres_petscspace_order 1 -snes_view -snes_error_if_not_converged -show_solution 0 -dm_mat_type is -ksp_type fetidp -ksp_rtol 1.0e-8 -ksp_fetidp_saddlepoint -fetidp_ksp_type cg -fetidp_fieldsplit_p_ksp_max_it 1 -fetidp_fieldsplit_p_ksp_type richardson -fetidp_fieldsplit_p_ksp_richardson_scale 2000 -fetidp_fieldsplit_p_pc_type none -ksp_fetidp_saddlepoint_flip 1 -fetidp_bddc_pc_bddc_vertex_size 3 -dim 3 -simplex 0 -fetidp_pc_discrete_harmonic -fetidp_harmonic_pc_factor_mat_solver_package cholmod -fetidp_harmonic_pc_type cholesky -petscpartitioner_type simple -fetidp_fieldsplit_lag_ksp_type preonly -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_solver_package umfpack -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_solver_package umfpack

832: TEST*/