Actual source code: ex62.c

petsc-3.7.3 2016-08-01
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[], 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[], 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[], 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[], 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[], 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[], 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[], 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: }

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

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

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

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

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

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

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

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

281: PetscErrorCode DMVecViewLocal(DM dm, Vec v, PetscViewer viewer)
282: {
283:   Vec            lv;
284:   PetscInt       p;
285:   PetscMPIInt    rank, numProcs;

289:   MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);
290:   MPI_Comm_size(PetscObjectComm((PetscObject)dm), &numProcs);
291:   DMGetLocalVector(dm, &lv);
292:   DMGlobalToLocalBegin(dm, v, INSERT_VALUES, lv);
293:   DMGlobalToLocalEnd(dm, v, INSERT_VALUES, lv);
294:   PetscPrintf(PETSC_COMM_WORLD, "Local function\n");
295:   for (p = 0; p < numProcs; ++p) {
296:     if (p == rank) {VecView(lv, PETSC_VIEWER_STDOUT_SELF);}
297:     PetscBarrier((PetscObject) dm);
298:   }
299:   DMRestoreLocalVector(dm, &lv);
300:   return(0);
301: }

305: PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
306: {
307:   PetscInt       dim             = user->dim;
308:   PetscBool      interpolate     = user->interpolate;
309:   PetscReal      refinementLimit = user->refinementLimit;
310:   const PetscInt cells[3]        = {3, 3, 3};

314:   PetscLogEventBegin(user->createMeshEvent,0,0,0,0);
315:   if (user->simplex) {DMPlexCreateBoxMesh(comm, dim, interpolate, dm);}
316:   else               {DMPlexCreateHexBoxMesh(comm, dim, cells, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, dm);}
317:   {
318:     DM refinedMesh     = NULL;
319:     DM distributedMesh = NULL;

321:     /* Refine mesh using a volume constraint */
322:     DMPlexSetRefinementLimit(*dm, refinementLimit);
323:     if (user->simplex) {DMRefine(*dm, comm, &refinedMesh);}
324:     if (refinedMesh) {
325:       DMDestroy(dm);
326:       *dm  = refinedMesh;
327:     }
328:     /* Setup test partitioning */
329:     if (user->testPartition) {
330:       PetscInt         triSizes_n2[2]       = {4, 4};
331:       PetscInt         triPoints_n2[8]      = {3, 5, 6, 7, 0, 1, 2, 4};
332:       PetscInt         triSizes_n3[3]       = {2, 3, 3};
333:       PetscInt         triPoints_n3[8]      = {3, 5, 1, 6, 7, 0, 2, 4};
334:       PetscInt         triSizes_n5[5]       = {1, 2, 2, 1, 2};
335:       PetscInt         triPoints_n5[8]      = {3, 5, 6, 4, 7, 0, 1, 2};
336:       PetscInt         triSizes_ref_n2[2]   = {8, 8};
337:       PetscInt         triPoints_ref_n2[16] = {1, 5, 6, 7, 10, 11, 14, 15, 0, 2, 3, 4, 8, 9, 12, 13};
338:       PetscInt         triSizes_ref_n3[3]   = {5, 6, 5};
339:       PetscInt         triPoints_ref_n3[16] = {1, 7, 10, 14, 15, 2, 6, 8, 11, 12, 13, 0, 3, 4, 5, 9};
340:       PetscInt         triSizes_ref_n5[5]   = {3, 4, 3, 3, 3};
341:       PetscInt         triPoints_ref_n5[16] = {1, 7, 10, 2, 11, 13, 14, 5, 6, 15, 0, 8, 9, 3, 4, 12};
342:       const PetscInt  *sizes = NULL;
343:       const PetscInt  *points = NULL;
344:       PetscPartitioner part;
345:       PetscInt         cEnd;
346:       PetscMPIInt      rank, numProcs;

348:       MPI_Comm_rank(comm, &rank);
349:       MPI_Comm_size(comm, &numProcs);
350:       DMPlexGetHeightStratum(*dm, 0, NULL, &cEnd);
351:       if (!rank) {
352:         if (dim == 2 && user->simplex && numProcs == 2 && cEnd == 8) {
353:            sizes = triSizes_n2; points = triPoints_n2;
354:         } else if (dim == 2 && user->simplex && numProcs == 3 && cEnd == 8) {
355:           sizes = triSizes_n3; points = triPoints_n3;
356:         } else if (dim == 2 && user->simplex && numProcs == 5 && cEnd == 8) {
357:           sizes = triSizes_n5; points = triPoints_n5;
358:         } else if (dim == 2 && user->simplex && numProcs == 2 && cEnd == 16) {
359:            sizes = triSizes_ref_n2; points = triPoints_ref_n2;
360:         } else if (dim == 2 && user->simplex && numProcs == 3 && cEnd == 16) {
361:           sizes = triSizes_ref_n3; points = triPoints_ref_n3;
362:         } else if (dim == 2 && user->simplex && numProcs == 5 && cEnd == 16) {
363:           sizes = triSizes_ref_n5; points = triPoints_ref_n5;
364:         } else SETERRQ(comm, PETSC_ERR_ARG_WRONG, "No stored partition matching run parameters");
365:       }
366:       DMPlexGetPartitioner(*dm, &part);
367:       PetscPartitionerSetType(part, PETSCPARTITIONERSHELL);
368:       PetscPartitionerShellSetPartition(part, numProcs, sizes, points);
369:     }
370:     /* Distribute mesh over processes */
371:     DMPlexDistribute(*dm, 0, NULL, &distributedMesh);
372:     if (distributedMesh) {
373:       DMDestroy(dm);
374:       *dm  = distributedMesh;
375:     }
376:   }
377:   DMSetFromOptions(*dm);
378:   DMViewFromOptions(*dm, NULL, "-dm_view");
379:   PetscLogEventEnd(user->createMeshEvent,0,0,0,0);
380:   return(0);
381: }

385: PetscErrorCode SetupProblem(DM dm, AppCtx *user)
386: {
387:   PetscDS        prob;

391:   DMGetDS(dm, &prob);
392:   PetscDSSetResidual(prob, 0, f0_u, f1_u);
393:   PetscDSSetResidual(prob, 1, f0_p, f1_p);
394:   PetscDSSetJacobian(prob, 0, 0, NULL, NULL,  NULL,  g3_uu);
395:   PetscDSSetJacobian(prob, 0, 1, NULL, NULL,  g2_up, NULL);
396:   PetscDSSetJacobian(prob, 1, 0, NULL, g1_pu, NULL,  NULL);
397:   switch (user->dim) {
398:   case 2:
399:     user->exactFuncs[0] = quadratic_u_2d;
400:     user->exactFuncs[1] = linear_p_2d;
401:     break;
402:   case 3:
403:     user->exactFuncs[0] = quadratic_u_3d;
404:     user->exactFuncs[1] = linear_p_3d;
405:     break;
406:   default:
407:     SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %d", user->dim);
408:   }
409:   return(0);
410: }

414: PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
415: {
416:   DM              cdm   = dm;
417:   const PetscInt  dim   = user->dim;
418:   const PetscInt  id    = 1;
419:   PetscFE         fe[2];
420:   PetscQuadrature q;
421:   PetscDS         prob;
422:   PetscInt        order;
423:   PetscErrorCode  ierr;

426:   /* Create finite element */
427:   PetscFECreateDefault(dm, dim, dim, user->simplex, "vel_", -1, &fe[0]);
428:   PetscObjectSetName((PetscObject) fe[0], "velocity");
429:   PetscFEGetQuadrature(fe[0], &q);
430:   PetscQuadratureGetOrder(q, &order);
431:   PetscFECreateDefault(dm, dim, 1, user->simplex, "pres_", order, &fe[1]);
432:   PetscObjectSetName((PetscObject) fe[1], "pressure");
433:   /* Set discretization and boundary conditions for each mesh */
434:   while (cdm) {
435:     DMGetDS(cdm, &prob);
436:     PetscDSSetDiscretization(prob, 0, (PetscObject) fe[0]);
437:     PetscDSSetDiscretization(prob, 1, (PetscObject) fe[1]);
438:     SetupProblem(cdm, user);
439:     DMAddBoundary(cdm, user->bcType == DIRICHLET ? PETSC_TRUE : PETSC_FALSE, "wall", user->bcType == NEUMANN ? "boundary" : "marker", 0, 0, NULL, (void (*)()) user->exactFuncs[0], 1, &id, user);
440:     DMGetCoarseDM(cdm, &cdm);
441:   }
442:   PetscFEDestroy(&fe[0]);
443:   PetscFEDestroy(&fe[1]);
444:   return(0);
445: }

449: PetscErrorCode CreatePressureNullSpace(DM dm, AppCtx *user, Vec *v, MatNullSpace *nullSpace)
450: {
451:   Vec              vec;
452:   PetscErrorCode (*funcs[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void* ctx) = {zero_vector, constant_p};
453:   PetscErrorCode   ierr;

456:   DMGetGlobalVector(dm, &vec);
457:   DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec);
458:   VecNormalize(vec, NULL);
459:   if (user->debug) {
460:     PetscPrintf(PetscObjectComm((PetscObject)dm), "Pressure Null Space\n");
461:     VecView(vec, PETSC_VIEWER_STDOUT_WORLD);
462:   }
463:   MatNullSpaceCreate(PetscObjectComm((PetscObject)dm), PETSC_FALSE, 1, &vec, nullSpace);
464:   if (v) {
465:     DMCreateGlobalVector(dm, v);
466:     VecCopy(vec, *v);
467:   }
468:   DMRestoreGlobalVector(dm, &vec);
469:   /* New style for field null spaces */
470:   {
471:     PetscObject  pressure;
472:     MatNullSpace nullSpacePres;

474:     DMGetField(dm, 1, &pressure);
475:     MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nullSpacePres);
476:     PetscObjectCompose(pressure, "nullspace", (PetscObject) nullSpacePres);
477:     MatNullSpaceDestroy(&nullSpacePres);
478:   }
479:   return(0);
480: }

484: int main(int argc, char **argv)
485: {
486:   SNES           snes;                 /* nonlinear solver */
487:   DM             dm;                   /* problem definition */
488:   Vec            u,r;                  /* solution, residual vectors */
489:   Mat            A,J;                  /* Jacobian matrix */
490:   MatNullSpace   nullSpace;            /* May be necessary for pressure */
491:   AppCtx         user;                 /* user-defined work context */
492:   PetscInt       its;                  /* iterations for convergence */
493:   PetscReal      error         = 0.0;  /* L_2 error in the solution */
494:   PetscReal      ferrors[2];

497:   PetscInitialize(&argc, &argv, NULL, help);
498:   ProcessOptions(PETSC_COMM_WORLD, &user);
499:   SNESCreate(PETSC_COMM_WORLD, &snes);
500:   CreateMesh(PETSC_COMM_WORLD, &user, &dm);
501:   SNESSetDM(snes, dm);
502:   DMSetApplicationContext(dm, &user);

504:   PetscMalloc(2 * sizeof(void (*)(const PetscReal[], PetscScalar *, void *)), &user.exactFuncs);
505:   SetupDiscretization(dm, &user);
506:   DMPlexCreateClosureIndex(dm, NULL);

508:   DMCreateGlobalVector(dm, &u);
509:   VecDuplicate(u, &r);

511:   DMSetMatType(dm,MATAIJ);
512:   DMCreateMatrix(dm, &J);
513:   A = J;
514:   CreatePressureNullSpace(dm, &user, NULL, &nullSpace);
515:   MatSetNullSpace(A, nullSpace);

517:   DMPlexSetSNESLocalFEM(dm,&user,&user,&user);
518:   SNESSetJacobian(snes, A, J, NULL, NULL);

520:   SNESSetFromOptions(snes);

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

527:     DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u);
528:     if (user.debug) {
529:       PetscPrintf(PETSC_COMM_WORLD, "Initial guess\n");
530:       VecView(u, PETSC_VIEWER_STDOUT_WORLD);
531:     }
532:     SNESSolve(snes, NULL, u);
533:     SNESGetIterationNumber(snes, &its);
534:     PetscPrintf(PETSC_COMM_WORLD, "Number of SNES iterations = %D\n", its);
535:     DMComputeL2Diff(dm, 0.0, user.exactFuncs, NULL, u, &error);
536:     DMComputeL2FieldDiff(dm, 0.0, user.exactFuncs, NULL, u, ferrors);
537:     PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: %.3g [%.3g, %.3g]\n", error, ferrors[0], ferrors[1]);
538:     if (user.showError) {
539:       Vec r;
540:       DMGetGlobalVector(dm, &r);
541:       DMProjectFunction(dm, 0.0, user.exactFuncs, NULL, INSERT_ALL_VALUES, r);
542:       VecAXPY(r, -1.0, u);
543:       PetscPrintf(PETSC_COMM_WORLD, "Solution Error\n");
544:       VecView(r, PETSC_VIEWER_STDOUT_WORLD);
545:       DMRestoreGlobalVector(dm, &r);
546:     }
547:     if (user.showSolution) {
548:       PetscPrintf(PETSC_COMM_WORLD, "Solution\n");
549:       VecChop(u, 3.0e-9);
550:       VecView(u, PETSC_VIEWER_STDOUT_WORLD);
551:     }
552:   } else {
553:     PetscReal res = 0.0;

555:     /* Check discretization error */
556:     PetscPrintf(PETSC_COMM_WORLD, "Initial guess\n");
557:     VecView(u, PETSC_VIEWER_STDOUT_WORLD);
558:     DMComputeL2Diff(dm, 0.0, user.exactFuncs, NULL, u, &error);
559:     if (error >= 1.0e-11) {PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: %g\n", error);}
560:     else                  {PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: < 1.0e-11\n");}
561:     /* Check residual */
562:     SNESComputeFunction(snes, u, r);
563:     PetscPrintf(PETSC_COMM_WORLD, "Initial Residual\n");
564:     VecChop(r, 1.0e-10);
565:     VecView(r, PETSC_VIEWER_STDOUT_WORLD);
566:     VecNorm(r, NORM_2, &res);
567:     PetscPrintf(PETSC_COMM_WORLD, "L_2 Residual: %g\n", res);
568:     /* Check Jacobian */
569:     {
570:       Vec          b;
571:       PetscBool    isNull;

573:       SNESComputeJacobian(snes, u, A, A);
574:       MatNullSpaceTest(nullSpace, J, &isNull);
575:       //if (!isNull) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "The null space calculated for the system operator is invalid.");
576:       VecDuplicate(u, &b);
577:       VecSet(r, 0.0);
578:       SNESComputeFunction(snes, r, b);
579:       MatMult(A, u, r);
580:       VecAXPY(r, 1.0, b);
581:       VecDestroy(&b);
582:       PetscPrintf(PETSC_COMM_WORLD, "Au - b = Au + F(0)\n");
583:       VecChop(r, 1.0e-10);
584:       VecView(r, PETSC_VIEWER_STDOUT_WORLD);
585:       VecNorm(r, NORM_2, &res);
586:       PetscPrintf(PETSC_COMM_WORLD, "Linear L_2 Residual: %g\n", res);
587:     }
588:   }
589:   VecViewFromOptions(u, NULL, "-sol_vec_view");

591:   MatNullSpaceDestroy(&nullSpace);
592:   if (A != J) {MatDestroy(&A);}
593:   MatDestroy(&J);
594:   VecDestroy(&u);
595:   VecDestroy(&r);
596:   SNESDestroy(&snes);
597:   DMDestroy(&dm);
598:   PetscFree(user.exactFuncs);
599:   PetscFinalize();
600:   return 0;
601: }