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