Actual source code: ex62.c

petsc-master 2014-12-17
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:   /* Problem definition */
 68:   BCType        bcType;
 69:   void       (**exactFuncs)(const PetscReal x[], PetscScalar *u, void *ctx);
 70: } AppCtx;

 72: void zero_scalar(const PetscReal coords[], PetscScalar *u, void *ctx)
 73: {
 74:   u[0] = 0.0;
 75: }
 76: void zero_vector(const PetscReal coords[], PetscScalar *u, void *ctx)
 77: {
 78:   PetscInt d;
 79:   for (d = 0; d < spatialDim; ++d) u[d] = 0.0;
 80: }

 82: /*
 83:   In 2D we use exact solution:

 85:     u = x^2 + y^2
 86:     v = 2 x^2 - 2xy
 87:     p = x + y - 1
 88:     f_x = f_y = 3

 90:   so that

 92:     -\Delta u + \nabla p + f = <-4, -4> + <1, 1> + <3, 3> = 0
 93:     \nabla \cdot u           = 2x - 2x                    = 0
 94: */
 95: void quadratic_u_2d(const PetscReal x[], PetscScalar *u, void *ctx)
 96: {
 97:   u[0] = x[0]*x[0] + x[1]*x[1];
 98:   u[1] = 2.0*x[0]*x[0] - 2.0*x[0]*x[1];
 99: }

101: void linear_p_2d(const PetscReal x[], PetscScalar *p, void *ctx)
102: {
103:   *p = x[0] + x[1] - 1.0;
104: }
105: void constant_p(const PetscReal x[], PetscScalar *p, void *ctx)
106: {
107:   *p = 1.0;
108: }

110: void f0_u(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar f0[])
111: {
112:   PetscInt c;
113:   for (c = 0; c < spatialDim; ++c) f0[c] = 3.0;
114: }

116: /* 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}
117:    u[Ncomp]          = {p} */
118: void f1_u(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar f1[])
119: {
120:   const PetscInt dim   = spatialDim;
121:   const PetscInt Ncomp = spatialDim;
122:   PetscInt       comp, d;

124:   for (comp = 0; comp < Ncomp; ++comp) {
125:     for (d = 0; d < dim; ++d) {
126:       /* f1[comp*dim+d] = 0.5*(gradU[comp*dim+d] + gradU[d*dim+comp]); */
127:       f1[comp*dim+d] = u_x[comp*dim+d];
128:     }
129:     f1[comp*dim+comp] -= u[Ncomp];
130:   }
131: }

133: /* 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} */
134: void f0_p(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar f0[])
135: {
136:   const PetscInt dim = spatialDim;
137:   PetscInt       d;

139:   f0[0] = 0.0;
140:   for (d = 0; d < dim; ++d) f0[0] += u_x[d*dim+d];
141: }

143: void f1_p(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar f1[])
144: {
145:   PetscInt d;
146:   for (d = 0; d < spatialDim; ++d) f1[d] = 0.0;
147: }

149: /* < q, \nabla\cdot u >
150:    NcompI = 1, NcompJ = dim */
151: void g1_pu(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar g1[])
152: {
153:   const PetscInt dim = spatialDim;
154:   PetscInt       d;

156:   for (d = 0; d < dim; ++d) g1[d*dim+d] = 1.0; /* \frac{\partial\phi^{u_d}}{\partial x_d} */
157: }

159: /* -< \nabla\cdot v, p >
160:     NcompI = dim, NcompJ = 1 */
161: void g2_up(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar g2[])
162: {
163:   const PetscInt dim = spatialDim;
164:   PetscInt       d;

166:   for (d = 0; d < dim; ++d) g2[d*dim+d] = -1.0; /* \frac{\partial\psi^{u_d}}{\partial x_d} */
167: }

169: /* < \nabla v, \nabla u + {\nabla u}^T >
170:    This just gives \nabla u, give the perdiagonal for the transpose */
171: void g3_uu(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar g3[])
172: {
173:   const PetscInt dim   = spatialDim;
174:   const PetscInt Ncomp = spatialDim;
175:   PetscInt       compI, d;

177:   for (compI = 0; compI < Ncomp; ++compI) {
178:     for (d = 0; d < dim; ++d) {
179:       g3[((compI*Ncomp+compI)*dim+d)*dim+d] = 1.0;
180:     }
181:   }
182: }

184: /*
185:   In 3D we use exact solution:

187:     u = x^2 + y^2
188:     v = y^2 + z^2
189:     w = x^2 + y^2 - 2(x+y)z
190:     p = x + y + z - 3/2
191:     f_x = f_y = f_z = 3

193:   so that

195:     -\Delta u + \nabla p + f = <-4, -4, -4> + <1, 1, 1> + <3, 3, 3> = 0
196:     \nabla \cdot u           = 2x + 2y - 2(x + y)                   = 0
197: */
198: void quadratic_u_3d(const PetscReal x[], PetscScalar *u, void *ctx)
199: {
200:   u[0] = x[0]*x[0] + x[1]*x[1];
201:   u[1] = x[1]*x[1] + x[2]*x[2];
202:   u[2] = x[0]*x[0] + x[1]*x[1] - 2.0*(x[0] + x[1])*x[2];
203: }

205: void linear_p_3d(const PetscReal x[], PetscScalar *p, void *ctx)
206: {
207:   *p = x[0] + x[1] + x[2] - 1.5;
208: }

212: PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
213: {
214:   const char    *bcTypes[2]  = {"neumann", "dirichlet"};
215:   const char    *runTypes[2] = {"full", "test"};
216:   PetscInt       bc, run;

220:   options->debug           = 0;
221:   options->runType         = RUN_FULL;
222:   options->dim             = 2;
223:   options->interpolate     = PETSC_FALSE;
224:   options->simplex         = PETSC_TRUE;
225:   options->refinementLimit = 0.0;
226:   options->bcType          = DIRICHLET;
227:   options->showInitial     = PETSC_FALSE;
228:   options->showSolution    = PETSC_TRUE;
229:   options->showError       = PETSC_FALSE;

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

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

238:   PetscOptionsInt("-dim", "The topological mesh dimension", "ex62.c", options->dim, &options->dim, NULL);
239:   spatialDim = options->dim;
240:   PetscOptionsBool("-interpolate", "Generate intermediate mesh elements", "ex62.c", options->interpolate, &options->interpolate, NULL);
241:   PetscOptionsBool("-simplex", "Use simplices or tensor product cells", "ex62.c", options->simplex, &options->simplex, NULL);
242:   PetscOptionsReal("-refinement_limit", "The largest allowable cell volume", "ex62.c", options->refinementLimit, &options->refinementLimit, NULL);
243:   bc   = options->bcType;
244:   PetscOptionsEList("-bc_type","Type of boundary condition","ex62.c",bcTypes,2,bcTypes[options->bcType],&bc,NULL);

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

248:   PetscOptionsBool("-show_initial", "Output the initial guess for verification", "ex62.c", options->showInitial, &options->showInitial, NULL);
249:   PetscOptionsBool("-show_solution", "Output the solution for verification", "ex62.c", options->showSolution, &options->showSolution, NULL);
250:   PetscOptionsBool("-show_error", "Output the error for verification", "ex62.c", options->showError, &options->showError, NULL);
251:   PetscOptionsEnd();

253:   PetscLogEventRegister("CreateMesh", DM_CLASSID, &options->createMeshEvent);
254:   return(0);
255: }

259: PetscErrorCode DMVecViewLocal(DM dm, Vec v, PetscViewer viewer)
260: {
261:   Vec            lv;
262:   PetscInt       p;
263:   PetscMPIInt    rank, numProcs;

267:   MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);
268:   MPI_Comm_size(PetscObjectComm((PetscObject)dm), &numProcs);
269:   DMGetLocalVector(dm, &lv);
270:   DMGlobalToLocalBegin(dm, v, INSERT_VALUES, lv);
271:   DMGlobalToLocalEnd(dm, v, INSERT_VALUES, lv);
272:   PetscPrintf(PETSC_COMM_WORLD, "Local function\n");
273:   for (p = 0; p < numProcs; ++p) {
274:     if (p == rank) {VecView(lv, PETSC_VIEWER_STDOUT_SELF);}
275:     PetscBarrier((PetscObject) dm);
276:   }
277:   DMRestoreLocalVector(dm, &lv);
278:   return(0);
279: }

283: PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
284: {
285:   PetscInt       dim             = user->dim;
286:   PetscBool      interpolate     = user->interpolate;
287:   PetscReal      refinementLimit = user->refinementLimit;
288:   const PetscInt cells[3]        = {3, 3, 3};

292:   PetscLogEventBegin(user->createMeshEvent,0,0,0,0);
293:   if (user->simplex) {DMPlexCreateBoxMesh(comm, dim, interpolate, dm);}
294:   else               {DMPlexCreateHexBoxMesh(comm, dim, cells, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, dm);}
295:   {
296:     DM refinedMesh     = NULL;
297:     DM distributedMesh = NULL;

299:     /* Refine mesh using a volume constraint */
300:     DMPlexSetRefinementLimit(*dm, refinementLimit);
301:     if (user->simplex) {DMRefine(*dm, comm, &refinedMesh);}
302:     if (refinedMesh) {
303:       DMDestroy(dm);
304:       *dm  = refinedMesh;
305:     }
306:     /* Distribute mesh over processes */
307:     DMPlexDistribute(*dm, 0, NULL, &distributedMesh);
308:     if (distributedMesh) {
309:       DMDestroy(dm);
310:       *dm  = distributedMesh;
311:     }
312:   }
313:   DMSetFromOptions(*dm);
314:   DMViewFromOptions(*dm, NULL, "-dm_view");
315:   PetscLogEventEnd(user->createMeshEvent,0,0,0,0);
316:   return(0);
317: }

321: PetscErrorCode SetupProblem(DM dm, AppCtx *user)
322: {
323:   PetscDS        prob;

327:   DMGetDS(dm, &prob);
328:   PetscDSSetResidual(prob, 0, f0_u, f1_u);
329:   PetscDSSetResidual(prob, 1, f0_p, f1_p);
330:   PetscDSSetJacobian(prob, 0, 0, NULL, NULL,  NULL,  g3_uu);
331:   PetscDSSetJacobian(prob, 0, 1, NULL, NULL,  g2_up, NULL);
332:   PetscDSSetJacobian(prob, 1, 0, NULL, g1_pu, NULL,  NULL);
333:   switch (user->dim) {
334:   case 2:
335:     user->exactFuncs[0] = quadratic_u_2d;
336:     user->exactFuncs[1] = linear_p_2d;
337:     break;
338:   case 3:
339:     user->exactFuncs[0] = quadratic_u_3d;
340:     user->exactFuncs[1] = linear_p_3d;
341:     break;
342:   default:
343:     SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %d", user->dim);
344:   }
345:   return(0);
346: }

350: PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
351: {
352:   DM              cdm   = dm;
353:   const PetscInt  dim   = user->dim;
354:   const PetscInt  id    = 1;
355:   PetscFE         fe[2];
356:   PetscQuadrature q;
357:   PetscDS         prob;
358:   PetscInt        order;
359:   PetscErrorCode  ierr;

362:   /* Create finite element */
363:   PetscFECreateDefault(dm, dim, dim, user->simplex, "vel_", -1, &fe[0]);
364:   PetscObjectSetName((PetscObject) fe[0], "velocity");
365:   PetscFEGetQuadrature(fe[0], &q);
366:   PetscQuadratureGetOrder(q, &order);
367:   PetscFECreateDefault(dm, dim, 1, user->simplex, "pres_", order, &fe[1]);
368:   PetscObjectSetName((PetscObject) fe[1], "pressure");
369:   /* Set discretization and boundary conditions for each mesh */
370:   while (cdm) {
371:     DMGetDS(cdm, &prob);
372:     PetscDSSetDiscretization(prob, 0, (PetscObject) fe[0]);
373:     PetscDSSetDiscretization(prob, 1, (PetscObject) fe[1]);
374:     SetupProblem(cdm, user);
375:     DMPlexAddBoundary(cdm, user->bcType == DIRICHLET, "wall", user->bcType == NEUMANN ? "boundary" : "marker", 0, user->exactFuncs[0], 1, &id, user);
376:     DMPlexGetCoarseDM(cdm, &cdm);
377:   }
378:   PetscFEDestroy(&fe[0]);
379:   PetscFEDestroy(&fe[1]);
380:   return(0);
381: }

385: PetscErrorCode CreatePressureNullSpace(DM dm, AppCtx *user, Vec *v, MatNullSpace *nullSpace)
386: {
387:   Vec            vec;
388:   void         (*funcs[2])(const PetscReal x[], PetscScalar *u, void* ctx) = {zero_vector, constant_p};

392:   DMGetGlobalVector(dm, &vec);
393:   DMPlexProjectFunction(dm, funcs, NULL, INSERT_ALL_VALUES, vec);
394:   VecNormalize(vec, NULL);
395:   if (user->debug) {
396:     PetscPrintf(PetscObjectComm((PetscObject)dm), "Pressure Null Space\n");
397:     VecView(vec, PETSC_VIEWER_STDOUT_WORLD);
398:   }
399:   MatNullSpaceCreate(PetscObjectComm((PetscObject)dm), PETSC_FALSE, 1, &vec, nullSpace);
400:   if (v) {
401:     DMCreateGlobalVector(dm, v);
402:     VecCopy(vec, *v);
403:   }
404:   DMRestoreGlobalVector(dm, &vec);
405:   /* New style for field null spaces */
406:   {
407:     PetscObject  pressure;
408:     MatNullSpace nullSpacePres;

410:     DMGetField(dm, 1, &pressure);
411:     MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nullSpacePres);
412:     PetscObjectCompose(pressure, "nullspace", (PetscObject) nullSpacePres);
413:     MatNullSpaceDestroy(&nullSpacePres);
414:   }
415:   return(0);
416: }

420: int main(int argc, char **argv)
421: {
422:   SNES           snes;                 /* nonlinear solver */
423:   DM             dm;                   /* problem definition */
424:   Vec            u,r;                  /* solution, residual vectors */
425:   Mat            A,J;                  /* Jacobian matrix */
426:   MatNullSpace   nullSpace;            /* May be necessary for pressure */
427:   AppCtx         user;                 /* user-defined work context */
428:   PetscInt       its;                  /* iterations for convergence */
429:   PetscReal      error         = 0.0;  /* L_2 error in the solution */
430:   PetscReal      ferrors[2];

433:   PetscInitialize(&argc, &argv, NULL, help);
434:   ProcessOptions(PETSC_COMM_WORLD, &user);
435:   SNESCreate(PETSC_COMM_WORLD, &snes);
436:   CreateMesh(PETSC_COMM_WORLD, &user, &dm);
437:   SNESSetDM(snes, dm);
438:   DMSetApplicationContext(dm, &user);

440:   PetscMalloc(2 * sizeof(void (*)(const PetscReal[], PetscScalar *, void *)), &user.exactFuncs);
441:   SetupDiscretization(dm, &user);
442:   DMPlexCreateClosureIndex(dm, NULL);

444:   DMCreateGlobalVector(dm, &u);
445:   VecDuplicate(u, &r);

447:   DMSetMatType(dm,MATAIJ);
448:   DMCreateMatrix(dm, &J);
449:   A = J;
450:   CreatePressureNullSpace(dm, &user, NULL, &nullSpace);
451:   MatSetNullSpace(J, nullSpace);
452:   if (A != J) {
453:     MatSetNullSpace(A, nullSpace);
454:   }

456:   DMSNESSetFunctionLocal(dm,  (PetscErrorCode (*)(DM,Vec,Vec,void*))DMPlexSNESComputeResidualFEM,&user);
457:   DMSNESSetJacobianLocal(dm,  (PetscErrorCode (*)(DM,Vec,Mat,Mat,void*))DMPlexSNESComputeJacobianFEM,&user);
458:   SNESSetJacobian(snes, A, J, NULL, NULL);

460:   SNESSetFromOptions(snes);

462:   DMPlexProjectFunction(dm, user.exactFuncs, NULL, INSERT_ALL_VALUES, u);
463:   if (user.showInitial) {DMVecViewLocal(dm, u, PETSC_VIEWER_STDOUT_SELF);}
464:   if (user.runType == RUN_FULL) {
465:     void (*initialGuess[2])(const PetscReal x[], PetscScalar *u, void* ctx) = {zero_vector, zero_scalar};

467:     DMPlexProjectFunction(dm, initialGuess, NULL, INSERT_VALUES, u);
468:     if (user.debug) {
469:       PetscPrintf(PETSC_COMM_WORLD, "Initial guess\n");
470:       VecView(u, PETSC_VIEWER_STDOUT_WORLD);
471:     }
472:     SNESSolve(snes, NULL, u);
473:     SNESGetIterationNumber(snes, &its);
474:     PetscPrintf(PETSC_COMM_WORLD, "Number of SNES iterations = %D\n", its);
475:     DMPlexComputeL2Diff(dm, user.exactFuncs, NULL, u, &error);
476:     DMPlexComputeL2FieldDiff(dm, user.exactFuncs, NULL, u, ferrors);
477:     PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: %.3g [%.3g, %.3g]\n", error, ferrors[0], ferrors[1]);
478:     if (user.showError) {
479:       Vec r;
480:       DMGetGlobalVector(dm, &r);
481:       DMPlexProjectFunction(dm, user.exactFuncs, NULL, INSERT_ALL_VALUES, r);
482:       VecAXPY(r, -1.0, u);
483:       PetscPrintf(PETSC_COMM_WORLD, "Solution Error\n");
484:       VecView(r, PETSC_VIEWER_STDOUT_WORLD);
485:       DMRestoreGlobalVector(dm, &r);
486:     }
487:     if (user.showSolution) {
488:       PetscPrintf(PETSC_COMM_WORLD, "Solution\n");
489:       VecChop(u, 3.0e-9);
490:       VecView(u, PETSC_VIEWER_STDOUT_WORLD);
491:     }
492:   } else {
493:     PetscReal res = 0.0;

495:     /* Check discretization error */
496:     PetscPrintf(PETSC_COMM_WORLD, "Initial guess\n");
497:     VecView(u, PETSC_VIEWER_STDOUT_WORLD);
498:     DMPlexComputeL2Diff(dm, user.exactFuncs, NULL, u, &error);
499:     if (error >= 1.0e-11) {PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: %g\n", error);}
500:     else                  {PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: < 1.0e-11\n");}
501:     /* Check residual */
502:     SNESComputeFunction(snes, u, r);
503:     PetscPrintf(PETSC_COMM_WORLD, "Initial Residual\n");
504:     VecChop(r, 1.0e-10);
505:     VecView(r, PETSC_VIEWER_STDOUT_WORLD);
506:     VecNorm(r, NORM_2, &res);
507:     PetscPrintf(PETSC_COMM_WORLD, "L_2 Residual: %g\n", res);
508:     /* Check Jacobian */
509:     {
510:       Vec          b;
511:       PetscBool    isNull;

513:       SNESComputeJacobian(snes, u, A, A);
514:       MatNullSpaceTest(nullSpace, J, &isNull);
515:       //if (!isNull) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "The null space calculated for the system operator is invalid.");
516:       VecDuplicate(u, &b);
517:       VecSet(r, 0.0);
518:       SNESComputeFunction(snes, r, b);
519:       MatMult(A, u, r);
520:       VecAXPY(r, 1.0, b);
521:       VecDestroy(&b);
522:       PetscPrintf(PETSC_COMM_WORLD, "Au - b = Au + F(0)\n");
523:       VecChop(r, 1.0e-10);
524:       VecView(r, PETSC_VIEWER_STDOUT_WORLD);
525:       VecNorm(r, NORM_2, &res);
526:       PetscPrintf(PETSC_COMM_WORLD, "Linear L_2 Residual: %g\n", res);
527:     }
528:   }
529:   VecViewFromOptions(u, NULL, "-sol_vec_view");

531:   MatNullSpaceDestroy(&nullSpace);
532:   if (A != J) {MatDestroy(&A);}
533:   MatDestroy(&J);
534:   VecDestroy(&u);
535:   VecDestroy(&r);
536:   SNESDestroy(&snes);
537:   DMDestroy(&dm);
538:   PetscFree(user.exactFuncs);
539:   PetscFinalize();
540:   return 0;
541: }