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