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