xref: /petsc/src/snes/utils/dmplexsnes.c (revision 909dfd52d29915732d44f9443658cd10cc8ebe94)
1 #include <petsc/private/dmpleximpl.h>   /*I "petscdmplex.h" I*/
2 #include <petsc/private/snesimpl.h>     /*I "petscsnes.h"   I*/
3 #include <petscds.h>
4 #include <petsc/private/petscimpl.h>
5 #include <petsc/private/petscfeimpl.h>
6 
7 static void pressure_Private(PetscInt dim, PetscInt Nf, PetscInt NfAux,
8                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
9                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
10                              PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar p[])
11 {
12   p[0] = u[uOff[1]];
13 }
14 
15 /*
16   SNESCorrectDiscretePressure_Private - Add a vector in the nullspace to make the continuum integral of the pressure field equal to zero. This is normally used only to evaluate convergence rates for the pressure accurately.
17 
18   Collective on SNES
19 
20   Input Parameters:
21 + snes      - The SNES
22 . pfield    - The field number for pressure
23 . nullspace - The pressure nullspace
24 . u         - The solution vector
25 - ctx       - An optional user context
26 
27   Output Parameter:
28 . u         - The solution with a continuum pressure integral of zero
29 
30   Notes:
31   If int(u) = a and int(n) = b, then int(u - a/b n) = a - a/b b = 0. We assume that the nullspace is a single vector given explicitly.
32 
33   Level: developer
34 
35 .seealso: SNESConvergedCorrectPressure()
36 */
37 static PetscErrorCode SNESCorrectDiscretePressure_Private(SNES snes, PetscInt pfield, MatNullSpace nullspace, Vec u, void *ctx)
38 {
39   DM             dm;
40   PetscDS        ds;
41   const Vec     *nullvecs;
42   PetscScalar    pintd, *intc, *intn;
43   MPI_Comm       comm;
44   PetscInt       Nf, Nv;
45   PetscErrorCode ierr;
46 
47   PetscFunctionBegin;
48   ierr = PetscObjectGetComm((PetscObject) snes, &comm);CHKERRQ(ierr);
49   ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr);
50   if (!dm) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Cannot compute test without a SNES DM");
51   if (!nullspace) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Cannot compute test without a Jacobian nullspace");
52   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
53   ierr = PetscDSSetObjective(ds, pfield, pressure_Private);CHKERRQ(ierr);
54   ierr = MatNullSpaceGetVecs(nullspace, NULL, &Nv, &nullvecs);CHKERRQ(ierr);
55   if (Nv != 1) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Can only handle a single null vector for pressure, not %D", Nv);
56   ierr = VecDot(nullvecs[0], u, &pintd);CHKERRQ(ierr);
57   if (PetscAbsScalar(pintd) > PETSC_SMALL) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Discrete integral of pressure: %g\n", (double) PetscRealPart(pintd));
58   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
59   ierr = PetscMalloc2(Nf, &intc, Nf, &intn);CHKERRQ(ierr);
60   ierr = DMPlexComputeIntegralFEM(dm, nullvecs[0], intn, ctx);CHKERRQ(ierr);
61   ierr = DMPlexComputeIntegralFEM(dm, u, intc, ctx);CHKERRQ(ierr);
62   ierr = VecAXPY(u, -intc[pfield]/intn[pfield], nullvecs[0]);CHKERRQ(ierr);
63 #if defined (PETSC_USE_DEBUG)
64   ierr = DMPlexComputeIntegralFEM(dm, u, intc, ctx);CHKERRQ(ierr);
65   if (PetscAbsScalar(intc[pfield]) > PETSC_SMALL) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Continuum integral of pressure after correction: %g\n", (double) PetscRealPart(intc[pfield]));
66 #endif
67   ierr = PetscFree2(intc, intn);CHKERRQ(ierr);
68   PetscFunctionReturn(0);
69 }
70 
71 /*@C
72    SNESConvergedCorrectPressure - Convergence test that adds a vector in the nullspace to make the continuum integral of the pressure field equal to zero. This is normally used only to evaluate convergence rates for the pressure accurately. The convergence test itself just mimics SNESConvergedDefault().
73 
74    Logically Collective on SNES
75 
76    Input Parameters:
77 +  snes - the SNES context
78 .  it - the iteration (0 indicates before any Newton steps)
79 .  xnorm - 2-norm of current iterate
80 .  snorm - 2-norm of current step
81 .  fnorm - 2-norm of function at current iterate
82 -  ctx   - Optional user context
83 
84    Output Parameter:
85 .  reason  - SNES_CONVERGED_ITERATING, SNES_CONVERGED_ITS, or SNES_DIVERGED_FNORM_NAN
86 
87    Notes:
88    In order to use this monitor, you must setup several PETSc structures. First fields must be added to the DM, and a PetscDS must be created with discretizations of those fields. We currently assume that the pressure field has index 1. The pressure field must have a nullspace, likely created using the DMSetNullSpaceConstructor() interface. Last we must be able to integrate the pressure over the domain, so the DM attached to the SNES must be a Plex at this time.
89 
90    Level: advanced
91 
92 .seealso: SNESConvergedDefault(), SNESSetConvergenceTest(), DMSetNullSpaceConstructor()
93 @*/
94 PetscErrorCode SNESConvergedCorrectPressure(SNES snes, PetscInt it, PetscReal xnorm, PetscReal gnorm, PetscReal f, SNESConvergedReason *reason, void *ctx)
95 {
96   PetscBool      monitorIntegral = PETSC_FALSE;
97   PetscErrorCode ierr;
98 
99   PetscFunctionBegin;
100   ierr = SNESConvergedDefault(snes, it, xnorm, gnorm, f, reason, ctx);CHKERRQ(ierr);
101   if (monitorIntegral) {
102     Mat          J;
103     Vec          u;
104     MatNullSpace nullspace;
105     const Vec   *nullvecs;
106     PetscScalar  pintd;
107 
108     ierr = SNESGetSolution(snes, &u);CHKERRQ(ierr);
109     ierr = SNESGetJacobian(snes, &J, NULL, NULL, NULL);CHKERRQ(ierr);
110     ierr = MatGetNullSpace(J, &nullspace);CHKERRQ(ierr);
111     ierr = MatNullSpaceGetVecs(nullspace, NULL, NULL, &nullvecs);CHKERRQ(ierr);
112     ierr = VecDot(nullvecs[0], u, &pintd);CHKERRQ(ierr);
113     ierr = PetscInfo1(snes, "SNES: Discrete integral of pressure: %g\n", (double) PetscRealPart(pintd));CHKERRQ(ierr);
114   }
115   if (*reason > 0) {
116     Mat          J;
117     Vec          u;
118     MatNullSpace nullspace;
119     PetscInt     pfield = 1;
120 
121     ierr = SNESGetSolution(snes, &u);CHKERRQ(ierr);
122     ierr = SNESGetJacobian(snes, &J, NULL, NULL, NULL);CHKERRQ(ierr);
123     ierr = MatGetNullSpace(J, &nullspace);CHKERRQ(ierr);
124     ierr = SNESCorrectDiscretePressure_Private(snes, pfield, nullspace, u, ctx);CHKERRQ(ierr);
125   }
126   PetscFunctionReturn(0);
127 }
128 
129 /************************** Interpolation *******************************/
130 
131 static PetscErrorCode DMSNESConvertPlex(DM dm, DM *plex, PetscBool copy)
132 {
133   PetscBool      isPlex;
134   PetscErrorCode ierr;
135 
136   PetscFunctionBegin;
137   ierr = PetscObjectTypeCompare((PetscObject) dm, DMPLEX, &isPlex);CHKERRQ(ierr);
138   if (isPlex) {
139     *plex = dm;
140     ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);
141   } else {
142     ierr = PetscObjectQuery((PetscObject) dm, "dm_plex", (PetscObject *) plex);CHKERRQ(ierr);
143     if (!*plex) {
144       ierr = DMConvert(dm,DMPLEX,plex);CHKERRQ(ierr);
145       ierr = PetscObjectCompose((PetscObject) dm, "dm_plex", (PetscObject) *plex);CHKERRQ(ierr);
146       if (copy) {
147         ierr = DMCopyDMSNES(dm, *plex);CHKERRQ(ierr);
148         ierr = DMCopyAuxiliaryVec(dm, *plex);CHKERRQ(ierr);
149       }
150     } else {
151       ierr = PetscObjectReference((PetscObject) *plex);CHKERRQ(ierr);
152     }
153   }
154   PetscFunctionReturn(0);
155 }
156 
157 /*@C
158   DMInterpolationCreate - Creates a DMInterpolationInfo context
159 
160   Collective
161 
162   Input Parameter:
163 . comm - the communicator
164 
165   Output Parameter:
166 . ctx - the context
167 
168   Level: beginner
169 
170 .seealso: DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationDestroy()
171 @*/
172 PetscErrorCode DMInterpolationCreate(MPI_Comm comm, DMInterpolationInfo *ctx)
173 {
174   PetscErrorCode ierr;
175 
176   PetscFunctionBegin;
177   PetscValidPointer(ctx, 2);
178   ierr = PetscNew(ctx);CHKERRQ(ierr);
179 
180   (*ctx)->comm   = comm;
181   (*ctx)->dim    = -1;
182   (*ctx)->nInput = 0;
183   (*ctx)->points = NULL;
184   (*ctx)->cells  = NULL;
185   (*ctx)->n      = -1;
186   (*ctx)->coords = NULL;
187   PetscFunctionReturn(0);
188 }
189 
190 /*@C
191   DMInterpolationSetDim - Sets the spatial dimension for the interpolation context
192 
193   Not collective
194 
195   Input Parameters:
196 + ctx - the context
197 - dim - the spatial dimension
198 
199   Level: intermediate
200 
201 .seealso: DMInterpolationGetDim(), DMInterpolationEvaluate(), DMInterpolationAddPoints()
202 @*/
203 PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo ctx, PetscInt dim)
204 {
205   PetscFunctionBegin;
206   if ((dim < 1) || (dim > 3)) SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension for points: %D", dim);
207   ctx->dim = dim;
208   PetscFunctionReturn(0);
209 }
210 
211 /*@C
212   DMInterpolationGetDim - Gets the spatial dimension for the interpolation context
213 
214   Not collective
215 
216   Input Parameter:
217 . ctx - the context
218 
219   Output Parameter:
220 . dim - the spatial dimension
221 
222   Level: intermediate
223 
224 .seealso: DMInterpolationSetDim(), DMInterpolationEvaluate(), DMInterpolationAddPoints()
225 @*/
226 PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo ctx, PetscInt *dim)
227 {
228   PetscFunctionBegin;
229   PetscValidIntPointer(dim, 2);
230   *dim = ctx->dim;
231   PetscFunctionReturn(0);
232 }
233 
234 /*@C
235   DMInterpolationSetDof - Sets the number of fields interpolated at a point for the interpolation context
236 
237   Not collective
238 
239   Input Parameters:
240 + ctx - the context
241 - dof - the number of fields
242 
243   Level: intermediate
244 
245 .seealso: DMInterpolationGetDof(), DMInterpolationEvaluate(), DMInterpolationAddPoints()
246 @*/
247 PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo ctx, PetscInt dof)
248 {
249   PetscFunctionBegin;
250   if (dof < 1) SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid number of components: %D", dof);
251   ctx->dof = dof;
252   PetscFunctionReturn(0);
253 }
254 
255 /*@C
256   DMInterpolationGetDof - Gets the number of fields interpolated at a point for the interpolation context
257 
258   Not collective
259 
260   Input Parameter:
261 . ctx - the context
262 
263   Output Parameter:
264 . dof - the number of fields
265 
266   Level: intermediate
267 
268 .seealso: DMInterpolationSetDof(), DMInterpolationEvaluate(), DMInterpolationAddPoints()
269 @*/
270 PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo ctx, PetscInt *dof)
271 {
272   PetscFunctionBegin;
273   PetscValidIntPointer(dof, 2);
274   *dof = ctx->dof;
275   PetscFunctionReturn(0);
276 }
277 
278 /*@C
279   DMInterpolationAddPoints - Add points at which we will interpolate the fields
280 
281   Not collective
282 
283   Input Parameters:
284 + ctx    - the context
285 . n      - the number of points
286 - points - the coordinates for each point, an array of size n * dim
287 
288   Note: The coordinate information is copied.
289 
290   Level: intermediate
291 
292 .seealso: DMInterpolationSetDim(), DMInterpolationEvaluate(), DMInterpolationCreate()
293 @*/
294 PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo ctx, PetscInt n, PetscReal points[])
295 {
296   PetscErrorCode ierr;
297 
298   PetscFunctionBegin;
299   if (ctx->dim < 0) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set");
300   if (ctx->points)  SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot add points multiple times yet");
301   ctx->nInput = n;
302 
303   ierr = PetscMalloc1(n*ctx->dim, &ctx->points);CHKERRQ(ierr);
304   ierr = PetscArraycpy(ctx->points, points, n*ctx->dim);CHKERRQ(ierr);
305   PetscFunctionReturn(0);
306 }
307 
308 /*@C
309   DMInterpolationSetUp - Compute spatial indices for point location during interpolation
310 
311   Collective on ctx
312 
313   Input Parameters:
314 + ctx - the context
315 . dm  - the DM for the function space used for interpolation
316 . redundantPoints - If PETSC_TRUE, all processes are passing in the same array of points. Otherwise, points need to be communicated among processes.
317 - ignoreOutsideDomain - If PETSC_TRUE, ignore points outside the domain, otherwise return an error
318 
319   Level: intermediate
320 
321 .seealso: DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
322 @*/
323 PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo ctx, DM dm, PetscBool redundantPoints, PetscBool ignoreOutsideDomain)
324 {
325   MPI_Comm          comm = ctx->comm;
326   PetscScalar       *a;
327   PetscInt          p, q, i;
328   PetscMPIInt       rank, size;
329   PetscErrorCode    ierr;
330   Vec               pointVec;
331   PetscSF           cellSF;
332   PetscLayout       layout;
333   PetscReal         *globalPoints;
334   PetscScalar       *globalPointsScalar;
335   const PetscInt    *ranges;
336   PetscMPIInt       *counts, *displs;
337   const PetscSFNode *foundCells;
338   const PetscInt    *foundPoints;
339   PetscMPIInt       *foundProcs, *globalProcs;
340   PetscInt          n, N, numFound;
341 
342   PetscFunctionBegin;
343   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
344   ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr);
345   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
346   if (ctx->dim < 0) SETERRQ(comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set");
347   /* Locate points */
348   n = ctx->nInput;
349   if (!redundantPoints) {
350     ierr = PetscLayoutCreate(comm, &layout);CHKERRQ(ierr);
351     ierr = PetscLayoutSetBlockSize(layout, 1);CHKERRQ(ierr);
352     ierr = PetscLayoutSetLocalSize(layout, n);CHKERRQ(ierr);
353     ierr = PetscLayoutSetUp(layout);CHKERRQ(ierr);
354     ierr = PetscLayoutGetSize(layout, &N);CHKERRQ(ierr);
355     /* Communicate all points to all processes */
356     ierr = PetscMalloc3(N*ctx->dim,&globalPoints,size,&counts,size,&displs);CHKERRQ(ierr);
357     ierr = PetscLayoutGetRanges(layout, &ranges);CHKERRQ(ierr);
358     for (p = 0; p < size; ++p) {
359       counts[p] = (ranges[p+1] - ranges[p])*ctx->dim;
360       displs[p] = ranges[p]*ctx->dim;
361     }
362     ierr = MPI_Allgatherv(ctx->points, n*ctx->dim, MPIU_REAL, globalPoints, counts, displs, MPIU_REAL, comm);CHKERRMPI(ierr);
363   } else {
364     N = n;
365     globalPoints = ctx->points;
366     counts = displs = NULL;
367     layout = NULL;
368   }
369 #if 0
370   ierr = PetscMalloc3(N,&foundCells,N,&foundProcs,N,&globalProcs);CHKERRQ(ierr);
371   /* foundCells[p] = m->locatePoint(&globalPoints[p*ctx->dim]); */
372 #else
373 #if defined(PETSC_USE_COMPLEX)
374   ierr = PetscMalloc1(N*ctx->dim,&globalPointsScalar);CHKERRQ(ierr);
375   for (i=0; i<N*ctx->dim; i++) globalPointsScalar[i] = globalPoints[i];
376 #else
377   globalPointsScalar = globalPoints;
378 #endif
379   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF, ctx->dim, N*ctx->dim, globalPointsScalar, &pointVec);CHKERRQ(ierr);
380   ierr = PetscMalloc2(N,&foundProcs,N,&globalProcs);CHKERRQ(ierr);
381   for (p = 0; p < N; ++p) {foundProcs[p] = size;}
382   cellSF = NULL;
383   ierr = DMLocatePoints(dm, pointVec, DM_POINTLOCATION_REMOVE, &cellSF);CHKERRQ(ierr);
384   ierr = PetscSFGetGraph(cellSF,NULL,&numFound,&foundPoints,&foundCells);CHKERRQ(ierr);
385 #endif
386   for (p = 0; p < numFound; ++p) {
387     if (foundCells[p].index >= 0) foundProcs[foundPoints ? foundPoints[p] : p] = rank;
388   }
389   /* Let the lowest rank process own each point */
390   ierr   = MPIU_Allreduce(foundProcs, globalProcs, N, MPI_INT, MPI_MIN, comm);CHKERRMPI(ierr);
391   ctx->n = 0;
392   for (p = 0; p < N; ++p) {
393     if (globalProcs[p] == size) {
394       if (!ignoreOutsideDomain) SETERRQ4(comm, PETSC_ERR_PLIB, "Point %d: %g %g %g not located in mesh", p, (double)globalPoints[p*ctx->dim+0], (double)(ctx->dim > 1 ? globalPoints[p*ctx->dim+1] : 0.0), (double)(ctx->dim > 2 ? globalPoints[p*ctx->dim+2] : 0.0));
395       else if (!rank) ++ctx->n;
396     } else if (globalProcs[p] == rank) ++ctx->n;
397   }
398   /* Create coordinates vector and array of owned cells */
399   ierr = PetscMalloc1(ctx->n, &ctx->cells);CHKERRQ(ierr);
400   ierr = VecCreate(comm, &ctx->coords);CHKERRQ(ierr);
401   ierr = VecSetSizes(ctx->coords, ctx->n*ctx->dim, PETSC_DECIDE);CHKERRQ(ierr);
402   ierr = VecSetBlockSize(ctx->coords, ctx->dim);CHKERRQ(ierr);
403   ierr = VecSetType(ctx->coords,VECSTANDARD);CHKERRQ(ierr);
404   ierr = VecGetArray(ctx->coords, &a);CHKERRQ(ierr);
405   for (p = 0, q = 0, i = 0; p < N; ++p) {
406     if (globalProcs[p] == rank) {
407       PetscInt d;
408 
409       for (d = 0; d < ctx->dim; ++d, ++i) a[i] = globalPoints[p*ctx->dim+d];
410       ctx->cells[q] = foundCells[q].index;
411       ++q;
412     }
413     if (globalProcs[p] == size && !rank) {
414       PetscInt d;
415 
416       for (d = 0; d < ctx->dim; ++d, ++i) a[i] = 0.;
417       ctx->cells[q] = -1;
418       ++q;
419     }
420   }
421   ierr = VecRestoreArray(ctx->coords, &a);CHKERRQ(ierr);
422 #if 0
423   ierr = PetscFree3(foundCells,foundProcs,globalProcs);CHKERRQ(ierr);
424 #else
425   ierr = PetscFree2(foundProcs,globalProcs);CHKERRQ(ierr);
426   ierr = PetscSFDestroy(&cellSF);CHKERRQ(ierr);
427   ierr = VecDestroy(&pointVec);CHKERRQ(ierr);
428 #endif
429   if ((void*)globalPointsScalar != (void*)globalPoints) {ierr = PetscFree(globalPointsScalar);CHKERRQ(ierr);}
430   if (!redundantPoints) {ierr = PetscFree3(globalPoints,counts,displs);CHKERRQ(ierr);}
431   ierr = PetscLayoutDestroy(&layout);CHKERRQ(ierr);
432   PetscFunctionReturn(0);
433 }
434 
435 /*@C
436   DMInterpolationGetCoordinates - Gets a Vec with the coordinates of each interpolation point
437 
438   Collective on ctx
439 
440   Input Parameter:
441 . ctx - the context
442 
443   Output Parameter:
444 . coordinates  - the coordinates of interpolation points
445 
446   Note: The local vector entries correspond to interpolation points lying on this process, according to the associated DM. This is a borrowed vector that the user should not destroy.
447 
448   Level: intermediate
449 
450 .seealso: DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
451 @*/
452 PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo ctx, Vec *coordinates)
453 {
454   PetscFunctionBegin;
455   PetscValidPointer(coordinates, 2);
456   if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
457   *coordinates = ctx->coords;
458   PetscFunctionReturn(0);
459 }
460 
461 /*@C
462   DMInterpolationGetVector - Gets a Vec which can hold all the interpolated field values
463 
464   Collective on ctx
465 
466   Input Parameter:
467 . ctx - the context
468 
469   Output Parameter:
470 . v  - a vector capable of holding the interpolated field values
471 
472   Note: This vector should be returned using DMInterpolationRestoreVector().
473 
474   Level: intermediate
475 
476 .seealso: DMInterpolationRestoreVector(), DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
477 @*/
478 PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo ctx, Vec *v)
479 {
480   PetscErrorCode ierr;
481 
482   PetscFunctionBegin;
483   PetscValidPointer(v, 2);
484   if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
485   ierr = VecCreate(ctx->comm, v);CHKERRQ(ierr);
486   ierr = VecSetSizes(*v, ctx->n*ctx->dof, PETSC_DECIDE);CHKERRQ(ierr);
487   ierr = VecSetBlockSize(*v, ctx->dof);CHKERRQ(ierr);
488   ierr = VecSetType(*v,VECSTANDARD);CHKERRQ(ierr);
489   PetscFunctionReturn(0);
490 }
491 
492 /*@C
493   DMInterpolationRestoreVector - Returns a Vec which can hold all the interpolated field values
494 
495   Collective on ctx
496 
497   Input Parameters:
498 + ctx - the context
499 - v  - a vector capable of holding the interpolated field values
500 
501   Level: intermediate
502 
503 .seealso: DMInterpolationGetVector(), DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
504 @*/
505 PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo ctx, Vec *v)
506 {
507   PetscErrorCode ierr;
508 
509   PetscFunctionBegin;
510   PetscValidPointer(v, 2);
511   if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
512   ierr = VecDestroy(v);CHKERRQ(ierr);
513   PetscFunctionReturn(0);
514 }
515 
516 PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Triangle_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
517 {
518   PetscReal      *v0, *J, *invJ, detJ;
519   const PetscScalar *coords;
520   PetscScalar    *a;
521   PetscInt       p;
522   PetscErrorCode ierr;
523 
524   PetscFunctionBegin;
525   ierr = PetscMalloc3(ctx->dim,&v0,ctx->dim*ctx->dim,&J,ctx->dim*ctx->dim,&invJ);CHKERRQ(ierr);
526   ierr = VecGetArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
527   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
528   for (p = 0; p < ctx->n; ++p) {
529     PetscInt     c = ctx->cells[p];
530     PetscScalar *x = NULL;
531     PetscReal    xi[4];
532     PetscInt     d, f, comp;
533 
534     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
535     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D", (double)detJ, c);
536     ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr);
537     for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] = x[0*ctx->dof+comp];
538 
539     for (d = 0; d < ctx->dim; ++d) {
540       xi[d] = 0.0;
541       for (f = 0; f < ctx->dim; ++f) xi[d] += invJ[d*ctx->dim+f]*0.5*PetscRealPart(coords[p*ctx->dim+f] - v0[f]);
542       for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] += PetscRealPart(x[(d+1)*ctx->dof+comp] - x[0*ctx->dof+comp])*xi[d];
543     }
544     ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr);
545   }
546   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
547   ierr = VecRestoreArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
548   ierr = PetscFree3(v0, J, invJ);CHKERRQ(ierr);
549   PetscFunctionReturn(0);
550 }
551 
552 PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Tetrahedron_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
553 {
554   PetscReal      *v0, *J, *invJ, detJ;
555   const PetscScalar *coords;
556   PetscScalar    *a;
557   PetscInt       p;
558   PetscErrorCode ierr;
559 
560   PetscFunctionBegin;
561   ierr = PetscMalloc3(ctx->dim,&v0,ctx->dim*ctx->dim,&J,ctx->dim*ctx->dim,&invJ);CHKERRQ(ierr);
562   ierr = VecGetArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
563   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
564   for (p = 0; p < ctx->n; ++p) {
565     PetscInt       c = ctx->cells[p];
566     const PetscInt order[3] = {2, 1, 3};
567     PetscScalar   *x = NULL;
568     PetscReal      xi[4];
569     PetscInt       d, f, comp;
570 
571     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
572     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D", (double)detJ, c);
573     ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr);
574     for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] = x[0*ctx->dof+comp];
575 
576     for (d = 0; d < ctx->dim; ++d) {
577       xi[d] = 0.0;
578       for (f = 0; f < ctx->dim; ++f) xi[d] += invJ[d*ctx->dim+f]*0.5*PetscRealPart(coords[p*ctx->dim+f] - v0[f]);
579       for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] += PetscRealPart(x[order[d]*ctx->dof+comp] - x[0*ctx->dof+comp])*xi[d];
580     }
581     ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr);
582   }
583   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
584   ierr = VecRestoreArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
585   ierr = PetscFree3(v0, J, invJ);CHKERRQ(ierr);
586   PetscFunctionReturn(0);
587 }
588 
589 PETSC_STATIC_INLINE PetscErrorCode QuadMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx)
590 {
591   const PetscScalar *vertices = (const PetscScalar*) ctx;
592   const PetscScalar x0        = vertices[0];
593   const PetscScalar y0        = vertices[1];
594   const PetscScalar x1        = vertices[2];
595   const PetscScalar y1        = vertices[3];
596   const PetscScalar x2        = vertices[4];
597   const PetscScalar y2        = vertices[5];
598   const PetscScalar x3        = vertices[6];
599   const PetscScalar y3        = vertices[7];
600   const PetscScalar f_1       = x1 - x0;
601   const PetscScalar g_1       = y1 - y0;
602   const PetscScalar f_3       = x3 - x0;
603   const PetscScalar g_3       = y3 - y0;
604   const PetscScalar f_01      = x2 - x1 - x3 + x0;
605   const PetscScalar g_01      = y2 - y1 - y3 + y0;
606   const PetscScalar *ref;
607   PetscScalar       *real;
608   PetscErrorCode    ierr;
609 
610   PetscFunctionBegin;
611   ierr = VecGetArrayRead(Xref,  &ref);CHKERRQ(ierr);
612   ierr = VecGetArray(Xreal, &real);CHKERRQ(ierr);
613   {
614     const PetscScalar p0 = ref[0];
615     const PetscScalar p1 = ref[1];
616 
617     real[0] = x0 + f_1 * p0 + f_3 * p1 + f_01 * p0 * p1;
618     real[1] = y0 + g_1 * p0 + g_3 * p1 + g_01 * p0 * p1;
619   }
620   ierr = PetscLogFlops(28);CHKERRQ(ierr);
621   ierr = VecRestoreArrayRead(Xref,  &ref);CHKERRQ(ierr);
622   ierr = VecRestoreArray(Xreal, &real);CHKERRQ(ierr);
623   PetscFunctionReturn(0);
624 }
625 
626 #include <petsc/private/dmimpl.h>
627 PETSC_STATIC_INLINE PetscErrorCode QuadJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx)
628 {
629   const PetscScalar *vertices = (const PetscScalar*) ctx;
630   const PetscScalar x0        = vertices[0];
631   const PetscScalar y0        = vertices[1];
632   const PetscScalar x1        = vertices[2];
633   const PetscScalar y1        = vertices[3];
634   const PetscScalar x2        = vertices[4];
635   const PetscScalar y2        = vertices[5];
636   const PetscScalar x3        = vertices[6];
637   const PetscScalar y3        = vertices[7];
638   const PetscScalar f_01      = x2 - x1 - x3 + x0;
639   const PetscScalar g_01      = y2 - y1 - y3 + y0;
640   const PetscScalar *ref;
641   PetscErrorCode    ierr;
642 
643   PetscFunctionBegin;
644   ierr = VecGetArrayRead(Xref,  &ref);CHKERRQ(ierr);
645   {
646     const PetscScalar x       = ref[0];
647     const PetscScalar y       = ref[1];
648     const PetscInt    rows[2] = {0, 1};
649     PetscScalar       values[4];
650 
651     values[0] = (x1 - x0 + f_01*y) * 0.5; values[1] = (x3 - x0 + f_01*x) * 0.5;
652     values[2] = (y1 - y0 + g_01*y) * 0.5; values[3] = (y3 - y0 + g_01*x) * 0.5;
653     ierr      = MatSetValues(J, 2, rows, 2, rows, values, INSERT_VALUES);CHKERRQ(ierr);
654   }
655   ierr = PetscLogFlops(30);CHKERRQ(ierr);
656   ierr = VecRestoreArrayRead(Xref,  &ref);CHKERRQ(ierr);
657   ierr = MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
658   ierr = MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
659   PetscFunctionReturn(0);
660 }
661 
662 PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Quad_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
663 {
664   DM             dmCoord;
665   PetscFE        fem = NULL;
666   SNES           snes;
667   KSP            ksp;
668   PC             pc;
669   Vec            coordsLocal, r, ref, real;
670   Mat            J;
671   PetscTabulation    T;
672   const PetscScalar *coords;
673   PetscScalar    *a;
674   PetscReal      xir[2];
675   PetscInt       Nf, p;
676   const PetscInt dof = ctx->dof;
677   PetscErrorCode ierr;
678 
679   PetscFunctionBegin;
680   ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr);
681   if (Nf) {ierr = DMGetField(dm, 0, NULL, (PetscObject *) &fem);CHKERRQ(ierr);}
682   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
683   ierr = DMGetCoordinateDM(dm, &dmCoord);CHKERRQ(ierr);
684   ierr = SNESCreate(PETSC_COMM_SELF, &snes);CHKERRQ(ierr);
685   ierr = SNESSetOptionsPrefix(snes, "quad_interp_");CHKERRQ(ierr);
686   ierr = VecCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr);
687   ierr = VecSetSizes(r, 2, 2);CHKERRQ(ierr);
688   ierr = VecSetType(r,dm->vectype);CHKERRQ(ierr);
689   ierr = VecDuplicate(r, &ref);CHKERRQ(ierr);
690   ierr = VecDuplicate(r, &real);CHKERRQ(ierr);
691   ierr = MatCreate(PETSC_COMM_SELF, &J);CHKERRQ(ierr);
692   ierr = MatSetSizes(J, 2, 2, 2, 2);CHKERRQ(ierr);
693   ierr = MatSetType(J, MATSEQDENSE);CHKERRQ(ierr);
694   ierr = MatSetUp(J);CHKERRQ(ierr);
695   ierr = SNESSetFunction(snes, r, QuadMap_Private, NULL);CHKERRQ(ierr);
696   ierr = SNESSetJacobian(snes, J, J, QuadJacobian_Private, NULL);CHKERRQ(ierr);
697   ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr);
698   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
699   ierr = PCSetType(pc, PCLU);CHKERRQ(ierr);
700   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
701 
702   ierr = VecGetArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
703   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
704   ierr = PetscFECreateTabulation(fem, 1, 1, xir, 0, &T);CHKERRQ(ierr);
705   for (p = 0; p < ctx->n; ++p) {
706     PetscScalar *x = NULL, *vertices = NULL;
707     PetscScalar *xi;
708     PetscInt     c = ctx->cells[p], comp, coordSize, xSize;
709 
710     /* Can make this do all points at once */
711     ierr = DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
712     if (4*2 != coordSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %D should be %d", coordSize, 4*2);
713     ierr   = DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
714     ierr   = SNESSetFunction(snes, NULL, NULL, (void*) vertices);CHKERRQ(ierr);
715     ierr   = SNESSetJacobian(snes, NULL, NULL, NULL, (void*) vertices);CHKERRQ(ierr);
716     ierr   = VecGetArray(real, &xi);CHKERRQ(ierr);
717     xi[0]  = coords[p*ctx->dim+0];
718     xi[1]  = coords[p*ctx->dim+1];
719     ierr   = VecRestoreArray(real, &xi);CHKERRQ(ierr);
720     ierr   = SNESSolve(snes, real, ref);CHKERRQ(ierr);
721     ierr   = VecGetArray(ref, &xi);CHKERRQ(ierr);
722     xir[0] = PetscRealPart(xi[0]);
723     xir[1] = PetscRealPart(xi[1]);
724     if (4*dof != xSize) {
725       PetscInt d;
726 
727       xir[0] = 2.0*xir[0] - 1.0; xir[1] = 2.0*xir[1] - 1.0;
728       ierr = PetscFEComputeTabulation(fem, 1, xir, 0, T);CHKERRQ(ierr);
729       for (comp = 0; comp < dof; ++comp) {
730         a[p*dof+comp] = 0.0;
731         for (d = 0; d < xSize/dof; ++d) {
732           a[p*dof+comp] += x[d*dof+comp]*T->T[0][d*dof+comp];
733         }
734       }
735     } else {
736       for (comp = 0; comp < dof; ++comp)
737         a[p*dof+comp] = x[0*dof+comp]*(1 - xir[0])*(1 - xir[1]) + x[1*dof+comp]*xir[0]*(1 - xir[1]) + x[2*dof+comp]*xir[0]*xir[1] + x[3*dof+comp]*(1 - xir[0])*xir[1];
738     }
739     ierr = VecRestoreArray(ref, &xi);CHKERRQ(ierr);
740     ierr = DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
741     ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
742   }
743   ierr = PetscTabulationDestroy(&T);CHKERRQ(ierr);
744   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
745   ierr = VecRestoreArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
746 
747   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
748   ierr = VecDestroy(&r);CHKERRQ(ierr);
749   ierr = VecDestroy(&ref);CHKERRQ(ierr);
750   ierr = VecDestroy(&real);CHKERRQ(ierr);
751   ierr = MatDestroy(&J);CHKERRQ(ierr);
752   PetscFunctionReturn(0);
753 }
754 
755 PETSC_STATIC_INLINE PetscErrorCode HexMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx)
756 {
757   const PetscScalar *vertices = (const PetscScalar*) ctx;
758   const PetscScalar x0        = vertices[0];
759   const PetscScalar y0        = vertices[1];
760   const PetscScalar z0        = vertices[2];
761   const PetscScalar x1        = vertices[9];
762   const PetscScalar y1        = vertices[10];
763   const PetscScalar z1        = vertices[11];
764   const PetscScalar x2        = vertices[6];
765   const PetscScalar y2        = vertices[7];
766   const PetscScalar z2        = vertices[8];
767   const PetscScalar x3        = vertices[3];
768   const PetscScalar y3        = vertices[4];
769   const PetscScalar z3        = vertices[5];
770   const PetscScalar x4        = vertices[12];
771   const PetscScalar y4        = vertices[13];
772   const PetscScalar z4        = vertices[14];
773   const PetscScalar x5        = vertices[15];
774   const PetscScalar y5        = vertices[16];
775   const PetscScalar z5        = vertices[17];
776   const PetscScalar x6        = vertices[18];
777   const PetscScalar y6        = vertices[19];
778   const PetscScalar z6        = vertices[20];
779   const PetscScalar x7        = vertices[21];
780   const PetscScalar y7        = vertices[22];
781   const PetscScalar z7        = vertices[23];
782   const PetscScalar f_1       = x1 - x0;
783   const PetscScalar g_1       = y1 - y0;
784   const PetscScalar h_1       = z1 - z0;
785   const PetscScalar f_3       = x3 - x0;
786   const PetscScalar g_3       = y3 - y0;
787   const PetscScalar h_3       = z3 - z0;
788   const PetscScalar f_4       = x4 - x0;
789   const PetscScalar g_4       = y4 - y0;
790   const PetscScalar h_4       = z4 - z0;
791   const PetscScalar f_01      = x2 - x1 - x3 + x0;
792   const PetscScalar g_01      = y2 - y1 - y3 + y0;
793   const PetscScalar h_01      = z2 - z1 - z3 + z0;
794   const PetscScalar f_12      = x7 - x3 - x4 + x0;
795   const PetscScalar g_12      = y7 - y3 - y4 + y0;
796   const PetscScalar h_12      = z7 - z3 - z4 + z0;
797   const PetscScalar f_02      = x5 - x1 - x4 + x0;
798   const PetscScalar g_02      = y5 - y1 - y4 + y0;
799   const PetscScalar h_02      = z5 - z1 - z4 + z0;
800   const PetscScalar f_012     = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7;
801   const PetscScalar g_012     = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7;
802   const PetscScalar h_012     = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7;
803   const PetscScalar *ref;
804   PetscScalar       *real;
805   PetscErrorCode    ierr;
806 
807   PetscFunctionBegin;
808   ierr = VecGetArrayRead(Xref,  &ref);CHKERRQ(ierr);
809   ierr = VecGetArray(Xreal, &real);CHKERRQ(ierr);
810   {
811     const PetscScalar p0 = ref[0];
812     const PetscScalar p1 = ref[1];
813     const PetscScalar p2 = ref[2];
814 
815     real[0] = x0 + f_1*p0 + f_3*p1 + f_4*p2 + f_01*p0*p1 + f_12*p1*p2 + f_02*p0*p2 + f_012*p0*p1*p2;
816     real[1] = y0 + g_1*p0 + g_3*p1 + g_4*p2 + g_01*p0*p1 + g_01*p0*p1 + g_12*p1*p2 + g_02*p0*p2 + g_012*p0*p1*p2;
817     real[2] = z0 + h_1*p0 + h_3*p1 + h_4*p2 + h_01*p0*p1 + h_01*p0*p1 + h_12*p1*p2 + h_02*p0*p2 + h_012*p0*p1*p2;
818   }
819   ierr = PetscLogFlops(114);CHKERRQ(ierr);
820   ierr = VecRestoreArrayRead(Xref,  &ref);CHKERRQ(ierr);
821   ierr = VecRestoreArray(Xreal, &real);CHKERRQ(ierr);
822   PetscFunctionReturn(0);
823 }
824 
825 PETSC_STATIC_INLINE PetscErrorCode HexJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx)
826 {
827   const PetscScalar *vertices = (const PetscScalar*) ctx;
828   const PetscScalar x0        = vertices[0];
829   const PetscScalar y0        = vertices[1];
830   const PetscScalar z0        = vertices[2];
831   const PetscScalar x1        = vertices[9];
832   const PetscScalar y1        = vertices[10];
833   const PetscScalar z1        = vertices[11];
834   const PetscScalar x2        = vertices[6];
835   const PetscScalar y2        = vertices[7];
836   const PetscScalar z2        = vertices[8];
837   const PetscScalar x3        = vertices[3];
838   const PetscScalar y3        = vertices[4];
839   const PetscScalar z3        = vertices[5];
840   const PetscScalar x4        = vertices[12];
841   const PetscScalar y4        = vertices[13];
842   const PetscScalar z4        = vertices[14];
843   const PetscScalar x5        = vertices[15];
844   const PetscScalar y5        = vertices[16];
845   const PetscScalar z5        = vertices[17];
846   const PetscScalar x6        = vertices[18];
847   const PetscScalar y6        = vertices[19];
848   const PetscScalar z6        = vertices[20];
849   const PetscScalar x7        = vertices[21];
850   const PetscScalar y7        = vertices[22];
851   const PetscScalar z7        = vertices[23];
852   const PetscScalar f_xy      = x2 - x1 - x3 + x0;
853   const PetscScalar g_xy      = y2 - y1 - y3 + y0;
854   const PetscScalar h_xy      = z2 - z1 - z3 + z0;
855   const PetscScalar f_yz      = x7 - x3 - x4 + x0;
856   const PetscScalar g_yz      = y7 - y3 - y4 + y0;
857   const PetscScalar h_yz      = z7 - z3 - z4 + z0;
858   const PetscScalar f_xz      = x5 - x1 - x4 + x0;
859   const PetscScalar g_xz      = y5 - y1 - y4 + y0;
860   const PetscScalar h_xz      = z5 - z1 - z4 + z0;
861   const PetscScalar f_xyz     = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7;
862   const PetscScalar g_xyz     = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7;
863   const PetscScalar h_xyz     = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7;
864   const PetscScalar *ref;
865   PetscErrorCode    ierr;
866 
867   PetscFunctionBegin;
868   ierr = VecGetArrayRead(Xref,  &ref);CHKERRQ(ierr);
869   {
870     const PetscScalar x       = ref[0];
871     const PetscScalar y       = ref[1];
872     const PetscScalar z       = ref[2];
873     const PetscInt    rows[3] = {0, 1, 2};
874     PetscScalar       values[9];
875 
876     values[0] = (x1 - x0 + f_xy*y + f_xz*z + f_xyz*y*z) / 2.0;
877     values[1] = (x3 - x0 + f_xy*x + f_yz*z + f_xyz*x*z) / 2.0;
878     values[2] = (x4 - x0 + f_yz*y + f_xz*x + f_xyz*x*y) / 2.0;
879     values[3] = (y1 - y0 + g_xy*y + g_xz*z + g_xyz*y*z) / 2.0;
880     values[4] = (y3 - y0 + g_xy*x + g_yz*z + g_xyz*x*z) / 2.0;
881     values[5] = (y4 - y0 + g_yz*y + g_xz*x + g_xyz*x*y) / 2.0;
882     values[6] = (z1 - z0 + h_xy*y + h_xz*z + h_xyz*y*z) / 2.0;
883     values[7] = (z3 - z0 + h_xy*x + h_yz*z + h_xyz*x*z) / 2.0;
884     values[8] = (z4 - z0 + h_yz*y + h_xz*x + h_xyz*x*y) / 2.0;
885 
886     ierr = MatSetValues(J, 3, rows, 3, rows, values, INSERT_VALUES);CHKERRQ(ierr);
887   }
888   ierr = PetscLogFlops(152);CHKERRQ(ierr);
889   ierr = VecRestoreArrayRead(Xref,  &ref);CHKERRQ(ierr);
890   ierr = MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
891   ierr = MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
892   PetscFunctionReturn(0);
893 }
894 
895 PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Hex_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
896 {
897   DM             dmCoord;
898   SNES           snes;
899   KSP            ksp;
900   PC             pc;
901   Vec            coordsLocal, r, ref, real;
902   Mat            J;
903   const PetscScalar *coords;
904   PetscScalar    *a;
905   PetscInt       p;
906   PetscErrorCode ierr;
907 
908   PetscFunctionBegin;
909   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
910   ierr = DMGetCoordinateDM(dm, &dmCoord);CHKERRQ(ierr);
911   ierr = SNESCreate(PETSC_COMM_SELF, &snes);CHKERRQ(ierr);
912   ierr = SNESSetOptionsPrefix(snes, "hex_interp_");CHKERRQ(ierr);
913   ierr = VecCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr);
914   ierr = VecSetSizes(r, 3, 3);CHKERRQ(ierr);
915   ierr = VecSetType(r,dm->vectype);CHKERRQ(ierr);
916   ierr = VecDuplicate(r, &ref);CHKERRQ(ierr);
917   ierr = VecDuplicate(r, &real);CHKERRQ(ierr);
918   ierr = MatCreate(PETSC_COMM_SELF, &J);CHKERRQ(ierr);
919   ierr = MatSetSizes(J, 3, 3, 3, 3);CHKERRQ(ierr);
920   ierr = MatSetType(J, MATSEQDENSE);CHKERRQ(ierr);
921   ierr = MatSetUp(J);CHKERRQ(ierr);
922   ierr = SNESSetFunction(snes, r, HexMap_Private, NULL);CHKERRQ(ierr);
923   ierr = SNESSetJacobian(snes, J, J, HexJacobian_Private, NULL);CHKERRQ(ierr);
924   ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr);
925   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
926   ierr = PCSetType(pc, PCLU);CHKERRQ(ierr);
927   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
928 
929   ierr = VecGetArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
930   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
931   for (p = 0; p < ctx->n; ++p) {
932     PetscScalar *x = NULL, *vertices = NULL;
933     PetscScalar *xi;
934     PetscReal    xir[3];
935     PetscInt     c = ctx->cells[p], comp, coordSize, xSize;
936 
937     /* Can make this do all points at once */
938     ierr = DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
939     if (8*3 != coordSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %D should be %d", coordSize, 8*3);
940     ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
941     if (8*ctx->dof != xSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %D should be %D", xSize, 8*ctx->dof);
942     ierr   = SNESSetFunction(snes, NULL, NULL, (void*) vertices);CHKERRQ(ierr);
943     ierr   = SNESSetJacobian(snes, NULL, NULL, NULL, (void*) vertices);CHKERRQ(ierr);
944     ierr   = VecGetArray(real, &xi);CHKERRQ(ierr);
945     xi[0]  = coords[p*ctx->dim+0];
946     xi[1]  = coords[p*ctx->dim+1];
947     xi[2]  = coords[p*ctx->dim+2];
948     ierr   = VecRestoreArray(real, &xi);CHKERRQ(ierr);
949     ierr   = SNESSolve(snes, real, ref);CHKERRQ(ierr);
950     ierr   = VecGetArray(ref, &xi);CHKERRQ(ierr);
951     xir[0] = PetscRealPart(xi[0]);
952     xir[1] = PetscRealPart(xi[1]);
953     xir[2] = PetscRealPart(xi[2]);
954     for (comp = 0; comp < ctx->dof; ++comp) {
955       a[p*ctx->dof+comp] =
956         x[0*ctx->dof+comp]*(1-xir[0])*(1-xir[1])*(1-xir[2]) +
957         x[3*ctx->dof+comp]*    xir[0]*(1-xir[1])*(1-xir[2]) +
958         x[2*ctx->dof+comp]*    xir[0]*    xir[1]*(1-xir[2]) +
959         x[1*ctx->dof+comp]*(1-xir[0])*    xir[1]*(1-xir[2]) +
960         x[4*ctx->dof+comp]*(1-xir[0])*(1-xir[1])*   xir[2] +
961         x[5*ctx->dof+comp]*    xir[0]*(1-xir[1])*   xir[2] +
962         x[6*ctx->dof+comp]*    xir[0]*    xir[1]*   xir[2] +
963         x[7*ctx->dof+comp]*(1-xir[0])*    xir[1]*   xir[2];
964     }
965     ierr = VecRestoreArray(ref, &xi);CHKERRQ(ierr);
966     ierr = DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
967     ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
968   }
969   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
970   ierr = VecRestoreArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
971 
972   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
973   ierr = VecDestroy(&r);CHKERRQ(ierr);
974   ierr = VecDestroy(&ref);CHKERRQ(ierr);
975   ierr = VecDestroy(&real);CHKERRQ(ierr);
976   ierr = MatDestroy(&J);CHKERRQ(ierr);
977   PetscFunctionReturn(0);
978 }
979 
980 /*@C
981   DMInterpolationEvaluate - Using the input from dm and x, calculates interpolated field values at the interpolation points.
982 
983   Input Parameters:
984 + ctx - The DMInterpolationInfo context
985 . dm  - The DM
986 - x   - The local vector containing the field to be interpolated
987 
988   Output Parameters:
989 . v   - The vector containing the interpolated values
990 
991   Note: A suitable v can be obtained using DMInterpolationGetVector().
992 
993   Level: beginner
994 
995 .seealso: DMInterpolationGetVector(), DMInterpolationAddPoints(), DMInterpolationCreate()
996 @*/
997 PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo ctx, DM dm, Vec x, Vec v)
998 {
999   PetscInt       n;
1000   PetscErrorCode ierr;
1001 
1002   PetscFunctionBegin;
1003   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1004   PetscValidHeaderSpecific(x, VEC_CLASSID, 3);
1005   PetscValidHeaderSpecific(v, VEC_CLASSID, 4);
1006   ierr = VecGetLocalSize(v, &n);CHKERRQ(ierr);
1007   if (n != ctx->n*ctx->dof) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid input vector size %D should be %D", n, ctx->n*ctx->dof);
1008   if (n) {
1009     PetscDS        ds;
1010     DMPolytopeType ct;
1011     PetscBool      done = PETSC_FALSE;
1012 
1013     ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
1014     if (ds) {
1015       const PetscScalar *coords;
1016       PetscScalar       *interpolant;
1017       PetscInt           cdim, d, p, Nf, field, c = 0;
1018 
1019       ierr = DMGetCoordinateDim(dm, &cdim);CHKERRQ(ierr);
1020       ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
1021       for (field = 0; field < Nf; ++field) {
1022         PetscTabulation T;
1023         PetscFE         fe;
1024         PetscClassId    id;
1025         PetscReal       xi[3];
1026         PetscInt        Nc, f, fc;
1027 
1028         ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr);
1029         ierr = PetscObjectGetClassId((PetscObject) fe, &id);CHKERRQ(ierr);
1030         if (id != PETSCFE_CLASSID) break;
1031         ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1032         ierr = VecGetArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
1033         ierr = VecGetArrayWrite(v, &interpolant);CHKERRQ(ierr);
1034         for (p = 0; p < ctx->n; ++p) {
1035           PetscScalar *xa = NULL;
1036           PetscReal    pcoords[3];
1037 
1038           if (ctx->cells[p] < 0) continue;
1039           for (d = 0; d < cdim; ++d) pcoords[d] = PetscRealPart(coords[p*cdim+d]);
1040           ierr = DMPlexCoordinatesToReference(dm, ctx->cells[p], 1, pcoords, xi);CHKERRQ(ierr);
1041           ierr = DMPlexVecGetClosure(dm, NULL, x, ctx->cells[p], NULL, &xa);CHKERRQ(ierr);
1042           ierr = PetscFECreateTabulation(fe, 1, 1, xi, 0, &T);CHKERRQ(ierr);
1043           {
1044             const PetscReal *basis = T->T[0];
1045             const PetscInt   Nb    = T->Nb;
1046             const PetscInt   Nc    = T->Nc;
1047             for (fc = 0; fc < Nc; ++fc) {
1048               interpolant[p*ctx->dof+c+fc] = 0.0;
1049               for (f = 0; f < Nb; ++f) {
1050                 interpolant[p*ctx->dof+c+fc] += xa[f]*basis[(0*Nb + f)*Nc + fc];
1051               }
1052             }
1053           }
1054           ierr = PetscTabulationDestroy(&T);CHKERRQ(ierr);
1055           ierr = DMPlexVecRestoreClosure(dm, NULL, x, ctx->cells[p], NULL, &xa);CHKERRQ(ierr);
1056         }
1057         ierr = VecRestoreArrayWrite(v, &interpolant);CHKERRQ(ierr);
1058         ierr = VecRestoreArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
1059         c += Nc;
1060       }
1061       if (field == Nf) {
1062         done = PETSC_TRUE;
1063         if (c != ctx->dof) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Total components %D != %D dof specified for interpolation", c, ctx->dof);
1064       }
1065     }
1066     if (!done) {
1067       /* TODO Check each cell individually */
1068       ierr = DMPlexGetCellType(dm, ctx->cells[0], &ct);CHKERRQ(ierr);
1069       switch (ct) {
1070         case DM_POLYTOPE_TRIANGLE:      ierr = DMInterpolate_Triangle_Private(ctx, dm, x, v);CHKERRQ(ierr);break;
1071         case DM_POLYTOPE_QUADRILATERAL: ierr = DMInterpolate_Quad_Private(ctx, dm, x, v);CHKERRQ(ierr);break;
1072         case DM_POLYTOPE_TETRAHEDRON:   ierr = DMInterpolate_Tetrahedron_Private(ctx, dm, x, v);CHKERRQ(ierr);break;
1073         case DM_POLYTOPE_HEXAHEDRON:    ierr = DMInterpolate_Hex_Private(ctx, dm, x, v);CHKERRQ(ierr);break;
1074         default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support fpr cell type %s", DMPolytopeTypes[ct]);
1075       }
1076     }
1077   }
1078   PetscFunctionReturn(0);
1079 }
1080 
1081 /*@C
1082   DMInterpolationDestroy - Destroys a DMInterpolationInfo context
1083 
1084   Collective on ctx
1085 
1086   Input Parameter:
1087 . ctx - the context
1088 
1089   Level: beginner
1090 
1091 .seealso: DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
1092 @*/
1093 PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *ctx)
1094 {
1095   PetscErrorCode ierr;
1096 
1097   PetscFunctionBegin;
1098   PetscValidPointer(ctx, 1);
1099   ierr = VecDestroy(&(*ctx)->coords);CHKERRQ(ierr);
1100   ierr = PetscFree((*ctx)->points);CHKERRQ(ierr);
1101   ierr = PetscFree((*ctx)->cells);CHKERRQ(ierr);
1102   ierr = PetscFree(*ctx);CHKERRQ(ierr);
1103   *ctx = NULL;
1104   PetscFunctionReturn(0);
1105 }
1106 
1107 /*@C
1108   SNESMonitorFields - Monitors the residual for each field separately
1109 
1110   Collective on SNES
1111 
1112   Input Parameters:
1113 + snes   - the SNES context
1114 . its    - iteration number
1115 . fgnorm - 2-norm of residual
1116 - vf  - PetscViewerAndFormat of type ASCII
1117 
1118   Notes:
1119   This routine prints the residual norm at each iteration.
1120 
1121   Level: intermediate
1122 
1123 .seealso: SNESMonitorSet(), SNESMonitorDefault()
1124 @*/
1125 PetscErrorCode SNESMonitorFields(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf)
1126 {
1127   PetscViewer        viewer = vf->viewer;
1128   Vec                res;
1129   DM                 dm;
1130   PetscSection       s;
1131   const PetscScalar *r;
1132   PetscReal         *lnorms, *norms;
1133   PetscInt           numFields, f, pStart, pEnd, p;
1134   PetscErrorCode     ierr;
1135 
1136   PetscFunctionBegin;
1137   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4);
1138   ierr = SNESGetFunction(snes, &res, NULL, NULL);CHKERRQ(ierr);
1139   ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr);
1140   ierr = DMGetLocalSection(dm, &s);CHKERRQ(ierr);
1141   ierr = PetscSectionGetNumFields(s, &numFields);CHKERRQ(ierr);
1142   ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr);
1143   ierr = PetscCalloc2(numFields, &lnorms, numFields, &norms);CHKERRQ(ierr);
1144   ierr = VecGetArrayRead(res, &r);CHKERRQ(ierr);
1145   for (p = pStart; p < pEnd; ++p) {
1146     for (f = 0; f < numFields; ++f) {
1147       PetscInt fdof, foff, d;
1148 
1149       ierr = PetscSectionGetFieldDof(s, p, f, &fdof);CHKERRQ(ierr);
1150       ierr = PetscSectionGetFieldOffset(s, p, f, &foff);CHKERRQ(ierr);
1151       for (d = 0; d < fdof; ++d) lnorms[f] += PetscRealPart(PetscSqr(r[foff+d]));
1152     }
1153   }
1154   ierr = VecRestoreArrayRead(res, &r);CHKERRQ(ierr);
1155   ierr = MPIU_Allreduce(lnorms, norms, numFields, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject) dm));CHKERRMPI(ierr);
1156   ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr);
1157   ierr = PetscViewerASCIIAddTab(viewer, ((PetscObject) snes)->tablevel);CHKERRQ(ierr);
1158   ierr = PetscViewerASCIIPrintf(viewer, "%3D SNES Function norm %14.12e [", its, (double) fgnorm);CHKERRQ(ierr);
1159   for (f = 0; f < numFields; ++f) {
1160     if (f > 0) {ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr);}
1161     ierr = PetscViewerASCIIPrintf(viewer, "%14.12e", (double) PetscSqrtReal(norms[f]));CHKERRQ(ierr);
1162   }
1163   ierr = PetscViewerASCIIPrintf(viewer, "]\n");CHKERRQ(ierr);
1164   ierr = PetscViewerASCIISubtractTab(viewer, ((PetscObject) snes)->tablevel);CHKERRQ(ierr);
1165   ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
1166   ierr = PetscFree2(lnorms, norms);CHKERRQ(ierr);
1167   PetscFunctionReturn(0);
1168 }
1169 
1170 /********************* Residual Computation **************************/
1171 
1172 PetscErrorCode DMPlexGetAllCells_Internal(DM plex, IS *cellIS)
1173 {
1174   PetscInt       depth;
1175   PetscErrorCode ierr;
1176 
1177   PetscFunctionBegin;
1178   ierr = DMPlexGetDepth(plex, &depth);CHKERRQ(ierr);
1179   ierr = DMGetStratumIS(plex, "dim", depth, cellIS);CHKERRQ(ierr);
1180   if (!*cellIS) {ierr = DMGetStratumIS(plex, "depth", depth, cellIS);CHKERRQ(ierr);}
1181   PetscFunctionReturn(0);
1182 }
1183 
1184 /*@
1185   DMPlexSNESComputeResidualFEM - Form the local residual F from the local input X using pointwise functions specified by the user
1186 
1187   Input Parameters:
1188 + dm - The mesh
1189 . X  - Local solution
1190 - user - The user context
1191 
1192   Output Parameter:
1193 . F  - Local output vector
1194 
1195   Level: developer
1196 
1197 .seealso: DMPlexComputeJacobianAction()
1198 @*/
1199 PetscErrorCode DMPlexSNESComputeResidualFEM(DM dm, Vec X, Vec F, void *user)
1200 {
1201   DM             plex;
1202   IS             allcellIS;
1203   PetscInt       Nds, s;
1204   PetscErrorCode ierr;
1205 
1206   PetscFunctionBegin;
1207   ierr = DMSNESConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr);
1208   ierr = DMPlexGetAllCells_Internal(plex, &allcellIS);CHKERRQ(ierr);
1209   ierr = DMGetNumDS(dm, &Nds);CHKERRQ(ierr);
1210   for (s = 0; s < Nds; ++s) {
1211     PetscDS          ds;
1212     IS               cellIS;
1213     PetscFormKey key;
1214 
1215     ierr = DMGetRegionNumDS(dm, s, &key.label, NULL, &ds);CHKERRQ(ierr);
1216     key.value = 0;
1217     key.field = 0;
1218     key.part  = 0;
1219     if (!key.label) {
1220       ierr = PetscObjectReference((PetscObject) allcellIS);CHKERRQ(ierr);
1221       cellIS = allcellIS;
1222     } else {
1223       IS pointIS;
1224 
1225       key.value = 1;
1226       ierr = DMLabelGetStratumIS(key.label, key.value, &pointIS);CHKERRQ(ierr);
1227       ierr = ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS);CHKERRQ(ierr);
1228       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
1229     }
1230     ierr = DMPlexComputeResidual_Internal(plex, key, cellIS, PETSC_MIN_REAL, X, NULL, 0.0, F, user);CHKERRQ(ierr);
1231     ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
1232   }
1233   ierr = ISDestroy(&allcellIS);CHKERRQ(ierr);
1234   ierr = DMDestroy(&plex);CHKERRQ(ierr);
1235   PetscFunctionReturn(0);
1236 }
1237 
1238 PetscErrorCode DMSNESComputeResidual(DM dm, Vec X, Vec F, void *user)
1239 {
1240   DM             plex;
1241   IS             allcellIS;
1242   PetscInt       Nds, s;
1243   PetscErrorCode ierr;
1244 
1245   PetscFunctionBegin;
1246   ierr = DMSNESConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr);
1247   ierr = DMPlexGetAllCells_Internal(plex, &allcellIS);CHKERRQ(ierr);
1248   ierr = DMGetNumDS(dm, &Nds);CHKERRQ(ierr);
1249   for (s = 0; s < Nds; ++s) {
1250     PetscDS ds;
1251     DMLabel label;
1252     IS      cellIS;
1253 
1254     ierr = DMGetRegionNumDS(dm, s, &label, NULL, &ds);CHKERRQ(ierr);
1255     {
1256       PetscWeakFormKind resmap[2] = {PETSC_WF_F0, PETSC_WF_F1};
1257       PetscWeakForm     wf;
1258       PetscInt          Nm = 2, m, Nk = 0, k, kp, off = 0;
1259       PetscFormKey *reskeys;
1260 
1261       /* Get unique residual keys */
1262       for (m = 0; m < Nm; ++m) {
1263         PetscInt Nkm;
1264         ierr = PetscHMapFormGetSize(ds->wf->form[resmap[m]], &Nkm);CHKERRQ(ierr);
1265         Nk  += Nkm;
1266       }
1267       ierr = PetscMalloc1(Nk, &reskeys);CHKERRQ(ierr);
1268       for (m = 0; m < Nm; ++m) {
1269         ierr = PetscHMapFormGetKeys(ds->wf->form[resmap[m]], &off, reskeys);CHKERRQ(ierr);
1270       }
1271       if (off != Nk) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of keys %D should be %D", off, Nk);
1272       ierr = PetscFormKeySort(Nk, reskeys);CHKERRQ(ierr);
1273       for (k = 0, kp = 1; kp < Nk; ++kp) {
1274         if ((reskeys[k].label != reskeys[kp].label) || (reskeys[k].value != reskeys[kp].value)) {
1275           ++k;
1276           if (kp != k) reskeys[k] = reskeys[kp];
1277         }
1278       }
1279       Nk = k;
1280 
1281       ierr = PetscDSGetWeakForm(ds, &wf);CHKERRQ(ierr);
1282       for (k = 0; k < Nk; ++k) {
1283         DMLabel  label = reskeys[k].label;
1284         PetscInt val   = reskeys[k].value;
1285 
1286         if (!label) {
1287           ierr = PetscObjectReference((PetscObject) allcellIS);CHKERRQ(ierr);
1288           cellIS = allcellIS;
1289         } else {
1290           IS pointIS;
1291 
1292           ierr = DMLabelGetStratumIS(label, val, &pointIS);CHKERRQ(ierr);
1293           ierr = ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS);CHKERRQ(ierr);
1294           ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
1295         }
1296         ierr = DMPlexComputeResidual_Internal(plex, reskeys[k], cellIS, PETSC_MIN_REAL, X, NULL, 0.0, F, user);CHKERRQ(ierr);
1297         ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
1298       }
1299       ierr = PetscFree(reskeys);CHKERRQ(ierr);
1300     }
1301   }
1302   ierr = ISDestroy(&allcellIS);CHKERRQ(ierr);
1303   ierr = DMDestroy(&plex);CHKERRQ(ierr);
1304   PetscFunctionReturn(0);
1305 }
1306 
1307 /*@
1308   DMPlexSNESComputeBoundaryFEM - Form the boundary values for the local input X
1309 
1310   Input Parameters:
1311 + dm - The mesh
1312 - user - The user context
1313 
1314   Output Parameter:
1315 . X  - Local solution
1316 
1317   Level: developer
1318 
1319 .seealso: DMPlexComputeJacobianAction()
1320 @*/
1321 PetscErrorCode DMPlexSNESComputeBoundaryFEM(DM dm, Vec X, void *user)
1322 {
1323   DM             plex;
1324   PetscErrorCode ierr;
1325 
1326   PetscFunctionBegin;
1327   ierr = DMSNESConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
1328   ierr = DMPlexInsertBoundaryValues(plex, PETSC_TRUE, X, PETSC_MIN_REAL, NULL, NULL, NULL);CHKERRQ(ierr);
1329   ierr = DMDestroy(&plex);CHKERRQ(ierr);
1330   PetscFunctionReturn(0);
1331 }
1332 
1333 /*@
1334   DMSNESComputeJacobianAction - Compute the action of the Jacobian J(X) on Y
1335 
1336   Input Parameters:
1337 + dm   - The DM
1338 . X    - Local solution vector
1339 . Y    - Local input vector
1340 - user - The user context
1341 
1342   Output Parameter:
1343 . F    - lcoal output vector
1344 
1345   Level: developer
1346 
1347   Notes:
1348   Users will typically use DMSNESCreateJacobianMF() followed by MatMult() instead of calling this routine directly.
1349 
1350 .seealso: DMSNESCreateJacobianMF()
1351 @*/
1352 PetscErrorCode DMSNESComputeJacobianAction(DM dm, Vec X, Vec Y, Vec F, void *user)
1353 {
1354   DM             plex;
1355   IS             allcellIS;
1356   PetscInt       Nds, s;
1357   PetscErrorCode ierr;
1358 
1359   PetscFunctionBegin;
1360   ierr = DMSNESConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr);
1361   ierr = DMPlexGetAllCells_Internal(plex, &allcellIS);CHKERRQ(ierr);
1362   ierr = DMGetNumDS(dm, &Nds);CHKERRQ(ierr);
1363   for (s = 0; s < Nds; ++s) {
1364     PetscDS ds;
1365     DMLabel label;
1366     IS      cellIS;
1367 
1368     ierr = DMGetRegionNumDS(dm, s, &label, NULL, &ds);CHKERRQ(ierr);
1369     {
1370       PetscWeakFormKind jacmap[4] = {PETSC_WF_G0, PETSC_WF_G1, PETSC_WF_G2, PETSC_WF_G3};
1371       PetscWeakForm     wf;
1372       PetscInt          Nm = 4, m, Nk = 0, k, kp, off = 0;
1373       PetscFormKey *jackeys;
1374 
1375       /* Get unique Jacobian keys */
1376       for (m = 0; m < Nm; ++m) {
1377         PetscInt Nkm;
1378         ierr = PetscHMapFormGetSize(ds->wf->form[jacmap[m]], &Nkm);CHKERRQ(ierr);
1379         Nk  += Nkm;
1380       }
1381       ierr = PetscMalloc1(Nk, &jackeys);CHKERRQ(ierr);
1382       for (m = 0; m < Nm; ++m) {
1383         ierr = PetscHMapFormGetKeys(ds->wf->form[jacmap[m]], &off, jackeys);CHKERRQ(ierr);
1384       }
1385       if (off != Nk) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of keys %D should be %D", off, Nk);
1386       ierr = PetscFormKeySort(Nk, jackeys);CHKERRQ(ierr);
1387       for (k = 0, kp = 1; kp < Nk; ++kp) {
1388         if ((jackeys[k].label != jackeys[kp].label) || (jackeys[k].value != jackeys[kp].value)) {
1389           ++k;
1390           if (kp != k) jackeys[k] = jackeys[kp];
1391         }
1392       }
1393       Nk = k;
1394 
1395       ierr = PetscDSGetWeakForm(ds, &wf);CHKERRQ(ierr);
1396       for (k = 0; k < Nk; ++k) {
1397         DMLabel  label = jackeys[k].label;
1398         PetscInt val   = jackeys[k].value;
1399 
1400         if (!label) {
1401           ierr = PetscObjectReference((PetscObject) allcellIS);CHKERRQ(ierr);
1402           cellIS = allcellIS;
1403         } else {
1404           IS pointIS;
1405 
1406           ierr = DMLabelGetStratumIS(label, val, &pointIS);CHKERRQ(ierr);
1407           ierr = ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS);CHKERRQ(ierr);
1408           ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
1409         }
1410         ierr = DMPlexComputeJacobian_Action_Internal(plex, jackeys[k], cellIS, 0.0, 0.0, X, NULL, Y, F, user);CHKERRQ(ierr);
1411         ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
1412       }
1413       ierr = PetscFree(jackeys);CHKERRQ(ierr);
1414     }
1415   }
1416   ierr = ISDestroy(&allcellIS);CHKERRQ(ierr);
1417   ierr = DMDestroy(&plex);CHKERRQ(ierr);
1418   PetscFunctionReturn(0);
1419 }
1420 
1421 /*@
1422   DMPlexSNESComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user.
1423 
1424   Input Parameters:
1425 + dm - The mesh
1426 . X  - Local input vector
1427 - user - The user context
1428 
1429   Output Parameter:
1430 . Jac  - Jacobian matrix
1431 
1432   Note:
1433   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
1434   like a GPU, or vectorize on a multicore machine.
1435 
1436   Level: developer
1437 
1438 .seealso: FormFunctionLocal()
1439 @*/
1440 PetscErrorCode DMPlexSNESComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP,void *user)
1441 {
1442   DM             plex;
1443   IS             allcellIS;
1444   PetscBool      hasJac, hasPrec;
1445   PetscInt       Nds, s;
1446   PetscErrorCode ierr;
1447 
1448   PetscFunctionBegin;
1449   ierr = DMSNESConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr);
1450   ierr = DMPlexGetAllCells_Internal(plex, &allcellIS);CHKERRQ(ierr);
1451   ierr = DMGetNumDS(dm, &Nds);CHKERRQ(ierr);
1452   for (s = 0; s < Nds; ++s) {
1453     PetscDS          ds;
1454     IS               cellIS;
1455     PetscFormKey key;
1456 
1457     ierr = DMGetRegionNumDS(dm, s, &key.label, NULL, &ds);CHKERRQ(ierr);
1458     key.value = 0;
1459     key.field = 0;
1460     key.part  = 0;
1461     if (!key.label) {
1462       ierr = PetscObjectReference((PetscObject) allcellIS);CHKERRQ(ierr);
1463       cellIS = allcellIS;
1464     } else {
1465       IS pointIS;
1466 
1467       key.value = 1;
1468       ierr = DMLabelGetStratumIS(key.label, key.value, &pointIS);CHKERRQ(ierr);
1469       ierr = ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS);CHKERRQ(ierr);
1470       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
1471     }
1472     if (!s) {
1473       ierr = PetscDSHasJacobian(ds, &hasJac);CHKERRQ(ierr);
1474       ierr = PetscDSHasJacobianPreconditioner(ds, &hasPrec);CHKERRQ(ierr);
1475       if (hasJac && hasPrec) {ierr = MatZeroEntries(Jac);CHKERRQ(ierr);}
1476       ierr = MatZeroEntries(JacP);CHKERRQ(ierr);
1477     }
1478     ierr = DMPlexComputeJacobian_Internal(plex, key, cellIS, 0.0, 0.0, X, NULL, Jac, JacP, user);CHKERRQ(ierr);
1479     ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
1480   }
1481   ierr = ISDestroy(&allcellIS);CHKERRQ(ierr);
1482   ierr = DMDestroy(&plex);CHKERRQ(ierr);
1483   PetscFunctionReturn(0);
1484 }
1485 
1486 struct _DMSNESJacobianMFCtx
1487 {
1488   DM    dm;
1489   Vec   X;
1490   void *ctx;
1491 };
1492 
1493 static PetscErrorCode DMSNESJacobianMF_Destroy_Private(Mat A)
1494 {
1495   struct _DMSNESJacobianMFCtx *ctx;
1496   PetscErrorCode               ierr;
1497 
1498   PetscFunctionBegin;
1499   ierr = MatShellGetContext(A, (void **) &ctx);CHKERRQ(ierr);
1500   ierr = MatShellSetContext(A, NULL);CHKERRQ(ierr);
1501   ierr = DMDestroy(&ctx->dm);CHKERRQ(ierr);
1502   ierr = VecDestroy(&ctx->X);CHKERRQ(ierr);
1503   ierr = PetscFree(ctx);CHKERRQ(ierr);
1504   PetscFunctionReturn(0);
1505 }
1506 
1507 static PetscErrorCode DMSNESJacobianMF_Mult_Private(Mat A, Vec Y, Vec Z)
1508 {
1509   struct _DMSNESJacobianMFCtx *ctx;
1510   PetscErrorCode               ierr;
1511 
1512   PetscFunctionBegin;
1513   ierr = MatShellGetContext(A, (void **) &ctx);CHKERRQ(ierr);
1514   ierr = DMSNESComputeJacobianAction(ctx->dm, ctx->X, Y, Z, ctx->ctx);CHKERRQ(ierr);
1515   PetscFunctionReturn(0);
1516 }
1517 
1518 /*@
1519   DMSNESCreateJacobianMF - Create a Mat which computes the action of the Jacobian matrix-free
1520 
1521   Collective on dm
1522 
1523   Input Parameters:
1524 + dm   - The DM
1525 . X    - The evaluation point for the Jacobian
1526 - user - A user context, or NULL
1527 
1528   Output Parameter:
1529 . J    - The Mat
1530 
1531   Level: advanced
1532 
1533   Notes:
1534   Vec X is kept in Mat J, so updating X then updates the evaluation point.
1535 
1536 .seealso: DMSNESComputeJacobianAction()
1537 @*/
1538 PetscErrorCode DMSNESCreateJacobianMF(DM dm, Vec X, void *user, Mat *J)
1539 {
1540   struct _DMSNESJacobianMFCtx *ctx;
1541   PetscInt                     n, N;
1542   PetscErrorCode               ierr;
1543 
1544   PetscFunctionBegin;
1545   ierr = MatCreate(PetscObjectComm((PetscObject) dm), J);CHKERRQ(ierr);
1546   ierr = MatSetType(*J, MATSHELL);CHKERRQ(ierr);
1547   ierr = VecGetLocalSize(X, &n);CHKERRQ(ierr);
1548   ierr = VecGetSize(X, &N);CHKERRQ(ierr);
1549   ierr = MatSetSizes(*J, n, n, N, N);CHKERRQ(ierr);
1550   ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);
1551   ierr = PetscObjectReference((PetscObject) X);CHKERRQ(ierr);
1552   ierr = PetscMalloc1(1, &ctx);CHKERRQ(ierr);
1553   ctx->dm  = dm;
1554   ctx->X   = X;
1555   ctx->ctx = user;
1556   ierr = MatShellSetContext(*J, ctx);CHKERRQ(ierr);
1557   ierr = MatShellSetOperation(*J, MATOP_DESTROY, (void (*)(void)) DMSNESJacobianMF_Destroy_Private);CHKERRQ(ierr);
1558   ierr = MatShellSetOperation(*J, MATOP_MULT,    (void (*)(void)) DMSNESJacobianMF_Mult_Private);CHKERRQ(ierr);
1559   PetscFunctionReturn(0);
1560 }
1561 
1562 /*
1563      MatComputeNeumannOverlap - Computes an unassembled (Neumann) local overlapping Mat in nonlinear context.
1564 
1565    Input Parameters:
1566 +     X - SNES linearization point
1567 .     ovl - index set of overlapping subdomains
1568 
1569    Output Parameter:
1570 .     J - unassembled (Neumann) local matrix
1571 
1572    Level: intermediate
1573 
1574 .seealso: DMCreateNeumannOverlap(), MATIS, PCHPDDMSetAuxiliaryMat()
1575 */
1576 static PetscErrorCode MatComputeNeumannOverlap_Plex(Mat J, PetscReal t, Vec X, Vec X_t, PetscReal s, IS ovl, void *ctx)
1577 {
1578   SNES           snes;
1579   Mat            pJ;
1580   DM             ovldm,origdm;
1581   DMSNES         sdm;
1582   PetscErrorCode (*bfun)(DM,Vec,void*);
1583   PetscErrorCode (*jfun)(DM,Vec,Mat,Mat,void*);
1584   void           *bctx,*jctx;
1585   PetscErrorCode ierr;
1586 
1587   PetscFunctionBegin;
1588   ierr = PetscObjectQuery((PetscObject)ovl,"_DM_Overlap_HPDDM_MATIS",(PetscObject*)&pJ);CHKERRQ(ierr);
1589   if (!pJ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing overlapping Mat");CHKERRQ(ierr);
1590   ierr = PetscObjectQuery((PetscObject)ovl,"_DM_Original_HPDDM",(PetscObject*)&origdm);CHKERRQ(ierr);
1591   if (!origdm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing original DM");CHKERRQ(ierr);
1592   ierr = MatGetDM(pJ,&ovldm);CHKERRQ(ierr);
1593   ierr = DMSNESGetBoundaryLocal(origdm,&bfun,&bctx);CHKERRQ(ierr);
1594   ierr = DMSNESSetBoundaryLocal(ovldm,bfun,bctx);CHKERRQ(ierr);
1595   ierr = DMSNESGetJacobianLocal(origdm,&jfun,&jctx);CHKERRQ(ierr);
1596   ierr = DMSNESSetJacobianLocal(ovldm,jfun,jctx);CHKERRQ(ierr);
1597   ierr = PetscObjectQuery((PetscObject)ovl,"_DM_Overlap_HPDDM_SNES",(PetscObject*)&snes);CHKERRQ(ierr);
1598   if (!snes) {
1599     ierr = SNESCreate(PetscObjectComm((PetscObject)ovl),&snes);CHKERRQ(ierr);
1600     ierr = SNESSetDM(snes,ovldm);CHKERRQ(ierr);
1601     ierr = PetscObjectCompose((PetscObject)ovl,"_DM_Overlap_HPDDM_SNES",(PetscObject)snes);CHKERRQ(ierr);
1602     ierr = PetscObjectDereference((PetscObject)snes);CHKERRQ(ierr);
1603   }
1604   ierr = DMGetDMSNES(ovldm,&sdm);CHKERRQ(ierr);
1605   ierr = VecLockReadPush(X);CHKERRQ(ierr);
1606   PetscStackPush("SNES user Jacobian function");
1607   ierr = (*sdm->ops->computejacobian)(snes,X,pJ,pJ,sdm->jacobianctx);CHKERRQ(ierr);
1608   PetscStackPop;
1609   ierr = VecLockReadPop(X);CHKERRQ(ierr);
1610   /* this is a no-hop, just in case we decide to change the placeholder for the local Neumann matrix */
1611   {
1612     Mat locpJ;
1613 
1614     ierr = MatISGetLocalMat(pJ,&locpJ);CHKERRQ(ierr);
1615     ierr = MatCopy(locpJ,J,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
1616   }
1617   PetscFunctionReturn(0);
1618 }
1619 
1620 /*@
1621   DMPlexSetSNESLocalFEM - Use DMPlex's internal FEM routines to compute SNES boundary values, residual, and Jacobian.
1622 
1623   Input Parameters:
1624 + dm - The DM object
1625 . boundaryctx - the user context that will be passed to pointwise evaluation of boundary values (see PetscDSAddBoundary())
1626 . residualctx - the user context that will be passed to pointwise evaluation of finite element residual computations (see PetscDSSetResidual())
1627 - jacobianctx - the user context that will be passed to pointwise evaluation of finite element Jacobian construction (see PetscDSSetJacobian())
1628 
1629   Level: developer
1630 @*/
1631 PetscErrorCode DMPlexSetSNESLocalFEM(DM dm, void *boundaryctx, void *residualctx, void *jacobianctx)
1632 {
1633   PetscErrorCode ierr;
1634 
1635   PetscFunctionBegin;
1636   ierr = DMSNESSetBoundaryLocal(dm,DMPlexSNESComputeBoundaryFEM,boundaryctx);CHKERRQ(ierr);
1637   ierr = DMSNESSetFunctionLocal(dm,DMPlexSNESComputeResidualFEM,residualctx);CHKERRQ(ierr);
1638   ierr = DMSNESSetJacobianLocal(dm,DMPlexSNESComputeJacobianFEM,jacobianctx);CHKERRQ(ierr);
1639   ierr = PetscObjectComposeFunction((PetscObject)dm,"MatComputeNeumannOverlap_C",MatComputeNeumannOverlap_Plex);CHKERRQ(ierr);
1640   PetscFunctionReturn(0);
1641 }
1642 
1643 /*@C
1644   DMSNESCheckDiscretization - Check the discretization error of the exact solution
1645 
1646   Input Parameters:
1647 + snes - the SNES object
1648 . dm   - the DM
1649 . t    - the time
1650 . u    - a DM vector
1651 - tol  - A tolerance for the check, or -1 to print the results instead
1652 
1653   Output Parameters:
1654 . error - An array which holds the discretization error in each field, or NULL
1655 
1656   Note: The user must call PetscDSSetExactSolution() beforehand
1657 
1658   Level: developer
1659 
1660 .seealso: DNSNESCheckFromOptions(), DMSNESCheckResidual(), DMSNESCheckJacobian(), PetscDSSetExactSolution()
1661 @*/
1662 PetscErrorCode DMSNESCheckDiscretization(SNES snes, DM dm, PetscReal t, Vec u, PetscReal tol, PetscReal error[])
1663 {
1664   PetscErrorCode (**exacts)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
1665   void            **ectxs;
1666   PetscReal        *err;
1667   MPI_Comm          comm;
1668   PetscInt          Nf, f;
1669   PetscErrorCode    ierr;
1670 
1671   PetscFunctionBegin;
1672   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1673   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1674   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
1675   if (error) PetscValidRealPointer(error, 6);
1676 
1677   ierr = DMComputeExactSolution(dm, t, u, NULL);CHKERRQ(ierr);
1678   ierr = VecViewFromOptions(u, NULL, "-vec_view");CHKERRQ(ierr);
1679 
1680   ierr = PetscObjectGetComm((PetscObject) snes, &comm);CHKERRQ(ierr);
1681   ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr);
1682   ierr = PetscCalloc3(Nf, &exacts, Nf, &ectxs, PetscMax(1, Nf), &err);CHKERRQ(ierr);
1683   {
1684     PetscInt Nds, s;
1685 
1686     ierr = DMGetNumDS(dm, &Nds);CHKERRQ(ierr);
1687     for (s = 0; s < Nds; ++s) {
1688       PetscDS         ds;
1689       DMLabel         label;
1690       IS              fieldIS;
1691       const PetscInt *fields;
1692       PetscInt        dsNf, f;
1693 
1694       ierr = DMGetRegionNumDS(dm, s, &label, &fieldIS, &ds);CHKERRQ(ierr);
1695       ierr = PetscDSGetNumFields(ds, &dsNf);CHKERRQ(ierr);
1696       ierr = ISGetIndices(fieldIS, &fields);CHKERRQ(ierr);
1697       for (f = 0; f < dsNf; ++f) {
1698         const PetscInt field = fields[f];
1699         ierr = PetscDSGetExactSolution(ds, field, &exacts[field], &ectxs[field]);CHKERRQ(ierr);
1700       }
1701       ierr = ISRestoreIndices(fieldIS, &fields);CHKERRQ(ierr);
1702     }
1703   }
1704   if (Nf > 1) {
1705     ierr = DMComputeL2FieldDiff(dm, t, exacts, ectxs, u, err);CHKERRQ(ierr);
1706     if (tol >= 0.0) {
1707       for (f = 0; f < Nf; ++f) {
1708         if (err[f] > tol) SETERRQ3(comm, PETSC_ERR_ARG_WRONG, "L_2 Error %g for field %D exceeds tolerance %g", (double) err[f], f, (double) tol);
1709       }
1710     } else if (error) {
1711       for (f = 0; f < Nf; ++f) error[f] = err[f];
1712     } else {
1713       ierr = PetscPrintf(comm, "L_2 Error: [");CHKERRQ(ierr);
1714       for (f = 0; f < Nf; ++f) {
1715         if (f) {ierr = PetscPrintf(comm, ", ");CHKERRQ(ierr);}
1716         ierr = PetscPrintf(comm, "%g", (double)err[f]);CHKERRQ(ierr);
1717       }
1718       ierr = PetscPrintf(comm, "]\n");CHKERRQ(ierr);
1719     }
1720   } else {
1721     ierr = DMComputeL2Diff(dm, t, exacts, ectxs, u, &err[0]);CHKERRQ(ierr);
1722     if (tol >= 0.0) {
1723       if (err[0] > tol) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "L_2 Error %g exceeds tolerance %g", (double) err[0], (double) tol);
1724     } else if (error) {
1725       error[0] = err[0];
1726     } else {
1727       ierr = PetscPrintf(comm, "L_2 Error: %g\n", (double) err[0]);CHKERRQ(ierr);
1728     }
1729   }
1730   ierr = PetscFree3(exacts, ectxs, err);CHKERRQ(ierr);
1731   PetscFunctionReturn(0);
1732 }
1733 
1734 /*@C
1735   DMSNESCheckResidual - Check the residual of the exact solution
1736 
1737   Input Parameters:
1738 + snes - the SNES object
1739 . dm   - the DM
1740 . u    - a DM vector
1741 - tol  - A tolerance for the check, or -1 to print the results instead
1742 
1743   Output Parameters:
1744 . residual - The residual norm of the exact solution, or NULL
1745 
1746   Level: developer
1747 
1748 .seealso: DNSNESCheckFromOptions(), DMSNESCheckDiscretization(), DMSNESCheckJacobian()
1749 @*/
1750 PetscErrorCode DMSNESCheckResidual(SNES snes, DM dm, Vec u, PetscReal tol, PetscReal *residual)
1751 {
1752   MPI_Comm       comm;
1753   Vec            r;
1754   PetscReal      res;
1755   PetscErrorCode ierr;
1756 
1757   PetscFunctionBegin;
1758   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1759   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1760   PetscValidHeaderSpecific(u, VEC_CLASSID, 3);
1761   if (residual) PetscValidRealPointer(residual, 5);
1762   ierr = PetscObjectGetComm((PetscObject) snes, &comm);CHKERRQ(ierr);
1763   ierr = DMComputeExactSolution(dm, 0.0, u, NULL);CHKERRQ(ierr);
1764   ierr = VecDuplicate(u, &r);CHKERRQ(ierr);
1765   ierr = SNESComputeFunction(snes, u, r);CHKERRQ(ierr);
1766   ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr);
1767   if (tol >= 0.0) {
1768     if (res > tol) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "L_2 Residual %g exceeds tolerance %g", (double) res, (double) tol);
1769   } else if (residual) {
1770     *residual = res;
1771   } else {
1772     ierr = PetscPrintf(comm, "L_2 Residual: %g\n", (double)res);CHKERRQ(ierr);
1773     ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr);
1774     ierr = PetscObjectSetName((PetscObject) r, "Initial Residual");CHKERRQ(ierr);
1775     ierr = PetscObjectSetOptionsPrefix((PetscObject)r,"res_");CHKERRQ(ierr);
1776     ierr = VecViewFromOptions(r, NULL, "-vec_view");CHKERRQ(ierr);
1777   }
1778   ierr = VecDestroy(&r);CHKERRQ(ierr);
1779   PetscFunctionReturn(0);
1780 }
1781 
1782 /*@C
1783   DMSNESCheckJacobian - Check the Jacobian of the exact solution against the residual using the Taylor Test
1784 
1785   Input Parameters:
1786 + snes - the SNES object
1787 . dm   - the DM
1788 . u    - a DM vector
1789 - tol  - A tolerance for the check, or -1 to print the results instead
1790 
1791   Output Parameters:
1792 + isLinear - Flag indicaing that the function looks linear, or NULL
1793 - convRate - The rate of convergence of the linear model, or NULL
1794 
1795   Level: developer
1796 
1797 .seealso: DNSNESCheckFromOptions(), DMSNESCheckDiscretization(), DMSNESCheckResidual()
1798 @*/
1799 PetscErrorCode DMSNESCheckJacobian(SNES snes, DM dm, Vec u, PetscReal tol, PetscBool *isLinear, PetscReal *convRate)
1800 {
1801   MPI_Comm       comm;
1802   PetscDS        ds;
1803   Mat            J, M;
1804   MatNullSpace   nullspace;
1805   PetscReal      slope, intercept;
1806   PetscBool      hasJac, hasPrec, isLin = PETSC_FALSE;
1807   PetscErrorCode ierr;
1808 
1809   PetscFunctionBegin;
1810   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1811   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1812   PetscValidHeaderSpecific(u, VEC_CLASSID, 3);
1813   if (isLinear) PetscValidBoolPointer(isLinear, 5);
1814   if (convRate) PetscValidRealPointer(convRate, 6);
1815   ierr = PetscObjectGetComm((PetscObject) snes, &comm);CHKERRQ(ierr);
1816   ierr = DMComputeExactSolution(dm, 0.0, u, NULL);CHKERRQ(ierr);
1817   /* Create and view matrices */
1818   ierr = DMCreateMatrix(dm, &J);CHKERRQ(ierr);
1819   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
1820   ierr = PetscDSHasJacobian(ds, &hasJac);CHKERRQ(ierr);
1821   ierr = PetscDSHasJacobianPreconditioner(ds, &hasPrec);CHKERRQ(ierr);
1822   if (hasJac && hasPrec) {
1823     ierr = DMCreateMatrix(dm, &M);CHKERRQ(ierr);
1824     ierr = SNESComputeJacobian(snes, u, J, M);CHKERRQ(ierr);
1825     ierr = PetscObjectSetName((PetscObject) M, "Preconditioning Matrix");CHKERRQ(ierr);
1826     ierr = PetscObjectSetOptionsPrefix((PetscObject) M, "jacpre_");CHKERRQ(ierr);
1827     ierr = MatViewFromOptions(M, NULL, "-mat_view");CHKERRQ(ierr);
1828     ierr = MatDestroy(&M);CHKERRQ(ierr);
1829   } else {
1830     ierr = SNESComputeJacobian(snes, u, J, J);CHKERRQ(ierr);
1831   }
1832   ierr = PetscObjectSetName((PetscObject) J, "Jacobian");CHKERRQ(ierr);
1833   ierr = PetscObjectSetOptionsPrefix((PetscObject) J, "jac_");CHKERRQ(ierr);
1834   ierr = MatViewFromOptions(J, NULL, "-mat_view");CHKERRQ(ierr);
1835   /* Check nullspace */
1836   ierr = MatGetNullSpace(J, &nullspace);CHKERRQ(ierr);
1837   if (nullspace) {
1838     PetscBool isNull;
1839     ierr = MatNullSpaceTest(nullspace, J, &isNull);CHKERRQ(ierr);
1840     if (!isNull) SETERRQ(comm, PETSC_ERR_PLIB, "The null space calculated for the system operator is invalid.");
1841   }
1842   /* Taylor test */
1843   {
1844     PetscRandom rand;
1845     Vec         du, uhat, r, rhat, df;
1846     PetscReal   h;
1847     PetscReal  *es, *hs, *errors;
1848     PetscReal   hMax = 1.0, hMin = 1e-6, hMult = 0.1;
1849     PetscInt    Nv, v;
1850 
1851     /* Choose a perturbation direction */
1852     ierr = PetscRandomCreate(comm, &rand);CHKERRQ(ierr);
1853     ierr = VecDuplicate(u, &du);CHKERRQ(ierr);
1854     ierr = VecSetRandom(du, rand);CHKERRQ(ierr);
1855     ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr);
1856     ierr = VecDuplicate(u, &df);CHKERRQ(ierr);
1857     ierr = MatMult(J, du, df);CHKERRQ(ierr);
1858     /* Evaluate residual at u, F(u), save in vector r */
1859     ierr = VecDuplicate(u, &r);CHKERRQ(ierr);
1860     ierr = SNESComputeFunction(snes, u, r);CHKERRQ(ierr);
1861     /* Look at the convergence of our Taylor approximation as we approach u */
1862     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv);
1863     ierr = PetscCalloc3(Nv, &es, Nv, &hs, Nv, &errors);CHKERRQ(ierr);
1864     ierr = VecDuplicate(u, &uhat);CHKERRQ(ierr);
1865     ierr = VecDuplicate(u, &rhat);CHKERRQ(ierr);
1866     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv) {
1867       ierr = VecWAXPY(uhat, h, du, u);CHKERRQ(ierr);
1868       /* F(\hat u) \approx F(u) + J(u) (uhat - u) = F(u) + h * J(u) du */
1869       ierr = SNESComputeFunction(snes, uhat, rhat);CHKERRQ(ierr);
1870       ierr = VecAXPBYPCZ(rhat, -1.0, -h, 1.0, r, df);CHKERRQ(ierr);
1871       ierr = VecNorm(rhat, NORM_2, &errors[Nv]);CHKERRQ(ierr);
1872 
1873       es[Nv] = PetscLog10Real(errors[Nv]);
1874       hs[Nv] = PetscLog10Real(h);
1875     }
1876     ierr = VecDestroy(&uhat);CHKERRQ(ierr);
1877     ierr = VecDestroy(&rhat);CHKERRQ(ierr);
1878     ierr = VecDestroy(&df);CHKERRQ(ierr);
1879     ierr = VecDestroy(&r);CHKERRQ(ierr);
1880     ierr = VecDestroy(&du);CHKERRQ(ierr);
1881     for (v = 0; v < Nv; ++v) {
1882       if ((tol >= 0) && (errors[v] > tol)) break;
1883       else if (errors[v] > PETSC_SMALL)    break;
1884     }
1885     if (v == Nv) isLin = PETSC_TRUE;
1886     ierr = PetscLinearRegression(Nv, hs, es, &slope, &intercept);CHKERRQ(ierr);
1887     ierr = PetscFree3(es, hs, errors);CHKERRQ(ierr);
1888     /* Slope should be about 2 */
1889     if (tol >= 0) {
1890       if (!isLin && PetscAbsReal(2 - slope) > tol) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Taylor approximation convergence rate should be 2, not %0.2f", (double) slope);
1891     } else if (isLinear || convRate) {
1892       if (isLinear) *isLinear = isLin;
1893       if (convRate) *convRate = slope;
1894     } else {
1895       if (!isLin) {ierr = PetscPrintf(comm, "Taylor approximation converging at order %3.2f\n", (double) slope);CHKERRQ(ierr);}
1896       else        {ierr = PetscPrintf(comm, "Function appears to be linear\n");CHKERRQ(ierr);}
1897     }
1898   }
1899   ierr = MatDestroy(&J);CHKERRQ(ierr);
1900   PetscFunctionReturn(0);
1901 }
1902 
1903 PetscErrorCode DMSNESCheck_Internal(SNES snes, DM dm, Vec u)
1904 {
1905   PetscErrorCode ierr;
1906 
1907   PetscFunctionBegin;
1908   ierr = DMSNESCheckDiscretization(snes, dm, 0.0, u, -1.0, NULL);CHKERRQ(ierr);
1909   ierr = DMSNESCheckResidual(snes, dm, u, -1.0, NULL);CHKERRQ(ierr);
1910   ierr = DMSNESCheckJacobian(snes, dm, u, -1.0, NULL, NULL);CHKERRQ(ierr);
1911   PetscFunctionReturn(0);
1912 }
1913 
1914 /*@C
1915   DMSNESCheckFromOptions - Check the residual and Jacobian functions using the exact solution by outputting some diagnostic information
1916 
1917   Input Parameters:
1918 + snes - the SNES object
1919 - u    - representative SNES vector
1920 
1921   Note: The user must call PetscDSSetExactSolution() beforehand
1922 
1923   Level: developer
1924 @*/
1925 PetscErrorCode DMSNESCheckFromOptions(SNES snes, Vec u)
1926 {
1927   DM             dm;
1928   Vec            sol;
1929   PetscBool      check;
1930   PetscErrorCode ierr;
1931 
1932   PetscFunctionBegin;
1933   ierr = PetscOptionsHasName(((PetscObject)snes)->options,((PetscObject)snes)->prefix, "-dmsnes_check", &check);CHKERRQ(ierr);
1934   if (!check) PetscFunctionReturn(0);
1935   ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr);
1936   ierr = VecDuplicate(u, &sol);CHKERRQ(ierr);
1937   ierr = SNESSetSolution(snes, sol);CHKERRQ(ierr);
1938   ierr = DMSNESCheck_Internal(snes, dm, sol);CHKERRQ(ierr);
1939   ierr = VecDestroy(&sol);CHKERRQ(ierr);
1940   PetscFunctionReturn(0);
1941 }
1942