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