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