xref: /petsc/src/snes/utils/dmplexsnes.c (revision fc8a9adeb7fcdc98711d755fa2dc544ddccf0f3e)
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 <petscblaslapack.h>
5 #include <petsc/private/petscimpl.h>
6 #include <petsc/private/petscfeimpl.h>
7 
8 /************************** Interpolation *******************************/
9 
10 static PetscErrorCode DMSNESConvertPlex(DM dm, DM *plex, PetscBool copy)
11 {
12   PetscBool      isPlex;
13   PetscErrorCode ierr;
14 
15   PetscFunctionBegin;
16   ierr = PetscObjectTypeCompare((PetscObject) dm, DMPLEX, &isPlex);CHKERRQ(ierr);
17   if (isPlex) {
18     *plex = dm;
19     ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);
20   } else {
21     ierr = PetscObjectQuery((PetscObject) dm, "dm_plex", (PetscObject *) plex);CHKERRQ(ierr);
22     if (!*plex) {
23       ierr = DMConvert(dm,DMPLEX,plex);CHKERRQ(ierr);
24       ierr = PetscObjectCompose((PetscObject) dm, "dm_plex", (PetscObject) *plex);CHKERRQ(ierr);
25       if (copy) {
26         PetscInt    i;
27         PetscObject obj;
28         const char *comps[3] = {"A","dmAux","dmCh"};
29 
30         ierr = DMCopyDMSNES(dm, *plex);CHKERRQ(ierr);
31         for (i = 0; i < 3; i++) {
32           ierr = PetscObjectQuery((PetscObject) dm, comps[i], &obj);CHKERRQ(ierr);
33           ierr = PetscObjectCompose((PetscObject) *plex, comps[i], obj);CHKERRQ(ierr);
34         }
35       }
36     } else {
37       ierr = PetscObjectReference((PetscObject) *plex);CHKERRQ(ierr);
38     }
39   }
40   PetscFunctionReturn(0);
41 }
42 
43 /*@C
44   DMInterpolationCreate - Creates a DMInterpolationInfo context
45 
46   Collective on comm
47 
48   Input Parameter:
49 . comm - the communicator
50 
51   Output Parameter:
52 . ctx - the context
53 
54   Level: beginner
55 
56 .seealso: DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationDestroy()
57 @*/
58 PetscErrorCode DMInterpolationCreate(MPI_Comm comm, DMInterpolationInfo *ctx)
59 {
60   PetscErrorCode ierr;
61 
62   PetscFunctionBegin;
63   PetscValidPointer(ctx, 2);
64   ierr = PetscNew(ctx);CHKERRQ(ierr);
65 
66   (*ctx)->comm   = comm;
67   (*ctx)->dim    = -1;
68   (*ctx)->nInput = 0;
69   (*ctx)->points = NULL;
70   (*ctx)->cells  = NULL;
71   (*ctx)->n      = -1;
72   (*ctx)->coords = NULL;
73   PetscFunctionReturn(0);
74 }
75 
76 /*@C
77   DMInterpolationSetDim - Sets the spatial dimension for the interpolation context
78 
79   Not collective
80 
81   Input Parameters:
82 + ctx - the context
83 - dim - the spatial dimension
84 
85   Level: intermediate
86 
87 .seealso: DMInterpolationGetDim(), DMInterpolationEvaluate(), DMInterpolationAddPoints()
88 @*/
89 PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo ctx, PetscInt dim)
90 {
91   PetscFunctionBegin;
92   if ((dim < 1) || (dim > 3)) SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension for points: %d", dim);
93   ctx->dim = dim;
94   PetscFunctionReturn(0);
95 }
96 
97 /*@C
98   DMInterpolationGetDim - Gets the spatial dimension for the interpolation context
99 
100   Not collective
101 
102   Input Parameter:
103 . ctx - the context
104 
105   Output Parameter:
106 . dim - the spatial dimension
107 
108   Level: intermediate
109 
110 .seealso: DMInterpolationSetDim(), DMInterpolationEvaluate(), DMInterpolationAddPoints()
111 @*/
112 PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo ctx, PetscInt *dim)
113 {
114   PetscFunctionBegin;
115   PetscValidIntPointer(dim, 2);
116   *dim = ctx->dim;
117   PetscFunctionReturn(0);
118 }
119 
120 /*@C
121   DMInterpolationSetDof - Sets the number of fields interpolated at a point for the interpolation context
122 
123   Not collective
124 
125   Input Parameters:
126 + ctx - the context
127 - dof - the number of fields
128 
129   Level: intermediate
130 
131 .seealso: DMInterpolationGetDof(), DMInterpolationEvaluate(), DMInterpolationAddPoints()
132 @*/
133 PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo ctx, PetscInt dof)
134 {
135   PetscFunctionBegin;
136   if (dof < 1) SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid number of components: %d", dof);
137   ctx->dof = dof;
138   PetscFunctionReturn(0);
139 }
140 
141 /*@C
142   DMInterpolationGetDof - Gets the number of fields interpolated at a point for the interpolation context
143 
144   Not collective
145 
146   Input Parameter:
147 . ctx - the context
148 
149   Output Parameter:
150 . dof - the number of fields
151 
152   Level: intermediate
153 
154 .seealso: DMInterpolationSetDof(), DMInterpolationEvaluate(), DMInterpolationAddPoints()
155 @*/
156 PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo ctx, PetscInt *dof)
157 {
158   PetscFunctionBegin;
159   PetscValidIntPointer(dof, 2);
160   *dof = ctx->dof;
161   PetscFunctionReturn(0);
162 }
163 
164 /*@C
165   DMInterpolationAddPoints - Add points at which we will interpolate the fields
166 
167   Not collective
168 
169   Input Parameters:
170 + ctx    - the context
171 . n      - the number of points
172 - points - the coordinates for each point, an array of size n * dim
173 
174   Note: The coordinate information is copied.
175 
176   Level: intermediate
177 
178 .seealso: DMInterpolationSetDim(), DMInterpolationEvaluate(), DMInterpolationCreate()
179 @*/
180 PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo ctx, PetscInt n, PetscReal points[])
181 {
182   PetscErrorCode ierr;
183 
184   PetscFunctionBegin;
185   if (ctx->dim < 0) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set");
186   if (ctx->points)  SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot add points multiple times yet");
187   ctx->nInput = n;
188 
189   ierr = PetscMalloc1(n*ctx->dim, &ctx->points);CHKERRQ(ierr);
190   ierr = PetscMemcpy(ctx->points, points, n*ctx->dim * sizeof(PetscReal));CHKERRQ(ierr);
191   PetscFunctionReturn(0);
192 }
193 
194 /*@C
195   DMInterpolationSetUp - Computea spatial indices that add in point location during interpolation
196 
197   Collective on ctx
198 
199   Input Parameters:
200 + ctx - the context
201 . dm  - the DM for the function space used for interpolation
202 - redundantPoints - If PETSC_TRUE, all processes are passing in the same array of points. Otherwise, points need to be communicated among processes.
203 
204   Level: intermediate
205 
206 .seealso: DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
207 @*/
208 PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo ctx, DM dm, PetscBool redundantPoints)
209 {
210   MPI_Comm          comm = ctx->comm;
211   PetscScalar       *a;
212   PetscInt          p, q, i;
213   PetscMPIInt       rank, size;
214   PetscErrorCode    ierr;
215   Vec               pointVec;
216   PetscSF           cellSF;
217   PetscLayout       layout;
218   PetscReal         *globalPoints;
219   PetscScalar       *globalPointsScalar;
220   const PetscInt    *ranges;
221   PetscMPIInt       *counts, *displs;
222   const PetscSFNode *foundCells;
223   const PetscInt    *foundPoints;
224   PetscMPIInt       *foundProcs, *globalProcs;
225   PetscInt          n, N, numFound;
226 
227   PetscFunctionBegin;
228   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
229   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
230   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
231   if (ctx->dim < 0) SETERRQ(comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set");
232   /* Locate points */
233   n = ctx->nInput;
234   if (!redundantPoints) {
235     ierr = PetscLayoutCreate(comm, &layout);CHKERRQ(ierr);
236     ierr = PetscLayoutSetBlockSize(layout, 1);CHKERRQ(ierr);
237     ierr = PetscLayoutSetLocalSize(layout, n);CHKERRQ(ierr);
238     ierr = PetscLayoutSetUp(layout);CHKERRQ(ierr);
239     ierr = PetscLayoutGetSize(layout, &N);CHKERRQ(ierr);
240     /* Communicate all points to all processes */
241     ierr = PetscMalloc3(N*ctx->dim,&globalPoints,size,&counts,size,&displs);CHKERRQ(ierr);
242     ierr = PetscLayoutGetRanges(layout, &ranges);CHKERRQ(ierr);
243     for (p = 0; p < size; ++p) {
244       counts[p] = (ranges[p+1] - ranges[p])*ctx->dim;
245       displs[p] = ranges[p]*ctx->dim;
246     }
247     ierr = MPI_Allgatherv(ctx->points, n*ctx->dim, MPIU_REAL, globalPoints, counts, displs, MPIU_REAL, comm);CHKERRQ(ierr);
248   } else {
249     N = n;
250     globalPoints = ctx->points;
251     counts = displs = NULL;
252     layout = NULL;
253   }
254 #if 0
255   ierr = PetscMalloc3(N,&foundCells,N,&foundProcs,N,&globalProcs);CHKERRQ(ierr);
256   /* foundCells[p] = m->locatePoint(&globalPoints[p*ctx->dim]); */
257 #else
258 #if defined(PETSC_USE_COMPLEX)
259   ierr = PetscMalloc1(N*ctx->dim,&globalPointsScalar);CHKERRQ(ierr);
260   for (i=0; i<N*ctx->dim; i++) globalPointsScalar[i] = globalPoints[i];
261 #else
262   globalPointsScalar = globalPoints;
263 #endif
264   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF, ctx->dim, N*ctx->dim, globalPointsScalar, &pointVec);CHKERRQ(ierr);
265   ierr = PetscMalloc2(N,&foundProcs,N,&globalProcs);CHKERRQ(ierr);
266   for (p = 0; p < N; ++p) {foundProcs[p] = size;}
267   cellSF = NULL;
268   ierr = DMLocatePoints(dm, pointVec, DM_POINTLOCATION_REMOVE, &cellSF);CHKERRQ(ierr);
269   ierr = PetscSFGetGraph(cellSF,NULL,&numFound,&foundPoints,&foundCells);CHKERRQ(ierr);
270 #endif
271   for (p = 0; p < numFound; ++p) {
272     if (foundCells[p].index >= 0) foundProcs[foundPoints ? foundPoints[p] : p] = rank;
273   }
274   /* Let the lowest rank process own each point */
275   ierr   = MPIU_Allreduce(foundProcs, globalProcs, N, MPI_INT, MPI_MIN, comm);CHKERRQ(ierr);
276   ctx->n = 0;
277   for (p = 0; p < N; ++p) {
278     if (globalProcs[p] == size) 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));
279     else if (globalProcs[p] == rank) ctx->n++;
280   }
281   /* Create coordinates vector and array of owned cells */
282   ierr = PetscMalloc1(ctx->n, &ctx->cells);CHKERRQ(ierr);
283   ierr = VecCreate(comm, &ctx->coords);CHKERRQ(ierr);
284   ierr = VecSetSizes(ctx->coords, ctx->n*ctx->dim, PETSC_DECIDE);CHKERRQ(ierr);
285   ierr = VecSetBlockSize(ctx->coords, ctx->dim);CHKERRQ(ierr);
286   ierr = VecSetType(ctx->coords,VECSTANDARD);CHKERRQ(ierr);
287   ierr = VecGetArray(ctx->coords, &a);CHKERRQ(ierr);
288   for (p = 0, q = 0, i = 0; p < N; ++p) {
289     if (globalProcs[p] == rank) {
290       PetscInt d;
291 
292       for (d = 0; d < ctx->dim; ++d, ++i) a[i] = globalPoints[p*ctx->dim+d];
293       ctx->cells[q] = foundCells[q].index;
294       ++q;
295     }
296   }
297   ierr = VecRestoreArray(ctx->coords, &a);CHKERRQ(ierr);
298 #if 0
299   ierr = PetscFree3(foundCells,foundProcs,globalProcs);CHKERRQ(ierr);
300 #else
301   ierr = PetscFree2(foundProcs,globalProcs);CHKERRQ(ierr);
302   ierr = PetscSFDestroy(&cellSF);CHKERRQ(ierr);
303   ierr = VecDestroy(&pointVec);CHKERRQ(ierr);
304 #endif
305   if ((void*)globalPointsScalar != (void*)globalPoints) {ierr = PetscFree(globalPointsScalar);CHKERRQ(ierr);}
306   if (!redundantPoints) {ierr = PetscFree3(globalPoints,counts,displs);CHKERRQ(ierr);}
307   ierr = PetscLayoutDestroy(&layout);CHKERRQ(ierr);
308   PetscFunctionReturn(0);
309 }
310 
311 /*@C
312   DMInterpolationGetCoordinates - Gets a Vec with the coordinates of each interpolation point
313 
314   Collective on ctx
315 
316   Input Parameter:
317 . ctx - the context
318 
319   Output Parameter:
320 . coordinates  - the coordinates of interpolation points
321 
322   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.
323 
324   Level: intermediate
325 
326 .seealso: DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
327 @*/
328 PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo ctx, Vec *coordinates)
329 {
330   PetscFunctionBegin;
331   PetscValidPointer(coordinates, 2);
332   if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
333   *coordinates = ctx->coords;
334   PetscFunctionReturn(0);
335 }
336 
337 /*@C
338   DMInterpolationGetVector - Gets a Vec which can hold all the interpolated field values
339 
340   Collective on ctx
341 
342   Input Parameter:
343 . ctx - the context
344 
345   Output Parameter:
346 . v  - a vector capable of holding the interpolated field values
347 
348   Note: This vector should be returned using DMInterpolationRestoreVector().
349 
350   Level: intermediate
351 
352 .seealso: DMInterpolationRestoreVector(), DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
353 @*/
354 PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo ctx, Vec *v)
355 {
356   PetscErrorCode ierr;
357 
358   PetscFunctionBegin;
359   PetscValidPointer(v, 2);
360   if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
361   ierr = VecCreate(ctx->comm, v);CHKERRQ(ierr);
362   ierr = VecSetSizes(*v, ctx->n*ctx->dof, PETSC_DECIDE);CHKERRQ(ierr);
363   ierr = VecSetBlockSize(*v, ctx->dof);CHKERRQ(ierr);
364   ierr = VecSetType(*v,VECSTANDARD);CHKERRQ(ierr);
365   PetscFunctionReturn(0);
366 }
367 
368 /*@C
369   DMInterpolationRestoreVector - Returns a Vec which can hold all the interpolated field values
370 
371   Collective on ctx
372 
373   Input Parameters:
374 + ctx - the context
375 - v  - a vector capable of holding the interpolated field values
376 
377   Level: intermediate
378 
379 .seealso: DMInterpolationGetVector(), DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
380 @*/
381 PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo ctx, Vec *v)
382 {
383   PetscErrorCode ierr;
384 
385   PetscFunctionBegin;
386   PetscValidPointer(v, 2);
387   if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
388   ierr = VecDestroy(v);CHKERRQ(ierr);
389   PetscFunctionReturn(0);
390 }
391 
392 PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Triangle_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
393 {
394   PetscReal      *v0, *J, *invJ, detJ;
395   const PetscScalar *coords;
396   PetscScalar    *a;
397   PetscInt       p;
398   PetscErrorCode ierr;
399 
400   PetscFunctionBegin;
401   ierr = PetscMalloc3(ctx->dim,&v0,ctx->dim*ctx->dim,&J,ctx->dim*ctx->dim,&invJ);CHKERRQ(ierr);
402   ierr = VecGetArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
403   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
404   for (p = 0; p < ctx->n; ++p) {
405     PetscInt     c = ctx->cells[p];
406     PetscScalar *x = NULL;
407     PetscReal    xi[4];
408     PetscInt     d, f, comp;
409 
410     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
411     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", (double)detJ, c);
412     ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr);
413     for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] = x[0*ctx->dof+comp];
414 
415     for (d = 0; d < ctx->dim; ++d) {
416       xi[d] = 0.0;
417       for (f = 0; f < ctx->dim; ++f) xi[d] += invJ[d*ctx->dim+f]*0.5*PetscRealPart(coords[p*ctx->dim+f] - v0[f]);
418       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];
419     }
420     ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr);
421   }
422   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
423   ierr = VecRestoreArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
424   ierr = PetscFree3(v0, J, invJ);CHKERRQ(ierr);
425   PetscFunctionReturn(0);
426 }
427 
428 PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Tetrahedron_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
429 {
430   PetscReal      *v0, *J, *invJ, detJ;
431   const PetscScalar *coords;
432   PetscScalar    *a;
433   PetscInt       p;
434   PetscErrorCode ierr;
435 
436   PetscFunctionBegin;
437   ierr = PetscMalloc3(ctx->dim,&v0,ctx->dim*ctx->dim,&J,ctx->dim*ctx->dim,&invJ);CHKERRQ(ierr);
438   ierr = VecGetArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
439   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
440   for (p = 0; p < ctx->n; ++p) {
441     PetscInt       c = ctx->cells[p];
442     const PetscInt order[3] = {2, 1, 3};
443     PetscScalar   *x = NULL;
444     PetscReal      xi[4];
445     PetscInt       d, f, comp;
446 
447     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
448     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", (double)detJ, c);
449     ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr);
450     for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] = x[0*ctx->dof+comp];
451 
452     for (d = 0; d < ctx->dim; ++d) {
453       xi[d] = 0.0;
454       for (f = 0; f < ctx->dim; ++f) xi[d] += invJ[d*ctx->dim+f]*0.5*PetscRealPart(coords[p*ctx->dim+f] - v0[f]);
455       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];
456     }
457     ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr);
458   }
459   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
460   ierr = VecRestoreArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
461   ierr = PetscFree3(v0, J, invJ);CHKERRQ(ierr);
462   PetscFunctionReturn(0);
463 }
464 
465 PETSC_STATIC_INLINE PetscErrorCode QuadMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx)
466 {
467   const PetscScalar *vertices = (const PetscScalar*) ctx;
468   const PetscScalar x0        = vertices[0];
469   const PetscScalar y0        = vertices[1];
470   const PetscScalar x1        = vertices[2];
471   const PetscScalar y1        = vertices[3];
472   const PetscScalar x2        = vertices[4];
473   const PetscScalar y2        = vertices[5];
474   const PetscScalar x3        = vertices[6];
475   const PetscScalar y3        = vertices[7];
476   const PetscScalar f_1       = x1 - x0;
477   const PetscScalar g_1       = y1 - y0;
478   const PetscScalar f_3       = x3 - x0;
479   const PetscScalar g_3       = y3 - y0;
480   const PetscScalar f_01      = x2 - x1 - x3 + x0;
481   const PetscScalar g_01      = y2 - y1 - y3 + y0;
482   const PetscScalar *ref;
483   PetscScalar       *real;
484   PetscErrorCode    ierr;
485 
486   PetscFunctionBegin;
487   ierr = VecGetArrayRead(Xref,  &ref);CHKERRQ(ierr);
488   ierr = VecGetArray(Xreal, &real);CHKERRQ(ierr);
489   {
490     const PetscScalar p0 = ref[0];
491     const PetscScalar p1 = ref[1];
492 
493     real[0] = x0 + f_1 * p0 + f_3 * p1 + f_01 * p0 * p1;
494     real[1] = y0 + g_1 * p0 + g_3 * p1 + g_01 * p0 * p1;
495   }
496   ierr = PetscLogFlops(28);CHKERRQ(ierr);
497   ierr = VecRestoreArrayRead(Xref,  &ref);CHKERRQ(ierr);
498   ierr = VecRestoreArray(Xreal, &real);CHKERRQ(ierr);
499   PetscFunctionReturn(0);
500 }
501 
502 #include <petsc/private/dmimpl.h>
503 PETSC_STATIC_INLINE PetscErrorCode QuadJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx)
504 {
505   const PetscScalar *vertices = (const PetscScalar*) ctx;
506   const PetscScalar x0        = vertices[0];
507   const PetscScalar y0        = vertices[1];
508   const PetscScalar x1        = vertices[2];
509   const PetscScalar y1        = vertices[3];
510   const PetscScalar x2        = vertices[4];
511   const PetscScalar y2        = vertices[5];
512   const PetscScalar x3        = vertices[6];
513   const PetscScalar y3        = vertices[7];
514   const PetscScalar f_01      = x2 - x1 - x3 + x0;
515   const PetscScalar g_01      = y2 - y1 - y3 + y0;
516   const PetscScalar *ref;
517   PetscErrorCode    ierr;
518 
519   PetscFunctionBegin;
520   ierr = VecGetArrayRead(Xref,  &ref);CHKERRQ(ierr);
521   {
522     const PetscScalar x       = ref[0];
523     const PetscScalar y       = ref[1];
524     const PetscInt    rows[2] = {0, 1};
525     PetscScalar       values[4];
526 
527     values[0] = (x1 - x0 + f_01*y) * 0.5; values[1] = (x3 - x0 + f_01*x) * 0.5;
528     values[2] = (y1 - y0 + g_01*y) * 0.5; values[3] = (y3 - y0 + g_01*x) * 0.5;
529     ierr      = MatSetValues(J, 2, rows, 2, rows, values, INSERT_VALUES);CHKERRQ(ierr);
530   }
531   ierr = PetscLogFlops(30);CHKERRQ(ierr);
532   ierr = VecRestoreArrayRead(Xref,  &ref);CHKERRQ(ierr);
533   ierr = MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
534   ierr = MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
535   PetscFunctionReturn(0);
536 }
537 
538 PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Quad_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
539 {
540   DM             dmCoord;
541   PetscFE        fem = NULL;
542   SNES           snes;
543   KSP            ksp;
544   PC             pc;
545   Vec            coordsLocal, r, ref, real;
546   Mat            J;
547   const PetscScalar *coords;
548   PetscScalar    *a;
549   PetscInt       Nf, p;
550   const PetscInt dof = ctx->dof;
551   PetscErrorCode ierr;
552 
553   PetscFunctionBegin;
554   ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr);
555   if (Nf) {ierr = DMGetField(dm, 0, NULL, (PetscObject *) &fem);CHKERRQ(ierr);}
556   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
557   ierr = DMGetCoordinateDM(dm, &dmCoord);CHKERRQ(ierr);
558   ierr = SNESCreate(PETSC_COMM_SELF, &snes);CHKERRQ(ierr);
559   ierr = SNESSetOptionsPrefix(snes, "quad_interp_");CHKERRQ(ierr);
560   ierr = VecCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr);
561   ierr = VecSetSizes(r, 2, 2);CHKERRQ(ierr);
562   ierr = VecSetType(r,dm->vectype);CHKERRQ(ierr);
563   ierr = VecDuplicate(r, &ref);CHKERRQ(ierr);
564   ierr = VecDuplicate(r, &real);CHKERRQ(ierr);
565   ierr = MatCreate(PETSC_COMM_SELF, &J);CHKERRQ(ierr);
566   ierr = MatSetSizes(J, 2, 2, 2, 2);CHKERRQ(ierr);
567   ierr = MatSetType(J, MATSEQDENSE);CHKERRQ(ierr);
568   ierr = MatSetUp(J);CHKERRQ(ierr);
569   ierr = SNESSetFunction(snes, r, QuadMap_Private, NULL);CHKERRQ(ierr);
570   ierr = SNESSetJacobian(snes, J, J, QuadJacobian_Private, NULL);CHKERRQ(ierr);
571   ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr);
572   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
573   ierr = PCSetType(pc, PCLU);CHKERRQ(ierr);
574   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
575 
576   ierr = VecGetArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
577   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
578   for (p = 0; p < ctx->n; ++p) {
579     PetscScalar *x = NULL, *vertices = NULL;
580     PetscScalar *xi;
581     PetscReal    xir[2];
582     PetscInt     c = ctx->cells[p], comp, coordSize, xSize;
583 
584     /* Can make this do all points at once */
585     ierr = DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
586     if (4*2 != coordSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", coordSize, 4*2);
587     ierr   = DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
588     ierr   = SNESSetFunction(snes, NULL, NULL, (void*) vertices);CHKERRQ(ierr);
589     ierr   = SNESSetJacobian(snes, NULL, NULL, NULL, (void*) vertices);CHKERRQ(ierr);
590     ierr   = VecGetArray(real, &xi);CHKERRQ(ierr);
591     xi[0]  = coords[p*ctx->dim+0];
592     xi[1]  = coords[p*ctx->dim+1];
593     ierr   = VecRestoreArray(real, &xi);CHKERRQ(ierr);
594     ierr   = SNESSolve(snes, real, ref);CHKERRQ(ierr);
595     ierr   = VecGetArray(ref, &xi);CHKERRQ(ierr);
596     xir[0] = PetscRealPart(xi[0]);
597     xir[1] = PetscRealPart(xi[1]);
598     if (4*dof != xSize) {
599       PetscReal *B;
600       PetscInt   d;
601 
602       xir[0] = 2.0*xir[0] - 1.0; xir[1] = 2.0*xir[1] - 1.0;
603       ierr = PetscFEGetTabulation(fem, 1, xir, &B, NULL, NULL);CHKERRQ(ierr);
604       for (comp = 0; comp < dof; ++comp) {
605         a[p*dof+comp] = 0.0;
606         for (d = 0; d < xSize/dof; ++d) {
607           a[p*dof+comp] += x[d*dof+comp]*B[d*dof+comp];
608         }
609       }
610       ierr = PetscFERestoreTabulation(fem, 1, xir, &B, NULL, NULL);CHKERRQ(ierr);
611     } else {
612       for (comp = 0; comp < dof; ++comp)
613         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];
614     }
615     ierr = VecRestoreArray(ref, &xi);CHKERRQ(ierr);
616     ierr = DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
617     ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
618   }
619   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
620   ierr = VecRestoreArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
621 
622   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
623   ierr = VecDestroy(&r);CHKERRQ(ierr);
624   ierr = VecDestroy(&ref);CHKERRQ(ierr);
625   ierr = VecDestroy(&real);CHKERRQ(ierr);
626   ierr = MatDestroy(&J);CHKERRQ(ierr);
627   PetscFunctionReturn(0);
628 }
629 
630 PETSC_STATIC_INLINE PetscErrorCode HexMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx)
631 {
632   const PetscScalar *vertices = (const PetscScalar*) ctx;
633   const PetscScalar x0        = vertices[0];
634   const PetscScalar y0        = vertices[1];
635   const PetscScalar z0        = vertices[2];
636   const PetscScalar x1        = vertices[9];
637   const PetscScalar y1        = vertices[10];
638   const PetscScalar z1        = vertices[11];
639   const PetscScalar x2        = vertices[6];
640   const PetscScalar y2        = vertices[7];
641   const PetscScalar z2        = vertices[8];
642   const PetscScalar x3        = vertices[3];
643   const PetscScalar y3        = vertices[4];
644   const PetscScalar z3        = vertices[5];
645   const PetscScalar x4        = vertices[12];
646   const PetscScalar y4        = vertices[13];
647   const PetscScalar z4        = vertices[14];
648   const PetscScalar x5        = vertices[15];
649   const PetscScalar y5        = vertices[16];
650   const PetscScalar z5        = vertices[17];
651   const PetscScalar x6        = vertices[18];
652   const PetscScalar y6        = vertices[19];
653   const PetscScalar z6        = vertices[20];
654   const PetscScalar x7        = vertices[21];
655   const PetscScalar y7        = vertices[22];
656   const PetscScalar z7        = vertices[23];
657   const PetscScalar f_1       = x1 - x0;
658   const PetscScalar g_1       = y1 - y0;
659   const PetscScalar h_1       = z1 - z0;
660   const PetscScalar f_3       = x3 - x0;
661   const PetscScalar g_3       = y3 - y0;
662   const PetscScalar h_3       = z3 - z0;
663   const PetscScalar f_4       = x4 - x0;
664   const PetscScalar g_4       = y4 - y0;
665   const PetscScalar h_4       = z4 - z0;
666   const PetscScalar f_01      = x2 - x1 - x3 + x0;
667   const PetscScalar g_01      = y2 - y1 - y3 + y0;
668   const PetscScalar h_01      = z2 - z1 - z3 + z0;
669   const PetscScalar f_12      = x7 - x3 - x4 + x0;
670   const PetscScalar g_12      = y7 - y3 - y4 + y0;
671   const PetscScalar h_12      = z7 - z3 - z4 + z0;
672   const PetscScalar f_02      = x5 - x1 - x4 + x0;
673   const PetscScalar g_02      = y5 - y1 - y4 + y0;
674   const PetscScalar h_02      = z5 - z1 - z4 + z0;
675   const PetscScalar f_012     = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7;
676   const PetscScalar g_012     = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7;
677   const PetscScalar h_012     = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7;
678   const PetscScalar *ref;
679   PetscScalar       *real;
680   PetscErrorCode    ierr;
681 
682   PetscFunctionBegin;
683   ierr = VecGetArrayRead(Xref,  &ref);CHKERRQ(ierr);
684   ierr = VecGetArray(Xreal, &real);CHKERRQ(ierr);
685   {
686     const PetscScalar p0 = ref[0];
687     const PetscScalar p1 = ref[1];
688     const PetscScalar p2 = ref[2];
689 
690     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;
691     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;
692     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;
693   }
694   ierr = PetscLogFlops(114);CHKERRQ(ierr);
695   ierr = VecRestoreArrayRead(Xref,  &ref);CHKERRQ(ierr);
696   ierr = VecRestoreArray(Xreal, &real);CHKERRQ(ierr);
697   PetscFunctionReturn(0);
698 }
699 
700 PETSC_STATIC_INLINE PetscErrorCode HexJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx)
701 {
702   const PetscScalar *vertices = (const PetscScalar*) ctx;
703   const PetscScalar x0        = vertices[0];
704   const PetscScalar y0        = vertices[1];
705   const PetscScalar z0        = vertices[2];
706   const PetscScalar x1        = vertices[9];
707   const PetscScalar y1        = vertices[10];
708   const PetscScalar z1        = vertices[11];
709   const PetscScalar x2        = vertices[6];
710   const PetscScalar y2        = vertices[7];
711   const PetscScalar z2        = vertices[8];
712   const PetscScalar x3        = vertices[3];
713   const PetscScalar y3        = vertices[4];
714   const PetscScalar z3        = vertices[5];
715   const PetscScalar x4        = vertices[12];
716   const PetscScalar y4        = vertices[13];
717   const PetscScalar z4        = vertices[14];
718   const PetscScalar x5        = vertices[15];
719   const PetscScalar y5        = vertices[16];
720   const PetscScalar z5        = vertices[17];
721   const PetscScalar x6        = vertices[18];
722   const PetscScalar y6        = vertices[19];
723   const PetscScalar z6        = vertices[20];
724   const PetscScalar x7        = vertices[21];
725   const PetscScalar y7        = vertices[22];
726   const PetscScalar z7        = vertices[23];
727   const PetscScalar f_xy      = x2 - x1 - x3 + x0;
728   const PetscScalar g_xy      = y2 - y1 - y3 + y0;
729   const PetscScalar h_xy      = z2 - z1 - z3 + z0;
730   const PetscScalar f_yz      = x7 - x3 - x4 + x0;
731   const PetscScalar g_yz      = y7 - y3 - y4 + y0;
732   const PetscScalar h_yz      = z7 - z3 - z4 + z0;
733   const PetscScalar f_xz      = x5 - x1 - x4 + x0;
734   const PetscScalar g_xz      = y5 - y1 - y4 + y0;
735   const PetscScalar h_xz      = z5 - z1 - z4 + z0;
736   const PetscScalar f_xyz     = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7;
737   const PetscScalar g_xyz     = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7;
738   const PetscScalar h_xyz     = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7;
739   const PetscScalar *ref;
740   PetscErrorCode    ierr;
741 
742   PetscFunctionBegin;
743   ierr = VecGetArrayRead(Xref,  &ref);CHKERRQ(ierr);
744   {
745     const PetscScalar x       = ref[0];
746     const PetscScalar y       = ref[1];
747     const PetscScalar z       = ref[2];
748     const PetscInt    rows[3] = {0, 1, 2};
749     PetscScalar       values[9];
750 
751     values[0] = (x1 - x0 + f_xy*y + f_xz*z + f_xyz*y*z) / 2.0;
752     values[1] = (x3 - x0 + f_xy*x + f_yz*z + f_xyz*x*z) / 2.0;
753     values[2] = (x4 - x0 + f_yz*y + f_xz*x + f_xyz*x*y) / 2.0;
754     values[3] = (y1 - y0 + g_xy*y + g_xz*z + g_xyz*y*z) / 2.0;
755     values[4] = (y3 - y0 + g_xy*x + g_yz*z + g_xyz*x*z) / 2.0;
756     values[5] = (y4 - y0 + g_yz*y + g_xz*x + g_xyz*x*y) / 2.0;
757     values[6] = (z1 - z0 + h_xy*y + h_xz*z + h_xyz*y*z) / 2.0;
758     values[7] = (z3 - z0 + h_xy*x + h_yz*z + h_xyz*x*z) / 2.0;
759     values[8] = (z4 - z0 + h_yz*y + h_xz*x + h_xyz*x*y) / 2.0;
760 
761     ierr = MatSetValues(J, 3, rows, 3, rows, values, INSERT_VALUES);CHKERRQ(ierr);
762   }
763   ierr = PetscLogFlops(152);CHKERRQ(ierr);
764   ierr = VecRestoreArrayRead(Xref,  &ref);CHKERRQ(ierr);
765   ierr = MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
766   ierr = MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
767   PetscFunctionReturn(0);
768 }
769 
770 PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Hex_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
771 {
772   DM             dmCoord;
773   SNES           snes;
774   KSP            ksp;
775   PC             pc;
776   Vec            coordsLocal, r, ref, real;
777   Mat            J;
778   const PetscScalar *coords;
779   PetscScalar    *a;
780   PetscInt       p;
781   PetscErrorCode ierr;
782 
783   PetscFunctionBegin;
784   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
785   ierr = DMGetCoordinateDM(dm, &dmCoord);CHKERRQ(ierr);
786   ierr = SNESCreate(PETSC_COMM_SELF, &snes);CHKERRQ(ierr);
787   ierr = SNESSetOptionsPrefix(snes, "hex_interp_");CHKERRQ(ierr);
788   ierr = VecCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr);
789   ierr = VecSetSizes(r, 3, 3);CHKERRQ(ierr);
790   ierr = VecSetType(r,dm->vectype);CHKERRQ(ierr);
791   ierr = VecDuplicate(r, &ref);CHKERRQ(ierr);
792   ierr = VecDuplicate(r, &real);CHKERRQ(ierr);
793   ierr = MatCreate(PETSC_COMM_SELF, &J);CHKERRQ(ierr);
794   ierr = MatSetSizes(J, 3, 3, 3, 3);CHKERRQ(ierr);
795   ierr = MatSetType(J, MATSEQDENSE);CHKERRQ(ierr);
796   ierr = MatSetUp(J);CHKERRQ(ierr);
797   ierr = SNESSetFunction(snes, r, HexMap_Private, NULL);CHKERRQ(ierr);
798   ierr = SNESSetJacobian(snes, J, J, HexJacobian_Private, NULL);CHKERRQ(ierr);
799   ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr);
800   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
801   ierr = PCSetType(pc, PCLU);CHKERRQ(ierr);
802   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
803 
804   ierr = VecGetArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
805   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
806   for (p = 0; p < ctx->n; ++p) {
807     PetscScalar *x = NULL, *vertices = NULL;
808     PetscScalar *xi;
809     PetscReal    xir[3];
810     PetscInt     c = ctx->cells[p], comp, coordSize, xSize;
811 
812     /* Can make this do all points at once */
813     ierr = DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
814     if (8*3 != coordSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", coordSize, 8*3);
815     ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
816     if (8*ctx->dof != xSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", xSize, 8*ctx->dof);
817     ierr   = SNESSetFunction(snes, NULL, NULL, (void*) vertices);CHKERRQ(ierr);
818     ierr   = SNESSetJacobian(snes, NULL, NULL, NULL, (void*) vertices);CHKERRQ(ierr);
819     ierr   = VecGetArray(real, &xi);CHKERRQ(ierr);
820     xi[0]  = coords[p*ctx->dim+0];
821     xi[1]  = coords[p*ctx->dim+1];
822     xi[2]  = coords[p*ctx->dim+2];
823     ierr   = VecRestoreArray(real, &xi);CHKERRQ(ierr);
824     ierr   = SNESSolve(snes, real, ref);CHKERRQ(ierr);
825     ierr   = VecGetArray(ref, &xi);CHKERRQ(ierr);
826     xir[0] = PetscRealPart(xi[0]);
827     xir[1] = PetscRealPart(xi[1]);
828     xir[2] = PetscRealPart(xi[2]);
829     for (comp = 0; comp < ctx->dof; ++comp) {
830       a[p*ctx->dof+comp] =
831         x[0*ctx->dof+comp]*(1-xir[0])*(1-xir[1])*(1-xir[2]) +
832         x[3*ctx->dof+comp]*    xir[0]*(1-xir[1])*(1-xir[2]) +
833         x[2*ctx->dof+comp]*    xir[0]*    xir[1]*(1-xir[2]) +
834         x[1*ctx->dof+comp]*(1-xir[0])*    xir[1]*(1-xir[2]) +
835         x[4*ctx->dof+comp]*(1-xir[0])*(1-xir[1])*   xir[2] +
836         x[5*ctx->dof+comp]*    xir[0]*(1-xir[1])*   xir[2] +
837         x[6*ctx->dof+comp]*    xir[0]*    xir[1]*   xir[2] +
838         x[7*ctx->dof+comp]*(1-xir[0])*    xir[1]*   xir[2];
839     }
840     ierr = VecRestoreArray(ref, &xi);CHKERRQ(ierr);
841     ierr = DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
842     ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
843   }
844   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
845   ierr = VecRestoreArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
846 
847   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
848   ierr = VecDestroy(&r);CHKERRQ(ierr);
849   ierr = VecDestroy(&ref);CHKERRQ(ierr);
850   ierr = VecDestroy(&real);CHKERRQ(ierr);
851   ierr = MatDestroy(&J);CHKERRQ(ierr);
852   PetscFunctionReturn(0);
853 }
854 
855 /*@C
856   DMInterpolationEvaluate - Using the input from dm and x, calculates interpolated field values at the interpolation points.
857 
858   Input Parameters:
859 + ctx - The DMInterpolationInfo context
860 . dm  - The DM
861 - x   - The local vector containing the field to be interpolated
862 
863   Output Parameters:
864 . v   - The vector containing the interpolated values
865 
866   Note: A suitable v can be obtained using DMInterpolationGetVector().
867 
868   Level: beginner
869 
870 .seealso: DMInterpolationGetVector(), DMInterpolationAddPoints(), DMInterpolationCreate()
871 @*/
872 PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo ctx, DM dm, Vec x, Vec v)
873 {
874   PetscInt       dim, coneSize, n;
875   PetscErrorCode ierr;
876 
877   PetscFunctionBegin;
878   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
879   PetscValidHeaderSpecific(x, VEC_CLASSID, 3);
880   PetscValidHeaderSpecific(v, VEC_CLASSID, 4);
881   ierr = VecGetLocalSize(v, &n);CHKERRQ(ierr);
882   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);
883   if (n) {
884     ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
885     ierr = DMPlexGetConeSize(dm, ctx->cells[0], &coneSize);CHKERRQ(ierr);
886     if (dim == 2) {
887       if (coneSize == 3) {
888         ierr = DMInterpolate_Triangle_Private(ctx, dm, x, v);CHKERRQ(ierr);
889       } else if (coneSize == 4) {
890         ierr = DMInterpolate_Quad_Private(ctx, dm, x, v);CHKERRQ(ierr);
891       } else SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %d for point interpolation", dim);
892     } else if (dim == 3) {
893       if (coneSize == 4) {
894         ierr = DMInterpolate_Tetrahedron_Private(ctx, dm, x, v);CHKERRQ(ierr);
895       } else {
896         ierr = DMInterpolate_Hex_Private(ctx, dm, x, v);CHKERRQ(ierr);
897       }
898     } else SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %d for point interpolation", dim);
899   }
900   PetscFunctionReturn(0);
901 }
902 
903 /*@C
904   DMInterpolationDestroy - Destroys a DMInterpolationInfo context
905 
906   Collective on ctx
907 
908   Input Parameter:
909 . ctx - the context
910 
911   Level: beginner
912 
913 .seealso: DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
914 @*/
915 PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *ctx)
916 {
917   PetscErrorCode ierr;
918 
919   PetscFunctionBegin;
920   PetscValidPointer(ctx, 2);
921   ierr = VecDestroy(&(*ctx)->coords);CHKERRQ(ierr);
922   ierr = PetscFree((*ctx)->points);CHKERRQ(ierr);
923   ierr = PetscFree((*ctx)->cells);CHKERRQ(ierr);
924   ierr = PetscFree(*ctx);CHKERRQ(ierr);
925   *ctx = NULL;
926   PetscFunctionReturn(0);
927 }
928 
929 /*@C
930   SNESMonitorFields - Monitors the residual for each field separately
931 
932   Collective on SNES
933 
934   Input Parameters:
935 + snes   - the SNES context
936 . its    - iteration number
937 . fgnorm - 2-norm of residual
938 - vf  - PetscViewerAndFormat of type ASCII
939 
940   Notes:
941   This routine prints the residual norm at each iteration.
942 
943   Level: intermediate
944 
945 .keywords: SNES, nonlinear, default, monitor, norm
946 .seealso: SNESMonitorSet(), SNESMonitorDefault()
947 @*/
948 PetscErrorCode SNESMonitorFields(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf)
949 {
950   PetscViewer        viewer = vf->viewer;
951   Vec                res;
952   DM                 dm;
953   PetscSection       s;
954   const PetscScalar *r;
955   PetscReal         *lnorms, *norms;
956   PetscInt           numFields, f, pStart, pEnd, p;
957   PetscErrorCode     ierr;
958 
959   PetscFunctionBegin;
960   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4);
961   ierr = SNESGetFunction(snes, &res, 0, 0);CHKERRQ(ierr);
962   ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr);
963   ierr = DMGetSection(dm, &s);CHKERRQ(ierr);
964   ierr = PetscSectionGetNumFields(s, &numFields);CHKERRQ(ierr);
965   ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr);
966   ierr = PetscCalloc2(numFields, &lnorms, numFields, &norms);CHKERRQ(ierr);
967   ierr = VecGetArrayRead(res, &r);CHKERRQ(ierr);
968   for (p = pStart; p < pEnd; ++p) {
969     for (f = 0; f < numFields; ++f) {
970       PetscInt fdof, foff, d;
971 
972       ierr = PetscSectionGetFieldDof(s, p, f, &fdof);CHKERRQ(ierr);
973       ierr = PetscSectionGetFieldOffset(s, p, f, &foff);CHKERRQ(ierr);
974       for (d = 0; d < fdof; ++d) lnorms[f] += PetscRealPart(PetscSqr(r[foff+d]));
975     }
976   }
977   ierr = VecRestoreArrayRead(res, &r);CHKERRQ(ierr);
978   ierr = MPIU_Allreduce(lnorms, norms, numFields, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject) dm));CHKERRQ(ierr);
979   ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr);
980   ierr = PetscViewerASCIIAddTab(viewer, ((PetscObject) snes)->tablevel);CHKERRQ(ierr);
981   ierr = PetscViewerASCIIPrintf(viewer, "%3D SNES Function norm %14.12e [", its, (double) fgnorm);CHKERRQ(ierr);
982   for (f = 0; f < numFields; ++f) {
983     if (f > 0) {ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr);}
984     ierr = PetscViewerASCIIPrintf(viewer, "%14.12e", (double) PetscSqrtReal(norms[f]));CHKERRQ(ierr);
985   }
986   ierr = PetscViewerASCIIPrintf(viewer, "]\n");CHKERRQ(ierr);
987   ierr = PetscViewerASCIISubtractTab(viewer, ((PetscObject) snes)->tablevel);CHKERRQ(ierr);
988   ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
989   ierr = PetscFree2(lnorms, norms);CHKERRQ(ierr);
990   PetscFunctionReturn(0);
991 }
992 
993 /********************* Residual Computation **************************/
994 
995 
996 /*@
997   DMPlexSNESGetGeometryFVM - Return precomputed geometric data
998 
999   Input Parameter:
1000 . dm - The DM
1001 
1002   Output Parameters:
1003 + facegeom - The values precomputed from face geometry
1004 . cellgeom - The values precomputed from cell geometry
1005 - minRadius - The minimum radius over the mesh of an inscribed sphere in a cell
1006 
1007   Level: developer
1008 
1009 .seealso: DMPlexTSSetRHSFunctionLocal()
1010 @*/
1011 PetscErrorCode DMPlexSNESGetGeometryFVM(DM dm, Vec *facegeom, Vec *cellgeom, PetscReal *minRadius)
1012 {
1013   DM             plex;
1014   PetscErrorCode ierr;
1015 
1016   PetscFunctionBegin;
1017   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1018   ierr = DMSNESConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
1019   ierr = DMPlexGetDataFVM(plex, NULL, cellgeom, facegeom, NULL);CHKERRQ(ierr);
1020   if (minRadius) {ierr = DMPlexGetMinRadius(plex, minRadius);CHKERRQ(ierr);}
1021   ierr = DMDestroy(&plex);CHKERRQ(ierr);
1022   PetscFunctionReturn(0);
1023 }
1024 
1025 /*@
1026   DMPlexSNESGetGradientDM - Return gradient data layout
1027 
1028   Input Parameters:
1029 + dm - The DM
1030 - fv - The PetscFV
1031 
1032   Output Parameter:
1033 . dmGrad - The layout for gradient values
1034 
1035   Level: developer
1036 
1037 .seealso: DMPlexSNESGetGeometryFVM()
1038 @*/
1039 PetscErrorCode DMPlexSNESGetGradientDM(DM dm, PetscFV fv, DM *dmGrad)
1040 {
1041   DM             plex;
1042   PetscBool      computeGradients;
1043   PetscErrorCode ierr;
1044 
1045   PetscFunctionBegin;
1046   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1047   PetscValidHeaderSpecific(fv,PETSCFV_CLASSID,2);
1048   PetscValidPointer(dmGrad,3);
1049   ierr = PetscFVGetComputeGradients(fv, &computeGradients);CHKERRQ(ierr);
1050   if (!computeGradients) {*dmGrad = NULL; PetscFunctionReturn(0);}
1051   ierr = DMSNESConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
1052   ierr = DMPlexGetDataFVM(plex, fv, NULL, NULL, dmGrad);CHKERRQ(ierr);
1053   ierr = DMDestroy(&plex);CHKERRQ(ierr);
1054   PetscFunctionReturn(0);
1055 }
1056 
1057 static PetscErrorCode DMPlexComputeBdResidual_Single_Internal(DM dm, PetscReal t, DMLabel label, PetscInt numValues, const PetscInt values[], PetscInt field, Vec locX, Vec locX_t, Vec locF, DMField coordField, IS facetIS)
1058 {
1059   DM_Plex         *mesh = (DM_Plex *) dm->data;
1060   DM               plex = NULL, plexA = NULL;
1061   PetscDS          prob, probAux = NULL;
1062   PetscSection     section, sectionAux = NULL;
1063   Vec              locA = NULL;
1064   PetscScalar     *u = NULL, *u_t = NULL, *a = NULL, *elemVec = NULL;
1065   PetscInt         v;
1066   PetscInt         totDim, totDimAux = 0;
1067   PetscErrorCode   ierr;
1068 
1069   PetscFunctionBegin;
1070   ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr);
1071   ierr = DMGetSection(dm, &section);CHKERRQ(ierr);
1072   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
1073   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
1074   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &locA);CHKERRQ(ierr);
1075   if (locA) {
1076     DM dmAux;
1077 
1078     ierr = VecGetDM(locA, &dmAux);CHKERRQ(ierr);
1079     ierr = DMConvert(dmAux, DMPLEX, &plexA);CHKERRQ(ierr);
1080     ierr = DMGetDS(plexA, &probAux);CHKERRQ(ierr);
1081     ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
1082     ierr = DMGetSection(plexA, &sectionAux);CHKERRQ(ierr);
1083   }
1084   for (v = 0; v < numValues; ++v) {
1085     PetscFEGeom    *fgeom;
1086     PetscInt        maxDegree;
1087     PetscQuadrature qGeom = NULL;
1088     IS              pointIS;
1089     const PetscInt *points;
1090     PetscInt        numFaces, face, Nq;
1091 
1092     ierr = DMLabelGetStratumIS(label, values[v], &pointIS);CHKERRQ(ierr);
1093     if (!pointIS) continue; /* No points with that id on this process */
1094     {
1095       IS isectIS;
1096 
1097       /* TODO: Special cases of ISIntersect where it is quick to check a priori if one is a superset of the other */
1098       ierr = ISIntersect_Caching_Internal(facetIS,pointIS,&isectIS);CHKERRQ(ierr);
1099       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
1100       pointIS = isectIS;
1101     }
1102     ierr = ISGetLocalSize(pointIS,&numFaces);CHKERRQ(ierr);
1103     ierr = ISGetIndices(pointIS,&points);CHKERRQ(ierr);
1104     ierr = PetscMalloc4(numFaces*totDim, &u, locX_t ? numFaces*totDim : 0, &u_t, numFaces*totDim, &elemVec, locA ? numFaces*totDimAux : 0, &a);CHKERRQ(ierr);
1105     ierr = DMFieldGetDegree(coordField,pointIS,NULL,&maxDegree);CHKERRQ(ierr);
1106     if (maxDegree <= 1) {
1107       ierr = DMFieldCreateDefaultQuadrature(coordField,pointIS,&qGeom);CHKERRQ(ierr);
1108     }
1109     if (!qGeom) {
1110       PetscFE fe;
1111 
1112       ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);CHKERRQ(ierr);
1113       ierr = PetscFEGetFaceQuadrature(fe, &qGeom);CHKERRQ(ierr);
1114       ierr = PetscObjectReference((PetscObject)qGeom);CHKERRQ(ierr);
1115     }
1116     ierr = PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);
1117     ierr = DMSNESGetFEGeom(coordField,pointIS,qGeom,PETSC_TRUE,&fgeom);CHKERRQ(ierr);
1118     for (face = 0; face < numFaces; ++face) {
1119       const PetscInt point = points[face], *support, *cone;
1120       PetscScalar   *x     = NULL;
1121       PetscInt       i, coneSize, faceLoc;
1122 
1123       ierr = DMPlexGetSupport(dm, point, &support);CHKERRQ(ierr);
1124       ierr = DMPlexGetConeSize(dm, support[0], &coneSize);CHKERRQ(ierr);
1125       ierr = DMPlexGetCone(dm, support[0], &cone);CHKERRQ(ierr);
1126       for (faceLoc = 0; faceLoc < coneSize; ++faceLoc) if (cone[faceLoc] == point) break;
1127       if (faceLoc == coneSize) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not find face %D in cone of support[0] %D", point, support[0]);
1128       fgeom->face[face][0] = faceLoc;
1129       ierr = DMPlexVecGetClosure(plex, section, locX, support[0], NULL, &x);CHKERRQ(ierr);
1130       for (i = 0; i < totDim; ++i) u[face*totDim+i] = x[i];
1131       ierr = DMPlexVecRestoreClosure(plex, section, locX, support[0], NULL, &x);CHKERRQ(ierr);
1132       if (locX_t) {
1133         ierr = DMPlexVecGetClosure(plex, section, locX_t, support[0], NULL, &x);CHKERRQ(ierr);
1134         for (i = 0; i < totDim; ++i) u_t[face*totDim+i] = x[i];
1135         ierr = DMPlexVecRestoreClosure(plex, section, locX_t, support[0], NULL, &x);CHKERRQ(ierr);
1136       }
1137       if (locA) {
1138         PetscInt subp;
1139 
1140         ierr = DMPlexGetAuxiliaryPoint(plex, plexA, support[0], &subp);CHKERRQ(ierr);
1141         ierr = DMPlexVecGetClosure(plexA, sectionAux, locA, subp, NULL, &x);CHKERRQ(ierr);
1142         for (i = 0; i < totDimAux; ++i) a[face*totDimAux+i] = x[i];
1143         ierr = DMPlexVecRestoreClosure(plexA, sectionAux, locA, subp, NULL, &x);CHKERRQ(ierr);
1144       }
1145     }
1146     ierr = PetscMemzero(elemVec, numFaces*totDim * sizeof(PetscScalar));CHKERRQ(ierr);
1147     {
1148       PetscFE         fe;
1149       PetscInt        Nb;
1150       PetscFEGeom     *chunkGeom = NULL;
1151       /* Conforming batches */
1152       PetscInt        numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1153       /* Remainder */
1154       PetscInt        Nr, offset;
1155 
1156       ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);CHKERRQ(ierr);
1157       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
1158       ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
1159       /* TODO: documentation is unclear about what is going on with these numbers: how should Nb / Nq factor in ? */
1160       blockSize = Nb;
1161       batchSize = numBlocks * blockSize;
1162       ierr =  PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1163       numChunks = numFaces / (numBatches*batchSize);
1164       Ne        = numChunks*numBatches*batchSize;
1165       Nr        = numFaces % (numBatches*batchSize);
1166       offset    = numFaces - Nr;
1167       ierr = PetscFEGeomGetChunk(fgeom,0,offset,&chunkGeom);CHKERRQ(ierr);
1168       ierr = PetscFEIntegrateBdResidual(prob, field, Ne, chunkGeom, u, u_t, probAux, a, t, elemVec);CHKERRQ(ierr);
1169       ierr = PetscFEGeomRestoreChunk(fgeom, 0, offset, &chunkGeom);CHKERRQ(ierr);
1170       ierr = PetscFEGeomGetChunk(fgeom,offset,numFaces,&chunkGeom);CHKERRQ(ierr);
1171       ierr = PetscFEIntegrateBdResidual(prob, field, Nr, chunkGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, a ? &a[offset*totDimAux] : NULL, t, &elemVec[offset*totDim]);CHKERRQ(ierr);
1172       ierr = PetscFEGeomRestoreChunk(fgeom,offset,numFaces,&chunkGeom);CHKERRQ(ierr);
1173     }
1174     for (face = 0; face < numFaces; ++face) {
1175       const PetscInt point = points[face], *support;
1176 
1177       if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "BdResidual", totDim, &elemVec[face*totDim]);CHKERRQ(ierr);}
1178       ierr = DMPlexGetSupport(plex, point, &support);CHKERRQ(ierr);
1179       ierr = DMPlexVecSetClosure(plex, NULL, locF, support[0], &elemVec[face*totDim], ADD_ALL_VALUES);CHKERRQ(ierr);
1180     }
1181     ierr = DMSNESRestoreFEGeom(coordField,pointIS,qGeom,PETSC_TRUE,&fgeom);CHKERRQ(ierr);
1182     ierr = PetscQuadratureDestroy(&qGeom);CHKERRQ(ierr);
1183     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
1184     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
1185     ierr = PetscFree4(u, u_t, elemVec, a);CHKERRQ(ierr);
1186   }
1187   if (plex)  {ierr = DMDestroy(&plex);CHKERRQ(ierr);}
1188   if (plexA) {ierr = DMDestroy(&plexA);CHKERRQ(ierr);}
1189   PetscFunctionReturn(0);
1190 }
1191 
1192 PetscErrorCode DMPlexComputeBdResidualSingle(DM dm, PetscReal t, DMLabel label, PetscInt numValues, const PetscInt values[], PetscInt field, Vec locX, Vec locX_t, Vec locF)
1193 {
1194   DMField        coordField;
1195   DMLabel        depthLabel;
1196   IS             facetIS;
1197   PetscInt       dim;
1198   PetscErrorCode ierr;
1199 
1200   PetscFunctionBegin;
1201   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1202   ierr = DMPlexGetDepthLabel(dm, &depthLabel);CHKERRQ(ierr);
1203   ierr = DMLabelGetStratumIS(depthLabel, dim-1, &facetIS);CHKERRQ(ierr);
1204   ierr = DMGetCoordinateField(dm, &coordField);CHKERRQ(ierr);
1205   ierr = DMPlexComputeBdResidual_Single_Internal(dm, t, label, numValues, values, field, locX, locX_t, locF, coordField, facetIS);CHKERRQ(ierr);
1206   PetscFunctionReturn(0);
1207 }
1208 
1209 PetscErrorCode DMPlexComputeBdResidual_Internal(DM dm, Vec locX, Vec locX_t, PetscReal t, Vec locF, void *user)
1210 {
1211   PetscDS        prob;
1212   PetscInt       numBd, bd;
1213   DMField        coordField = NULL;
1214   IS             facetIS    = NULL;
1215   DMLabel        depthLabel;
1216   PetscInt       dim;
1217   PetscErrorCode ierr;
1218 
1219   PetscFunctionBegin;
1220   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
1221   ierr = DMPlexGetDepthLabel(dm, &depthLabel);CHKERRQ(ierr);
1222   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1223   ierr = DMLabelGetStratumIS(depthLabel,dim - 1,&facetIS);CHKERRQ(ierr);
1224   ierr = PetscDSGetNumBoundary(prob, &numBd);CHKERRQ(ierr);
1225   for (bd = 0; bd < numBd; ++bd) {
1226     DMBoundaryConditionType type;
1227     const char             *bdLabel;
1228     DMLabel                 label;
1229     const PetscInt         *values;
1230     PetscInt                field, numValues;
1231     PetscObject             obj;
1232     PetscClassId            id;
1233 
1234     ierr = PetscDSGetBoundary(prob, bd, &type, NULL, &bdLabel, &field, NULL, NULL, NULL, &numValues, &values, NULL);CHKERRQ(ierr);
1235     ierr = PetscDSGetDiscretization(prob, field, &obj);CHKERRQ(ierr);
1236     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1237     if ((id != PETSCFE_CLASSID) || (type & DM_BC_ESSENTIAL)) continue;
1238     if (!facetIS) {
1239       DMLabel  depthLabel;
1240       PetscInt dim;
1241 
1242       ierr = DMPlexGetDepthLabel(dm, &depthLabel);CHKERRQ(ierr);
1243       ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1244       ierr = DMLabelGetStratumIS(depthLabel, dim - 1, &facetIS);CHKERRQ(ierr);
1245     }
1246     ierr = DMGetCoordinateField(dm, &coordField);CHKERRQ(ierr);
1247     ierr = DMGetLabel(dm, bdLabel, &label);CHKERRQ(ierr);
1248     ierr = DMPlexComputeBdResidual_Single_Internal(dm, t, label, numValues, values, field, locX, locX_t, locF, coordField, facetIS);CHKERRQ(ierr);
1249   }
1250   ierr = ISDestroy(&facetIS);CHKERRQ(ierr);
1251   PetscFunctionReturn(0);
1252 }
1253 
1254 PetscErrorCode DMPlexComputeResidual_Internal(DM dm, IS cellIS, PetscReal time, Vec locX, Vec locX_t, PetscReal t, Vec locF, void *user)
1255 {
1256   DM_Plex         *mesh       = (DM_Plex *) dm->data;
1257   const char      *name       = "Residual";
1258   DM               dmAux      = NULL;
1259   DM               dmGrad     = NULL;
1260   DMLabel          ghostLabel = NULL;
1261   PetscDS          prob       = NULL;
1262   PetscDS          probAux    = NULL;
1263   PetscSection     section    = NULL;
1264   PetscBool        useFEM     = PETSC_FALSE;
1265   PetscBool        useFVM     = PETSC_FALSE;
1266   PetscBool        isImplicit = (locX_t || time == PETSC_MIN_REAL) ? PETSC_TRUE : PETSC_FALSE;
1267   PetscFV          fvm        = NULL;
1268   PetscFVCellGeom *cgeomFVM   = NULL;
1269   PetscFVFaceGeom *fgeomFVM   = NULL;
1270   DMField          coordField = NULL;
1271   Vec              locA, cellGeometryFVM = NULL, faceGeometryFVM = NULL, grad, locGrad = NULL;
1272   PetscScalar     *u = NULL, *u_t, *a, *uL, *uR;
1273   IS               chunkIS;
1274   const PetscInt  *cells;
1275   PetscInt         cStart, cEnd, numCells;
1276   PetscInt         Nf, f, totDim, totDimAux, numChunks, cellChunkSize, faceChunkSize, chunk, fStart, fEnd;
1277   PetscInt         maxDegree = PETSC_MAX_INT;
1278   PetscQuadrature  affineQuad = NULL, *quads = NULL;
1279   PetscFEGeom     *affineGeom = NULL, **geoms = NULL;
1280   PetscErrorCode   ierr;
1281 
1282   PetscFunctionBegin;
1283   ierr = PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
1284   /* TODO The places where we have to use isFE are probably the member functions for the PetscDisc class */
1285   /* TODO The FVM geometry is over-manipulated. Make the precalc functions return exactly what we need */
1286   /* FEM+FVM */
1287   ierr = ISGetPointRange(cellIS, &cStart, &cEnd, &cells);CHKERRQ(ierr);
1288   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
1289   /* 1: Get sizes from dm and dmAux */
1290   ierr = DMGetSection(dm, &section);CHKERRQ(ierr);
1291   ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr);
1292   ierr = DMGetCellDS(dm, cStart, &prob);CHKERRQ(ierr);
1293   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
1294   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
1295   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &locA);CHKERRQ(ierr);
1296   if (locA) {
1297     PetscInt subcell;
1298     ierr = DMPlexGetAuxiliaryPoint(dm, dmAux, cStart, &subcell);CHKERRQ(ierr);
1299     ierr = VecGetDM(locA, &dmAux);CHKERRQ(ierr);
1300     ierr = DMGetCellDS(dmAux, subcell, &probAux);CHKERRQ(ierr);
1301     ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
1302   }
1303   /* 2: Get geometric data */
1304   for (f = 0; f < Nf; ++f) {
1305     PetscObject  obj;
1306     PetscClassId id;
1307     PetscBool    fimp;
1308 
1309     ierr = PetscDSGetImplicit(prob, f, &fimp);CHKERRQ(ierr);
1310     if (isImplicit != fimp) continue;
1311     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
1312     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1313     if (id == PETSCFE_CLASSID) {useFEM = PETSC_TRUE;}
1314     if (id == PETSCFV_CLASSID) {useFVM = PETSC_TRUE; fvm = (PetscFV) obj;}
1315   }
1316   if (useFEM) {
1317     ierr = DMGetCoordinateField(dm, &coordField);CHKERRQ(ierr);
1318     ierr = DMFieldGetDegree(coordField,cellIS,NULL,&maxDegree);CHKERRQ(ierr);
1319     if (maxDegree <= 1) {
1320       ierr = DMFieldCreateDefaultQuadrature(coordField,cellIS,&affineQuad);CHKERRQ(ierr);
1321       if (affineQuad) {
1322         ierr = DMSNESGetFEGeom(coordField,cellIS,affineQuad,PETSC_FALSE,&affineGeom);CHKERRQ(ierr);
1323       }
1324     } else {
1325       ierr = PetscCalloc2(Nf,&quads,Nf,&geoms);CHKERRQ(ierr);
1326       for (f = 0; f < Nf; ++f) {
1327         PetscObject  obj;
1328         PetscClassId id;
1329         PetscBool    fimp;
1330 
1331         ierr = PetscDSGetImplicit(prob, f, &fimp);CHKERRQ(ierr);
1332         if (isImplicit != fimp) continue;
1333         ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
1334         ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1335         if (id == PETSCFE_CLASSID) {
1336           PetscFE fe = (PetscFE) obj;
1337 
1338           ierr = PetscFEGetQuadrature(fe, &quads[f]);CHKERRQ(ierr);
1339           ierr = PetscObjectReference((PetscObject)quads[f]);CHKERRQ(ierr);
1340           ierr = DMSNESGetFEGeom(coordField,cellIS,quads[f],PETSC_FALSE,&geoms[f]);CHKERRQ(ierr);
1341         }
1342       }
1343     }
1344   }
1345   if (useFVM) {
1346     ierr = DMPlexSNESGetGeometryFVM(dm, &faceGeometryFVM, &cellGeometryFVM, NULL);CHKERRQ(ierr);
1347     ierr = VecGetArrayRead(faceGeometryFVM, (const PetscScalar **) &fgeomFVM);CHKERRQ(ierr);
1348     ierr = VecGetArrayRead(cellGeometryFVM, (const PetscScalar **) &cgeomFVM);CHKERRQ(ierr);
1349     /* Reconstruct and limit cell gradients */
1350     ierr = DMPlexSNESGetGradientDM(dm, fvm, &dmGrad);CHKERRQ(ierr);
1351     if (dmGrad) {
1352       ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
1353       ierr = DMGetGlobalVector(dmGrad, &grad);CHKERRQ(ierr);
1354       ierr = DMPlexReconstructGradients_Internal(dm, fvm, fStart, fEnd, faceGeometryFVM, cellGeometryFVM, locX, grad);CHKERRQ(ierr);
1355       /* Communicate gradient values */
1356       ierr = DMGetLocalVector(dmGrad, &locGrad);CHKERRQ(ierr);
1357       ierr = DMGlobalToLocalBegin(dmGrad, grad, INSERT_VALUES, locGrad);CHKERRQ(ierr);
1358       ierr = DMGlobalToLocalEnd(dmGrad, grad, INSERT_VALUES, locGrad);CHKERRQ(ierr);
1359       ierr = DMRestoreGlobalVector(dmGrad, &grad);CHKERRQ(ierr);
1360     }
1361     /* Handle non-essential (e.g. outflow) boundary values */
1362     ierr = DMPlexInsertBoundaryValues(dm, PETSC_FALSE, locX, time, faceGeometryFVM, cellGeometryFVM, locGrad);CHKERRQ(ierr);
1363   }
1364   /* Loop over chunks */
1365   if (useFEM) {ierr = ISCreate(PETSC_COMM_SELF, &chunkIS);CHKERRQ(ierr);}
1366   numCells      = cEnd - cStart;
1367   numChunks     = 1;
1368   cellChunkSize = numCells/numChunks;
1369   faceChunkSize = (fEnd - fStart)/numChunks;
1370   numChunks     = PetscMin(1,numCells);
1371   for (chunk = 0; chunk < numChunks; ++chunk) {
1372     PetscScalar     *elemVec, *fluxL, *fluxR;
1373     PetscReal       *vol;
1374     PetscFVFaceGeom *fgeom;
1375     PetscInt         cS = cStart+chunk*cellChunkSize, cE = PetscMin(cS+cellChunkSize, cEnd), numCells = cE - cS, c;
1376     PetscInt         fS = fStart+chunk*faceChunkSize, fE = PetscMin(fS+faceChunkSize, fEnd), numFaces = 0, face;
1377 
1378     /* Extract field coefficients */
1379     if (useFEM) {
1380       ierr = ISGetPointSubrange(chunkIS, cS, cE, cells);CHKERRQ(ierr);
1381       ierr = DMPlexGetCellFields(dm, chunkIS, locX, locX_t, locA, &u, &u_t, &a);CHKERRQ(ierr);
1382       ierr = DMGetWorkArray(dm, numCells*totDim, MPIU_SCALAR, &elemVec);CHKERRQ(ierr);
1383       ierr = PetscMemzero(elemVec, numCells*totDim * sizeof(PetscScalar));CHKERRQ(ierr);
1384     }
1385     if (useFVM) {
1386       ierr = DMPlexGetFaceFields(dm, fS, fE, locX, locX_t, faceGeometryFVM, cellGeometryFVM, locGrad, &numFaces, &uL, &uR);CHKERRQ(ierr);
1387       ierr = DMPlexGetFaceGeometry(dm, fS, fE, faceGeometryFVM, cellGeometryFVM, &numFaces, &fgeom, &vol);CHKERRQ(ierr);
1388       ierr = DMGetWorkArray(dm, numFaces*totDim, MPIU_SCALAR, &fluxL);CHKERRQ(ierr);
1389       ierr = DMGetWorkArray(dm, numFaces*totDim, MPIU_SCALAR, &fluxR);CHKERRQ(ierr);
1390       ierr = PetscMemzero(fluxL, numFaces*totDim * sizeof(PetscScalar));CHKERRQ(ierr);
1391       ierr = PetscMemzero(fluxR, numFaces*totDim * sizeof(PetscScalar));CHKERRQ(ierr);
1392     }
1393     /* TODO We will interlace both our field coefficients (u, u_t, uL, uR, etc.) and our output (elemVec, fL, fR). I think this works */
1394     /* Loop over fields */
1395     for (f = 0; f < Nf; ++f) {
1396       PetscObject  obj;
1397       PetscClassId id;
1398       PetscBool    fimp;
1399       PetscInt     numChunks, numBatches, batchSize, numBlocks, blockSize, Ne, Nr, offset;
1400 
1401       ierr = PetscDSGetImplicit(prob, f, &fimp);CHKERRQ(ierr);
1402       if (isImplicit != fimp) continue;
1403       ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
1404       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1405       if (id == PETSCFE_CLASSID) {
1406         PetscFE         fe = (PetscFE) obj;
1407         PetscFEGeom    *geom = affineGeom ? affineGeom : geoms[f];
1408         PetscFEGeom    *chunkGeom = NULL;
1409         PetscQuadrature quad = affineQuad ? affineQuad : quads[f];
1410         PetscInt        Nq, Nb;
1411 
1412         ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
1413         ierr = PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);
1414         ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
1415         blockSize = Nb;
1416         batchSize = numBlocks * blockSize;
1417         ierr      = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1418         numChunks = numCells / (numBatches*batchSize);
1419         Ne        = numChunks*numBatches*batchSize;
1420         Nr        = numCells % (numBatches*batchSize);
1421         offset    = numCells - Nr;
1422         /* Integrate FE residual to get elemVec (need fields at quadrature points) */
1423         /*   For FV, I think we use a P0 basis and the cell coefficients (for subdivided cells, we can tweak the basis tabulation to be the indicator function) */
1424         ierr = PetscFEGeomGetChunk(geom,0,offset,&chunkGeom);CHKERRQ(ierr);
1425         ierr = PetscFEIntegrateResidual(prob, f, Ne, chunkGeom, u, u_t, probAux, a, t, elemVec);CHKERRQ(ierr);
1426         ierr = PetscFEGeomGetChunk(geom,offset,numCells,&chunkGeom);CHKERRQ(ierr);
1427         ierr = PetscFEIntegrateResidual(prob, f, Nr, chunkGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], t, &elemVec[offset*totDim]);CHKERRQ(ierr);
1428         ierr = PetscFEGeomRestoreChunk(geom,offset,numCells,&chunkGeom);CHKERRQ(ierr);
1429       } else if (id == PETSCFV_CLASSID) {
1430         PetscFV fv = (PetscFV) obj;
1431 
1432         Ne = numFaces;
1433         /* Riemann solve over faces (need fields at face centroids) */
1434         /*   We need to evaluate FE fields at those coordinates */
1435         ierr = PetscFVIntegrateRHSFunction(fv, prob, f, Ne, fgeom, vol, uL, uR, fluxL, fluxR);CHKERRQ(ierr);
1436       } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
1437     }
1438     /* Loop over domain */
1439     if (useFEM) {
1440       /* Add elemVec to locX */
1441       for (c = cS; c < cE; ++c) {
1442         const PetscInt cell = cells ? cells[c] : c;
1443         const PetscInt cind = c - cStart;
1444 
1445         if (mesh->printFEM > 1) {ierr = DMPrintCellVector(cell, name, totDim, &elemVec[cind*totDim]);CHKERRQ(ierr);}
1446         if (ghostLabel) {
1447           PetscInt ghostVal;
1448 
1449           ierr = DMLabelGetValue(ghostLabel,cell,&ghostVal);CHKERRQ(ierr);
1450           if (ghostVal > 0) continue;
1451         }
1452         ierr = DMPlexVecSetClosure(dm, section, locF, cell, &elemVec[cind*totDim], ADD_ALL_VALUES);CHKERRQ(ierr);
1453       }
1454     }
1455     if (useFVM) {
1456       PetscScalar *fa;
1457       PetscInt     iface;
1458 
1459       ierr = VecGetArray(locF, &fa);CHKERRQ(ierr);
1460       for (f = 0; f < Nf; ++f) {
1461         PetscFV      fv;
1462         PetscObject  obj;
1463         PetscClassId id;
1464         PetscInt     foff, pdim;
1465 
1466         ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
1467         ierr = PetscDSGetFieldOffset(prob, f, &foff);CHKERRQ(ierr);
1468         ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1469         if (id != PETSCFV_CLASSID) continue;
1470         fv   = (PetscFV) obj;
1471         ierr = PetscFVGetNumComponents(fv, &pdim);CHKERRQ(ierr);
1472         /* Accumulate fluxes to cells */
1473         for (face = fS, iface = 0; face < fE; ++face) {
1474           const PetscInt *scells;
1475           PetscScalar    *fL = NULL, *fR = NULL;
1476           PetscInt        ghost, d, nsupp, nchild;
1477 
1478           ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr);
1479           ierr = DMPlexGetSupportSize(dm, face, &nsupp);CHKERRQ(ierr);
1480           ierr = DMPlexGetTreeChildren(dm, face, &nchild, NULL);CHKERRQ(ierr);
1481           if (ghost >= 0 || nsupp > 2 || nchild > 0) continue;
1482           ierr = DMPlexGetSupport(dm, face, &scells);CHKERRQ(ierr);
1483           ierr = DMLabelGetValue(ghostLabel,scells[0],&ghost);CHKERRQ(ierr);
1484           if (ghost <= 0) {ierr = DMPlexPointLocalFieldRef(dm, scells[0], f, fa, &fL);CHKERRQ(ierr);}
1485           ierr = DMLabelGetValue(ghostLabel,scells[1],&ghost);CHKERRQ(ierr);
1486           if (ghost <= 0) {ierr = DMPlexPointLocalFieldRef(dm, scells[1], f, fa, &fR);CHKERRQ(ierr);}
1487           for (d = 0; d < pdim; ++d) {
1488             if (fL) fL[d] -= fluxL[iface*totDim+foff+d];
1489             if (fR) fR[d] += fluxR[iface*totDim+foff+d];
1490           }
1491           ++iface;
1492         }
1493       }
1494       ierr = VecRestoreArray(locF, &fa);CHKERRQ(ierr);
1495     }
1496     /* Handle time derivative */
1497     if (locX_t) {
1498       PetscScalar *x_t, *fa;
1499 
1500       ierr = VecGetArray(locF, &fa);CHKERRQ(ierr);
1501       ierr = VecGetArray(locX_t, &x_t);CHKERRQ(ierr);
1502       for (f = 0; f < Nf; ++f) {
1503         PetscFV      fv;
1504         PetscObject  obj;
1505         PetscClassId id;
1506         PetscInt     pdim, d;
1507 
1508         ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
1509         ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1510         if (id != PETSCFV_CLASSID) continue;
1511         fv   = (PetscFV) obj;
1512         ierr = PetscFVGetNumComponents(fv, &pdim);CHKERRQ(ierr);
1513         for (c = cS; c < cE; ++c) {
1514           const PetscInt cell = cells ? cells[c] : c;
1515           PetscScalar   *u_t, *r;
1516 
1517           if (ghostLabel) {
1518             PetscInt ghostVal;
1519 
1520             ierr = DMLabelGetValue(ghostLabel, cell, &ghostVal);CHKERRQ(ierr);
1521             if (ghostVal > 0) continue;
1522           }
1523           ierr = DMPlexPointLocalFieldRead(dm, cell, f, x_t, &u_t);CHKERRQ(ierr);
1524           ierr = DMPlexPointLocalFieldRef(dm, cell, f, fa, &r);CHKERRQ(ierr);
1525           for (d = 0; d < pdim; ++d) r[d] += u_t[d];
1526         }
1527       }
1528       ierr = VecRestoreArray(locX_t, &x_t);CHKERRQ(ierr);
1529       ierr = VecRestoreArray(locF, &fa);CHKERRQ(ierr);
1530     }
1531     if (useFEM) {
1532       ierr = DMPlexRestoreCellFields(dm, chunkIS, locX, locX_t, locA, &u, &u_t, &a);CHKERRQ(ierr);
1533       ierr = DMRestoreWorkArray(dm, numCells*totDim, MPIU_SCALAR, &elemVec);CHKERRQ(ierr);
1534     }
1535     if (useFVM) {
1536       ierr = DMPlexRestoreFaceFields(dm, fS, fE, locX, locX_t, faceGeometryFVM, cellGeometryFVM, locGrad, &numFaces, &uL, &uR);CHKERRQ(ierr);
1537       ierr = DMPlexRestoreFaceGeometry(dm, fS, fE, faceGeometryFVM, cellGeometryFVM, &numFaces, &fgeom, &vol);CHKERRQ(ierr);
1538       ierr = DMRestoreWorkArray(dm, numFaces*totDim, MPIU_SCALAR, &fluxL);CHKERRQ(ierr);
1539       ierr = DMRestoreWorkArray(dm, numFaces*totDim, MPIU_SCALAR, &fluxR);CHKERRQ(ierr);
1540       if (dmGrad) {ierr = DMRestoreLocalVector(dmGrad, &locGrad);CHKERRQ(ierr);}
1541     }
1542   }
1543   if (useFEM) {ierr = ISDestroy(&chunkIS);CHKERRQ(ierr);}
1544   ierr = ISRestorePointRange(cellIS, &cStart, &cEnd, &cells);CHKERRQ(ierr);
1545 
1546   if (useFEM) {
1547     ierr = DMPlexComputeBdResidual_Internal(dm, locX, locX_t, t, locF, user);CHKERRQ(ierr);
1548 
1549     if (maxDegree <= 1) {
1550       ierr = DMSNESRestoreFEGeom(coordField,cellIS,affineQuad,PETSC_FALSE,&affineGeom);CHKERRQ(ierr);
1551       ierr = PetscQuadratureDestroy(&affineQuad);CHKERRQ(ierr);
1552     } else {
1553       for (f = 0; f < Nf; ++f) {
1554         ierr = DMSNESRestoreFEGeom(coordField,cellIS,quads[f],PETSC_FALSE,&geoms[f]);CHKERRQ(ierr);
1555         ierr = PetscQuadratureDestroy(&quads[f]);CHKERRQ(ierr);
1556       }
1557       ierr = PetscFree2(quads,geoms);CHKERRQ(ierr);
1558     }
1559   }
1560 
1561   /* FEM */
1562   /* 1: Get sizes from dm and dmAux */
1563   /* 2: Get geometric data */
1564   /* 3: Handle boundary values */
1565   /* 4: Loop over domain */
1566   /*   Extract coefficients */
1567   /* Loop over fields */
1568   /*   Set tiling for FE*/
1569   /*   Integrate FE residual to get elemVec */
1570   /*     Loop over subdomain */
1571   /*       Loop over quad points */
1572   /*         Transform coords to real space */
1573   /*         Evaluate field and aux fields at point */
1574   /*         Evaluate residual at point */
1575   /*         Transform residual to real space */
1576   /*       Add residual to elemVec */
1577   /* Loop over domain */
1578   /*   Add elemVec to locX */
1579 
1580   /* FVM */
1581   /* Get geometric data */
1582   /* If using gradients */
1583   /*   Compute gradient data */
1584   /*   Loop over domain faces */
1585   /*     Count computational faces */
1586   /*     Reconstruct cell gradient */
1587   /*   Loop over domain cells */
1588   /*     Limit cell gradients */
1589   /* Handle boundary values */
1590   /* Loop over domain faces */
1591   /*   Read out field, centroid, normal, volume for each side of face */
1592   /* Riemann solve over faces */
1593   /* Loop over domain faces */
1594   /*   Accumulate fluxes to cells */
1595   /* TODO Change printFEM to printDisc here */
1596   if (mesh->printFEM) {
1597     Vec         locFbc;
1598     PetscInt    pStart, pEnd, p, maxDof;
1599     PetscScalar *zeroes;
1600 
1601     ierr = VecDuplicate(locF,&locFbc);CHKERRQ(ierr);
1602     ierr = VecCopy(locF,locFbc);CHKERRQ(ierr);
1603     ierr = PetscSectionGetChart(section,&pStart,&pEnd);CHKERRQ(ierr);
1604     ierr = PetscSectionGetMaxDof(section,&maxDof);CHKERRQ(ierr);
1605     ierr = PetscCalloc1(maxDof,&zeroes);CHKERRQ(ierr);
1606     for (p = pStart; p < pEnd; p++) {
1607       ierr = VecSetValuesSection(locFbc,section,p,zeroes,INSERT_BC_VALUES);CHKERRQ(ierr);
1608     }
1609     ierr = PetscFree(zeroes);CHKERRQ(ierr);
1610     ierr = DMPrintLocalVec(dm, name, mesh->printTol, locFbc);CHKERRQ(ierr);
1611     ierr = VecDestroy(&locFbc);CHKERRQ(ierr);
1612   }
1613   ierr = PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
1614   PetscFunctionReturn(0);
1615 }
1616 
1617 /*@
1618   DMPlexSNESComputeResidualFEM - Form the local residual F from the local input X using pointwise functions specified by the user
1619 
1620   Input Parameters:
1621 + dm - The mesh
1622 . X  - Local solution
1623 - user - The user context
1624 
1625   Output Parameter:
1626 . F  - Local output vector
1627 
1628   Level: developer
1629 
1630 .seealso: DMPlexComputeJacobianAction()
1631 @*/
1632 PetscErrorCode DMPlexSNESComputeResidualFEM(DM dm, Vec X, Vec F, void *user)
1633 {
1634   DM             plex;
1635   IS             cellIS;
1636   PetscInt       depth;
1637   PetscErrorCode ierr;
1638 
1639   PetscFunctionBegin;
1640   ierr = DMSNESConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
1641   ierr = DMPlexGetDepth(plex, &depth);CHKERRQ(ierr);
1642   ierr = DMGetStratumIS(plex, "dim", depth, &cellIS);CHKERRQ(ierr);
1643   if (!cellIS) {
1644     ierr = DMGetStratumIS(plex, "depth", depth, &cellIS);CHKERRQ(ierr);
1645   }
1646   ierr = DMPlexComputeResidual_Internal(plex, cellIS, PETSC_MIN_REAL, X, NULL, 0.0, F, user);CHKERRQ(ierr);
1647   ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
1648   ierr = DMDestroy(&plex);CHKERRQ(ierr);
1649   PetscFunctionReturn(0);
1650 }
1651 
1652 /*@
1653   DMPlexSNESComputeBoundaryFEM - Form the boundary values for the local input X
1654 
1655   Input Parameters:
1656 + dm - The mesh
1657 - user - The user context
1658 
1659   Output Parameter:
1660 . X  - Local solution
1661 
1662   Level: developer
1663 
1664 .seealso: DMPlexComputeJacobianAction()
1665 @*/
1666 PetscErrorCode DMPlexSNESComputeBoundaryFEM(DM dm, Vec X, void *user)
1667 {
1668   DM             plex;
1669   PetscErrorCode ierr;
1670 
1671   PetscFunctionBegin;
1672   ierr = DMSNESConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
1673   ierr = DMPlexInsertBoundaryValues(plex, PETSC_TRUE, X, PETSC_MIN_REAL, NULL, NULL, NULL);CHKERRQ(ierr);
1674   ierr = DMDestroy(&plex);CHKERRQ(ierr);
1675   PetscFunctionReturn(0);
1676 }
1677 
1678 PetscErrorCode DMPlexComputeBdJacobian_Single_Internal(DM dm, PetscReal t, DMLabel label, PetscInt numValues, const PetscInt values[], PetscInt fieldI, Vec locX, Vec locX_t, PetscReal X_tShift, Mat Jac, Mat JacP, DMField coordField, IS facetIS)
1679 {
1680   DM_Plex       *mesh = (DM_Plex *) dm->data;
1681   DM             plex = NULL, plexA = NULL, tdm;
1682   PetscDS        prob, probAux = NULL;
1683   PetscSection   section, sectionAux = NULL;
1684   PetscSection   globalSection, subSection = NULL;
1685   Vec            locA = NULL, tv;
1686   PetscScalar   *u = NULL, *u_t = NULL, *a = NULL, *elemMat = NULL;
1687   PetscInt       v;
1688   PetscInt       Nf, totDim, totDimAux = 0;
1689   PetscBool      isMatISP, transform;
1690   PetscErrorCode ierr;
1691 
1692   PetscFunctionBegin;
1693   ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr);
1694   ierr = DMHasBasisTransform(dm, &transform);CHKERRQ(ierr);
1695   ierr = DMGetBasisTransformDM_Internal(dm, &tdm);CHKERRQ(ierr);
1696   ierr = DMGetBasisTransformVec_Internal(dm, &tv);CHKERRQ(ierr);
1697   ierr = DMGetSection(dm, &section);CHKERRQ(ierr);
1698   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
1699   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
1700   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
1701   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &locA);CHKERRQ(ierr);
1702   if (locA) {
1703     DM dmAux;
1704 
1705     ierr = VecGetDM(locA, &dmAux);CHKERRQ(ierr);
1706     ierr = DMConvert(dmAux, DMPLEX, &plexA);CHKERRQ(ierr);
1707     ierr = DMGetDS(plexA, &probAux);CHKERRQ(ierr);
1708     ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
1709     ierr = DMGetSection(plexA, &sectionAux);CHKERRQ(ierr);
1710   }
1711 
1712   ierr = PetscObjectTypeCompare((PetscObject) JacP, MATIS, &isMatISP);CHKERRQ(ierr);
1713   ierr = DMGetGlobalSection(dm, &globalSection);CHKERRQ(ierr);
1714   if (isMatISP) {ierr = DMPlexGetSubdomainSection(dm, &subSection);CHKERRQ(ierr);}
1715   for (v = 0; v < numValues; ++v) {
1716     PetscFEGeom    *fgeom;
1717     PetscInt        maxDegree;
1718     PetscQuadrature qGeom = NULL;
1719     IS              pointIS;
1720     const PetscInt *points;
1721     PetscInt        numFaces, face, Nq;
1722 
1723     ierr = DMLabelGetStratumIS(label, values[v], &pointIS);CHKERRQ(ierr);
1724     if (!pointIS) continue; /* No points with that id on this process */
1725     {
1726       IS isectIS;
1727 
1728       /* TODO: Special cases of ISIntersect where it is quick to check a prior if one is a superset of the other */
1729       ierr = ISIntersect_Caching_Internal(facetIS,pointIS,&isectIS);CHKERRQ(ierr);
1730       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
1731       pointIS = isectIS;
1732     }
1733     ierr = ISGetLocalSize(pointIS, &numFaces);CHKERRQ(ierr);
1734     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
1735     ierr = PetscMalloc4(numFaces*totDim, &u, locX_t ? numFaces*totDim : 0, &u_t, numFaces*totDim*totDim, &elemMat, locA ? numFaces*totDimAux : 0, &a);CHKERRQ(ierr);
1736     ierr = DMFieldGetDegree(coordField,pointIS,NULL,&maxDegree);CHKERRQ(ierr);
1737     if (maxDegree <= 1) {
1738       ierr = DMFieldCreateDefaultQuadrature(coordField,pointIS,&qGeom);CHKERRQ(ierr);
1739     }
1740     if (!qGeom) {
1741       PetscFE fe;
1742 
1743       ierr = PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);CHKERRQ(ierr);
1744       ierr = PetscFEGetFaceQuadrature(fe, &qGeom);CHKERRQ(ierr);
1745       ierr = PetscObjectReference((PetscObject)qGeom);CHKERRQ(ierr);
1746     }
1747     ierr = PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);
1748     ierr = DMSNESGetFEGeom(coordField,pointIS,qGeom,PETSC_TRUE,&fgeom);CHKERRQ(ierr);
1749     for (face = 0; face < numFaces; ++face) {
1750       const PetscInt point = points[face], *support, *cone;
1751       PetscScalar   *x     = NULL;
1752       PetscInt       i, coneSize, faceLoc;
1753 
1754       ierr = DMPlexGetSupport(dm, point, &support);CHKERRQ(ierr);
1755       ierr = DMPlexGetConeSize(dm, support[0], &coneSize);CHKERRQ(ierr);
1756       ierr = DMPlexGetCone(dm, support[0], &cone);CHKERRQ(ierr);
1757       for (faceLoc = 0; faceLoc < coneSize; ++faceLoc) if (cone[faceLoc] == point) break;
1758       if (faceLoc == coneSize) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not find face %d in cone of support[0] %d", point, support[0]);
1759       fgeom->face[face][0] = faceLoc;
1760       ierr = DMPlexVecGetClosure(plex, section, locX, support[0], NULL, &x);CHKERRQ(ierr);
1761       for (i = 0; i < totDim; ++i) u[face*totDim+i] = x[i];
1762       ierr = DMPlexVecRestoreClosure(plex, section, locX, support[0], NULL, &x);CHKERRQ(ierr);
1763       if (locX_t) {
1764         ierr = DMPlexVecGetClosure(plex, section, locX_t, support[0], NULL, &x);CHKERRQ(ierr);
1765         for (i = 0; i < totDim; ++i) u_t[face*totDim+i] = x[i];
1766         ierr = DMPlexVecRestoreClosure(plex, section, locX_t, support[0], NULL, &x);CHKERRQ(ierr);
1767       }
1768       if (locA) {
1769         PetscInt subp;
1770         ierr = DMPlexGetSubpoint(plexA, support[0], &subp);CHKERRQ(ierr);
1771         ierr = DMPlexVecGetClosure(plexA, sectionAux, locA, subp, NULL, &x);CHKERRQ(ierr);
1772         for (i = 0; i < totDimAux; ++i) a[face*totDimAux+i] = x[i];
1773         ierr = DMPlexVecRestoreClosure(plexA, sectionAux, locA, subp, NULL, &x);CHKERRQ(ierr);
1774       }
1775     }
1776     ierr = PetscMemzero(elemMat, numFaces*totDim*totDim * sizeof(PetscScalar));CHKERRQ(ierr);
1777     {
1778       PetscFE         fe;
1779       PetscInt        Nb;
1780       /* Conforming batches */
1781       PetscInt        numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1782       /* Remainder */
1783       PetscFEGeom    *chunkGeom = NULL;
1784       PetscInt        fieldJ, Nr, offset;
1785 
1786       ierr = PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);CHKERRQ(ierr);
1787       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
1788       ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
1789       blockSize = Nb;
1790       batchSize = numBlocks * blockSize;
1791       ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1792       numChunks = numFaces / (numBatches*batchSize);
1793       Ne        = numChunks*numBatches*batchSize;
1794       Nr        = numFaces % (numBatches*batchSize);
1795       offset    = numFaces - Nr;
1796       ierr = PetscFEGeomGetChunk(fgeom,0,offset,&chunkGeom);CHKERRQ(ierr);
1797       for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
1798         ierr = PetscFEIntegrateBdJacobian(prob, fieldI, fieldJ, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMat);CHKERRQ(ierr);
1799       }
1800       ierr = PetscFEGeomGetChunk(fgeom,offset,numFaces,&chunkGeom);CHKERRQ(ierr);
1801       for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
1802         ierr = PetscFEIntegrateBdJacobian(prob, fieldI, fieldJ, Nr, chunkGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, a ? &a[offset*totDimAux] : NULL, t, X_tShift, &elemMat[offset*totDim*totDim]);CHKERRQ(ierr);
1803       }
1804       ierr = PetscFEGeomRestoreChunk(fgeom,offset,numFaces,&chunkGeom);CHKERRQ(ierr);
1805     }
1806     for (face = 0; face < numFaces; ++face) {
1807       const PetscInt point = points[face], *support;
1808 
1809       /* Transform to global basis before insertion in Jacobian */
1810       ierr = DMPlexGetSupport(plex, point, &support);CHKERRQ(ierr);
1811       if (transform) {ierr = DMPlexBasisTransformPointTensor_Internal(dm, tdm, tv, support[0], PETSC_TRUE, totDim, &elemMat[face*totDim*totDim]);CHKERRQ(ierr);}
1812       if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(point, "BdJacobian", totDim, totDim, &elemMat[face*totDim*totDim]);CHKERRQ(ierr);}
1813       if (!isMatISP) {
1814         ierr = DMPlexMatSetClosure(plex, section, globalSection, JacP, support[0], &elemMat[face*totDim*totDim], ADD_VALUES);CHKERRQ(ierr);
1815       } else {
1816         Mat lJ;
1817 
1818         ierr = MatISGetLocalMat(JacP, &lJ);CHKERRQ(ierr);
1819         ierr = DMPlexMatSetClosure(plex, section, subSection, lJ, support[0], &elemMat[face*totDim*totDim], ADD_VALUES);CHKERRQ(ierr);
1820       }
1821     }
1822     ierr = DMSNESRestoreFEGeom(coordField,pointIS,qGeom,PETSC_TRUE,&fgeom);CHKERRQ(ierr);
1823     ierr = PetscQuadratureDestroy(&qGeom);CHKERRQ(ierr);
1824     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
1825     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
1826     ierr = PetscFree4(u, u_t, elemMat, a);CHKERRQ(ierr);
1827   }
1828   if (plex)  {ierr = DMDestroy(&plex);CHKERRQ(ierr);}
1829   if (plexA) {ierr = DMDestroy(&plexA);CHKERRQ(ierr);}
1830   PetscFunctionReturn(0);
1831 }
1832 
1833 PetscErrorCode DMPlexComputeBdJacobianSingle(DM dm, PetscReal t, DMLabel label, PetscInt numValues, const PetscInt values[], PetscInt field, Vec locX, Vec locX_t, PetscReal X_tShift, Mat Jac, Mat JacP)
1834 {
1835   DMField        coordField;
1836   DMLabel        depthLabel;
1837   IS             facetIS;
1838   PetscInt       dim;
1839   PetscErrorCode ierr;
1840 
1841   PetscFunctionBegin;
1842   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1843   ierr = DMPlexGetDepthLabel(dm, &depthLabel);CHKERRQ(ierr);
1844   ierr = DMLabelGetStratumIS(depthLabel, dim-1, &facetIS);CHKERRQ(ierr);
1845   ierr = DMGetCoordinateField(dm, &coordField);CHKERRQ(ierr);
1846   ierr = DMPlexComputeBdJacobian_Single_Internal(dm, t, label, numValues, values, field, locX, locX_t, X_tShift, Jac, JacP, coordField, facetIS);CHKERRQ(ierr);
1847   PetscFunctionReturn(0);
1848 }
1849 
1850 PetscErrorCode DMPlexComputeBdJacobian_Internal(DM dm, Vec locX, Vec locX_t, PetscReal t, PetscReal X_tShift, Mat Jac, Mat JacP, void *user)
1851 {
1852   PetscDS          prob;
1853   PetscInt         dim, numBd, bd;
1854   DMLabel          depthLabel;
1855   DMField          coordField = NULL;
1856   IS               facetIS;
1857   PetscErrorCode   ierr;
1858 
1859   PetscFunctionBegin;
1860   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
1861   ierr = DMPlexGetDepthLabel(dm, &depthLabel);CHKERRQ(ierr);
1862   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1863   ierr = DMLabelGetStratumIS(depthLabel, dim-1, &facetIS);CHKERRQ(ierr);
1864   ierr = PetscDSGetNumBoundary(prob, &numBd);CHKERRQ(ierr);
1865   ierr = DMGetCoordinateField(dm, &coordField);CHKERRQ(ierr);
1866   for (bd = 0; bd < numBd; ++bd) {
1867     DMBoundaryConditionType type;
1868     const char             *bdLabel;
1869     DMLabel                 label;
1870     const PetscInt         *values;
1871     PetscInt                fieldI, numValues;
1872     PetscObject             obj;
1873     PetscClassId            id;
1874 
1875     ierr = PetscDSGetBoundary(prob, bd, &type, NULL, &bdLabel, &fieldI, NULL, NULL, NULL, &numValues, &values, NULL);CHKERRQ(ierr);
1876     ierr = PetscDSGetDiscretization(prob, fieldI, &obj);CHKERRQ(ierr);
1877     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1878     if ((id != PETSCFE_CLASSID) || (type & DM_BC_ESSENTIAL)) continue;
1879     ierr = DMGetLabel(dm, bdLabel, &label);CHKERRQ(ierr);
1880     ierr = DMPlexComputeBdJacobian_Single_Internal(dm, t, label, numValues, values, fieldI, locX, locX_t, X_tShift, Jac, JacP, coordField, facetIS);CHKERRQ(ierr);
1881   }
1882   ierr = ISDestroy(&facetIS);CHKERRQ(ierr);
1883   PetscFunctionReturn(0);
1884 }
1885 
1886 PetscErrorCode DMPlexComputeJacobian_Internal(DM dm, IS cellIS, PetscReal t, PetscReal X_tShift, Vec X, Vec X_t, Mat Jac, Mat JacP,void *user)
1887 {
1888   DM_Plex        *mesh  = (DM_Plex *) dm->data;
1889   const char     *name  = "Jacobian";
1890   DM              dmAux, plex, tdm;
1891   Vec             A, tv;
1892   DMField         coordField;
1893   PetscDS         prob, probAux = NULL;
1894   PetscSection    section, globalSection, subSection, sectionAux;
1895   PetscScalar    *elemMat, *elemMatP, *elemMatD, *u, *u_t, *a = NULL;
1896   const PetscInt *cells;
1897   PetscInt        Nf, fieldI, fieldJ;
1898   PetscInt        totDim, totDimAux, cStart, cEnd, numCells, c;
1899   PetscBool       isMatIS, isMatISP, hasJac, hasPrec, hasDyn, hasFV = PETSC_FALSE, transform;
1900   PetscErrorCode  ierr;
1901 
1902   PetscFunctionBegin;
1903   ierr = PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
1904   ierr = DMHasBasisTransform(dm, &transform);CHKERRQ(ierr);
1905   ierr = DMGetBasisTransformDM_Internal(dm, &tdm);CHKERRQ(ierr);
1906   ierr = DMGetBasisTransformVec_Internal(dm, &tv);CHKERRQ(ierr);
1907   ierr = DMGetSection(dm, &section);CHKERRQ(ierr);
1908   ierr = PetscObjectTypeCompare((PetscObject) JacP, MATIS, &isMatISP);CHKERRQ(ierr);
1909   ierr = DMGetGlobalSection(dm, &globalSection);CHKERRQ(ierr);
1910   if (isMatISP) {ierr = DMPlexGetSubdomainSection(dm, &subSection);CHKERRQ(ierr);}
1911   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
1912   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
1913   ierr = PetscDSHasJacobian(prob, &hasJac);CHKERRQ(ierr);
1914   ierr = PetscDSHasJacobianPreconditioner(prob, &hasPrec);CHKERRQ(ierr);
1915   ierr = PetscDSHasDynamicJacobian(prob, &hasDyn);CHKERRQ(ierr);
1916   hasDyn = hasDyn && (X_tShift != 0.0) ? PETSC_TRUE : PETSC_FALSE;
1917   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
1918   ierr = ISGetLocalSize(cellIS, &numCells);CHKERRQ(ierr);
1919   ierr = ISGetPointRange(cellIS, &cStart, &cEnd, &cells);CHKERRQ(ierr);
1920   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
1921   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
1922   if (dmAux) {
1923     ierr = DMConvert(dmAux, DMPLEX, &plex);CHKERRQ(ierr);
1924     ierr = DMGetSection(plex, &sectionAux);CHKERRQ(ierr);
1925     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
1926     ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
1927   }
1928   ierr = PetscMalloc5(numCells*totDim,&u,X_t ? numCells*totDim : 0,&u_t,hasJac ? numCells*totDim*totDim : 0,&elemMat,hasPrec ? numCells*totDim*totDim : 0, &elemMatP,hasDyn ? numCells*totDim*totDim : 0, &elemMatD);CHKERRQ(ierr);
1929   if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);}
1930   ierr = DMGetCoordinateField(dm, &coordField);CHKERRQ(ierr);
1931   for (c = cStart; c < cEnd; ++c) {
1932     const PetscInt cell = cells ? cells[c] : c;
1933     const PetscInt cind = c - cStart;
1934     PetscScalar   *x = NULL,  *x_t = NULL;
1935     PetscInt       i;
1936 
1937     ierr = DMPlexVecGetClosure(dm, section, X, cell, NULL, &x);CHKERRQ(ierr);
1938     for (i = 0; i < totDim; ++i) u[cind*totDim+i] = x[i];
1939     ierr = DMPlexVecRestoreClosure(dm, section, X, cell, NULL, &x);CHKERRQ(ierr);
1940     if (X_t) {
1941       ierr = DMPlexVecGetClosure(dm, section, X_t, cell, NULL, &x_t);CHKERRQ(ierr);
1942       for (i = 0; i < totDim; ++i) u_t[cind*totDim+i] = x_t[i];
1943       ierr = DMPlexVecRestoreClosure(dm, section, X_t, cell, NULL, &x_t);CHKERRQ(ierr);
1944     }
1945     if (dmAux) {
1946       PetscInt subcell;
1947       ierr = DMPlexGetAuxiliaryPoint(dm, dmAux, cell, &subcell);CHKERRQ(ierr);
1948       ierr = DMPlexVecGetClosure(plex, sectionAux, A, subcell, NULL, &x);CHKERRQ(ierr);
1949       for (i = 0; i < totDimAux; ++i) a[cind*totDimAux+i] = x[i];
1950       ierr = DMPlexVecRestoreClosure(plex, sectionAux, A, subcell, NULL, &x);CHKERRQ(ierr);
1951     }
1952   }
1953   if (hasJac)  {ierr = PetscMemzero(elemMat,  numCells*totDim*totDim * sizeof(PetscScalar));CHKERRQ(ierr);}
1954   if (hasPrec) {ierr = PetscMemzero(elemMatP, numCells*totDim*totDim * sizeof(PetscScalar));CHKERRQ(ierr);}
1955   if (hasDyn)  {ierr = PetscMemzero(elemMatD, numCells*totDim*totDim * sizeof(PetscScalar));CHKERRQ(ierr);}
1956   for (fieldI = 0; fieldI < Nf; ++fieldI) {
1957     PetscClassId    id;
1958     PetscFE         fe;
1959     PetscQuadrature qGeom = NULL;
1960     PetscInt        Nb;
1961     /* Conforming batches */
1962     PetscInt        numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1963     /* Remainder */
1964     PetscInt        Nr, offset, Nq;
1965     PetscInt        maxDegree;
1966     PetscFEGeom     *cgeomFEM, *chunkGeom = NULL, *remGeom = NULL;
1967 
1968     ierr = PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);CHKERRQ(ierr);
1969     ierr = PetscObjectGetClassId((PetscObject) fe, &id);CHKERRQ(ierr);
1970     if (id == PETSCFV_CLASSID) {hasFV = PETSC_TRUE; continue;}
1971     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
1972     ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
1973     ierr = DMFieldGetDegree(coordField,cellIS,NULL,&maxDegree);CHKERRQ(ierr);
1974     if (maxDegree <= 1) {
1975       ierr = DMFieldCreateDefaultQuadrature(coordField,cellIS,&qGeom);CHKERRQ(ierr);
1976     }
1977     if (!qGeom) {
1978       ierr = PetscFEGetQuadrature(fe,&qGeom);CHKERRQ(ierr);
1979       ierr = PetscObjectReference((PetscObject)qGeom);CHKERRQ(ierr);
1980     }
1981     ierr = PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);
1982     ierr = DMSNESGetFEGeom(coordField,cellIS,qGeom,PETSC_FALSE,&cgeomFEM);CHKERRQ(ierr);
1983     blockSize = Nb;
1984     batchSize = numBlocks * blockSize;
1985     ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1986     numChunks = numCells / (numBatches*batchSize);
1987     Ne        = numChunks*numBatches*batchSize;
1988     Nr        = numCells % (numBatches*batchSize);
1989     offset    = numCells - Nr;
1990     ierr = PetscFEGeomGetChunk(cgeomFEM,0,offset,&chunkGeom);CHKERRQ(ierr);
1991     ierr = PetscFEGeomGetChunk(cgeomFEM,offset,numCells,&remGeom);CHKERRQ(ierr);
1992     for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
1993       if (hasJac) {
1994         ierr = PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN, fieldI, fieldJ, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMat);CHKERRQ(ierr);
1995         ierr = PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN, fieldI, fieldJ, Nr, remGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], t, X_tShift, &elemMat[offset*totDim*totDim]);CHKERRQ(ierr);
1996       }
1997       if (hasPrec) {
1998         ierr = PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_PRE, fieldI, fieldJ, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMatP);CHKERRQ(ierr);
1999         ierr = PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_PRE, fieldI, fieldJ, Nr, remGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], t, X_tShift, &elemMatP[offset*totDim*totDim]);CHKERRQ(ierr);
2000       }
2001       if (hasDyn) {
2002         ierr = PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_DYN, fieldI, fieldJ, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMatD);CHKERRQ(ierr);
2003         ierr = PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_DYN, fieldI, fieldJ, Nr, remGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], t, X_tShift, &elemMatD[offset*totDim*totDim]);CHKERRQ(ierr);
2004       }
2005     }
2006     ierr = PetscFEGeomRestoreChunk(cgeomFEM,offset,numCells,&remGeom);CHKERRQ(ierr);
2007     ierr = PetscFEGeomRestoreChunk(cgeomFEM,0,offset,&chunkGeom);CHKERRQ(ierr);
2008     ierr = DMSNESRestoreFEGeom(coordField,cellIS,qGeom,PETSC_FALSE,&cgeomFEM);CHKERRQ(ierr);
2009     ierr = PetscQuadratureDestroy(&qGeom);CHKERRQ(ierr);
2010   }
2011   /*   Add contribution from X_t */
2012   if (hasDyn) {for (c = 0; c < numCells*totDim*totDim; ++c) elemMat[c] += X_tShift*elemMatD[c];}
2013   if (hasFV) {
2014     PetscClassId id;
2015     PetscFV      fv;
2016     PetscInt     offsetI, NcI, NbI = 1, fc, f;
2017 
2018     for (fieldI = 0; fieldI < Nf; ++fieldI) {
2019       ierr = PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fv);CHKERRQ(ierr);
2020       ierr = PetscDSGetFieldOffset(prob, fieldI, &offsetI);CHKERRQ(ierr);
2021       ierr = PetscObjectGetClassId((PetscObject) fv, &id);CHKERRQ(ierr);
2022       if (id != PETSCFV_CLASSID) continue;
2023       /* Put in the identity */
2024       ierr = PetscFVGetNumComponents(fv, &NcI);CHKERRQ(ierr);
2025       for (c = cStart; c < cEnd; ++c) {
2026         const PetscInt cind    = c - cStart;
2027         const PetscInt eOffset = cind*totDim*totDim;
2028         for (fc = 0; fc < NcI; ++fc) {
2029           for (f = 0; f < NbI; ++f) {
2030             const PetscInt i = offsetI + f*NcI+fc;
2031             if (hasPrec) {
2032               if (hasJac) {elemMat[eOffset+i*totDim+i] = 1.0;}
2033               elemMatP[eOffset+i*totDim+i] = 1.0;
2034             } else {elemMat[eOffset+i*totDim+i] = 1.0;}
2035           }
2036         }
2037       }
2038     }
2039     /* No allocated space for FV stuff, so ignore the zero entries */
2040     ierr = MatSetOption(JacP, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE);CHKERRQ(ierr);
2041   }
2042   /* Insert values into matrix */
2043   isMatIS = PETSC_FALSE;
2044   if (hasPrec && hasJac) {
2045     ierr = PetscObjectTypeCompare((PetscObject) JacP, MATIS, &isMatIS);CHKERRQ(ierr);
2046   }
2047   if (isMatIS && !subSection) {
2048     ierr = DMPlexGetSubdomainSection(dm, &subSection);CHKERRQ(ierr);
2049   }
2050   for (c = cStart; c < cEnd; ++c) {
2051     const PetscInt cell = cells ? cells[c] : c;
2052     const PetscInt cind = c - cStart;
2053 
2054     /* Transform to global basis before insertion in Jacobian */
2055     if (transform) {ierr = DMPlexBasisTransformPointTensor_Internal(dm, tdm, tv, cell, PETSC_TRUE, totDim, &elemMat[cind*totDim*totDim]);CHKERRQ(ierr);}
2056     if (hasPrec) {
2057       if (hasJac) {
2058         if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(cell, name, totDim, totDim, &elemMat[cind*totDim*totDim]);CHKERRQ(ierr);}
2059         if (!isMatIS) {
2060           ierr = DMPlexMatSetClosure(dm, section, globalSection, Jac, cell, &elemMat[cind*totDim*totDim], ADD_VALUES);CHKERRQ(ierr);
2061         } else {
2062           Mat lJ;
2063 
2064           ierr = MatISGetLocalMat(Jac,&lJ);CHKERRQ(ierr);
2065           ierr = DMPlexMatSetClosure(dm, section, subSection, lJ, cell, &elemMat[cind*totDim*totDim], ADD_VALUES);CHKERRQ(ierr);
2066         }
2067       }
2068       if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(cell, name, totDim, totDim, &elemMatP[cind*totDim*totDim]);CHKERRQ(ierr);}
2069       if (!isMatISP) {
2070         ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, cell, &elemMatP[cind*totDim*totDim], ADD_VALUES);CHKERRQ(ierr);
2071       } else {
2072         Mat lJ;
2073 
2074         ierr = MatISGetLocalMat(JacP,&lJ);CHKERRQ(ierr);
2075         ierr = DMPlexMatSetClosure(dm, section, subSection, lJ, cell, &elemMatP[cind*totDim*totDim], ADD_VALUES);CHKERRQ(ierr);
2076       }
2077     } else {
2078       if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(cell, name, totDim, totDim, &elemMat[cind*totDim*totDim]);CHKERRQ(ierr);}
2079       if (!isMatISP) {
2080         ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, cell, &elemMat[cind*totDim*totDim], ADD_VALUES);CHKERRQ(ierr);
2081       } else {
2082         Mat lJ;
2083 
2084         ierr = MatISGetLocalMat(JacP,&lJ);CHKERRQ(ierr);
2085         ierr = DMPlexMatSetClosure(dm, section, subSection, lJ, cell, &elemMat[cind*totDim*totDim], ADD_VALUES);CHKERRQ(ierr);
2086       }
2087     }
2088   }
2089   ierr = ISRestorePointRange(cellIS, &cStart, &cEnd, &cells);CHKERRQ(ierr);
2090   if (hasFV) {ierr = MatSetOption(JacP, MAT_IGNORE_ZERO_ENTRIES, PETSC_FALSE);CHKERRQ(ierr);}
2091   ierr = PetscFree5(u,u_t,elemMat,elemMatP,elemMatD);CHKERRQ(ierr);
2092   if (dmAux) {
2093     ierr = PetscFree(a);CHKERRQ(ierr);
2094     ierr = DMDestroy(&plex);CHKERRQ(ierr);
2095   }
2096   /* Compute boundary integrals */
2097   ierr = DMPlexComputeBdJacobian_Internal(dm, X, X_t, t, X_tShift, Jac, JacP, user);CHKERRQ(ierr);
2098   /* Assemble matrix */
2099   if (hasJac && hasPrec) {
2100     ierr = MatAssemblyBegin(Jac, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2101     ierr = MatAssemblyEnd(Jac, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2102   }
2103   ierr = MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2104   ierr = MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2105   ierr = PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
2106   PetscFunctionReturn(0);
2107 }
2108 
2109 /*@
2110   DMPlexComputeJacobianAction - Form the local portion of the Jacobian action Z = J(X) Y at the local solution X using pointwise functions specified by the user.
2111 
2112   Input Parameters:
2113 + dm - The mesh
2114 . cellIS -
2115 . t  - The time
2116 . X_tShift - The multiplier for the Jacobian with repsect to X_t
2117 . X  - Local solution vector
2118 . X_t  - Time-derivative of the local solution vector
2119 . Y  - Local input vector
2120 - user - The user context
2121 
2122   Output Parameter:
2123 . Z - Local output vector
2124 
2125   Note:
2126   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
2127   like a GPU, or vectorize on a multicore machine.
2128 
2129   Level: developer
2130 
2131 .seealso: FormFunctionLocal()
2132 @*/
2133 PetscErrorCode DMPlexComputeJacobianAction(DM dm, IS cellIS, PetscReal t, PetscReal X_tShift, Vec X, Vec X_t, Vec Y, Vec Z, void *user)
2134 {
2135   DM_Plex          *mesh  = (DM_Plex *) dm->data;
2136   const char       *name  = "Jacobian";
2137   DM                dmAux, plex, plexAux = NULL;
2138   Vec               A;
2139   PetscDS           prob, probAux = NULL;
2140   PetscQuadrature   quad;
2141   PetscSection      section, globalSection, sectionAux;
2142   PetscScalar      *elemMat, *elemMatD, *u, *u_t, *a = NULL, *y, *z;
2143   PetscInt          Nf, fieldI, fieldJ;
2144   PetscInt          totDim, totDimAux = 0;
2145   const PetscInt   *cells;
2146   PetscInt          cStart, cEnd, numCells, c;
2147   PetscBool         hasDyn;
2148   DMField           coordField;
2149   PetscErrorCode    ierr;
2150 
2151   PetscFunctionBegin;
2152   ierr = PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
2153   ierr = DMSNESConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr);
2154   if (!cellIS) {
2155     PetscInt depth;
2156 
2157     ierr = DMPlexGetDepth(plex, &depth);CHKERRQ(ierr);
2158     ierr = DMGetStratumIS(plex, "dim", depth, &cellIS);CHKERRQ(ierr);
2159     if (!cellIS) {ierr = DMGetStratumIS(plex, "depth", depth, &cellIS);CHKERRQ(ierr);}
2160   } else {
2161     ierr = PetscObjectReference((PetscObject) cellIS);CHKERRQ(ierr);
2162   }
2163   ierr = DMGetSection(dm, &section);CHKERRQ(ierr);
2164   ierr = DMGetGlobalSection(dm, &globalSection);CHKERRQ(ierr);
2165   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
2166   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
2167   ierr = PetscDSHasDynamicJacobian(prob, &hasDyn);CHKERRQ(ierr);
2168   hasDyn = hasDyn && (X_tShift != 0.0) ? PETSC_TRUE : PETSC_FALSE;
2169   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
2170   ierr = ISGetLocalSize(cellIS, &numCells);CHKERRQ(ierr);
2171   ierr = ISGetPointRange(cellIS, &cStart, &cEnd, &cells);CHKERRQ(ierr);
2172   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
2173   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
2174   if (dmAux) {
2175     ierr = DMConvert(dmAux, DMPLEX, &plexAux);CHKERRQ(ierr);
2176     ierr = DMGetSection(plexAux, &sectionAux);CHKERRQ(ierr);
2177     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
2178     ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
2179   }
2180   ierr = VecSet(Z, 0.0);CHKERRQ(ierr);
2181   ierr = PetscMalloc6(numCells*totDim,&u,X_t ? numCells*totDim : 0,&u_t,numCells*totDim*totDim,&elemMat,hasDyn ? numCells*totDim*totDim : 0, &elemMatD,numCells*totDim,&y,totDim,&z);CHKERRQ(ierr);
2182   if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);}
2183   ierr = DMGetCoordinateField(dm, &coordField);CHKERRQ(ierr);
2184   for (c = cStart; c < cEnd; ++c) {
2185     const PetscInt cell = cells ? cells[c] : c;
2186     const PetscInt cind = c - cStart;
2187     PetscScalar   *x = NULL,  *x_t = NULL;
2188     PetscInt       i;
2189 
2190     ierr = DMPlexVecGetClosure(dm, section, X, cell, NULL, &x);CHKERRQ(ierr);
2191     for (i = 0; i < totDim; ++i) u[cind*totDim+i] = x[i];
2192     ierr = DMPlexVecRestoreClosure(dm, section, X, cell, NULL, &x);CHKERRQ(ierr);
2193     if (X_t) {
2194       ierr = DMPlexVecGetClosure(dm, section, X_t, cell, NULL, &x_t);CHKERRQ(ierr);
2195       for (i = 0; i < totDim; ++i) u_t[cind*totDim+i] = x_t[i];
2196       ierr = DMPlexVecRestoreClosure(dm, section, X_t, cell, NULL, &x_t);CHKERRQ(ierr);
2197     }
2198     if (dmAux) {
2199       PetscInt subcell;
2200       ierr = DMPlexGetAuxiliaryPoint(dm, dmAux, cell, &subcell);CHKERRQ(ierr);
2201       ierr = DMPlexVecGetClosure(plexAux, sectionAux, A, subcell, NULL, &x);CHKERRQ(ierr);
2202       for (i = 0; i < totDimAux; ++i) a[cind*totDimAux+i] = x[i];
2203       ierr = DMPlexVecRestoreClosure(plexAux, sectionAux, A, subcell, NULL, &x);CHKERRQ(ierr);
2204     }
2205     ierr = DMPlexVecGetClosure(dm, section, Y, cell, NULL, &x);CHKERRQ(ierr);
2206     for (i = 0; i < totDim; ++i) y[cind*totDim+i] = x[i];
2207     ierr = DMPlexVecRestoreClosure(dm, section, Y, cell, NULL, &x);CHKERRQ(ierr);
2208   }
2209   ierr = PetscMemzero(elemMat, numCells*totDim*totDim * sizeof(PetscScalar));CHKERRQ(ierr);
2210   if (hasDyn)  {ierr = PetscMemzero(elemMatD, numCells*totDim*totDim * sizeof(PetscScalar));CHKERRQ(ierr);}
2211   for (fieldI = 0; fieldI < Nf; ++fieldI) {
2212     PetscFE  fe;
2213     PetscInt Nb;
2214     /* Conforming batches */
2215     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
2216     /* Remainder */
2217     PetscInt Nr, offset, Nq;
2218     PetscQuadrature qGeom = NULL;
2219     PetscInt    maxDegree;
2220     PetscFEGeom *cgeomFEM, *chunkGeom = NULL, *remGeom = NULL;
2221 
2222     ierr = PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);CHKERRQ(ierr);
2223     ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
2224     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
2225     ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
2226     ierr = DMFieldGetDegree(coordField,cellIS,NULL,&maxDegree);CHKERRQ(ierr);
2227     if (maxDegree <= 1) {ierr = DMFieldCreateDefaultQuadrature(coordField,cellIS,&qGeom);CHKERRQ(ierr);}
2228     if (!qGeom) {
2229       ierr = PetscFEGetQuadrature(fe,&qGeom);CHKERRQ(ierr);
2230       ierr = PetscObjectReference((PetscObject)qGeom);CHKERRQ(ierr);
2231     }
2232     ierr = PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);
2233     ierr = DMSNESGetFEGeom(coordField,cellIS,qGeom,PETSC_FALSE,&cgeomFEM);CHKERRQ(ierr);
2234     blockSize = Nb;
2235     batchSize = numBlocks * blockSize;
2236     ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
2237     numChunks = numCells / (numBatches*batchSize);
2238     Ne        = numChunks*numBatches*batchSize;
2239     Nr        = numCells % (numBatches*batchSize);
2240     offset    = numCells - Nr;
2241     ierr = PetscFEGeomGetChunk(cgeomFEM,0,offset,&chunkGeom);CHKERRQ(ierr);
2242     ierr = PetscFEGeomGetChunk(cgeomFEM,offset,numCells,&remGeom);CHKERRQ(ierr);
2243     for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
2244       ierr = PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN, fieldI, fieldJ, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMat);CHKERRQ(ierr);
2245       ierr = PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN, fieldI, fieldJ, Nr, remGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], t, X_tShift, &elemMat[offset*totDim*totDim]);CHKERRQ(ierr);
2246       if (hasDyn) {
2247         ierr = PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_DYN, fieldI, fieldJ, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMatD);CHKERRQ(ierr);
2248         ierr = PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_DYN, fieldI, fieldJ, Nr, remGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], t, X_tShift, &elemMatD[offset*totDim*totDim]);CHKERRQ(ierr);
2249       }
2250     }
2251     ierr = PetscFEGeomRestoreChunk(cgeomFEM,offset,numCells,&remGeom);CHKERRQ(ierr);
2252     ierr = PetscFEGeomRestoreChunk(cgeomFEM,0,offset,&chunkGeom);CHKERRQ(ierr);
2253     ierr = DMSNESRestoreFEGeom(coordField,cellIS,qGeom,PETSC_FALSE,&cgeomFEM);CHKERRQ(ierr);
2254     ierr = PetscQuadratureDestroy(&qGeom);CHKERRQ(ierr);
2255   }
2256   if (hasDyn) {
2257     for (c = 0; c < numCells*totDim*totDim; ++c) elemMat[c] += X_tShift*elemMatD[c];
2258   }
2259   for (c = cStart; c < cEnd; ++c) {
2260     const PetscInt     cell = cells ? cells[c] : c;
2261     const PetscInt     cind = c - cStart;
2262     const PetscBLASInt M = totDim, one = 1;
2263     const PetscScalar  a = 1.0, b = 0.0;
2264 
2265     PetscStackCallBLAS("BLASgemv", BLASgemv_("N", &M, &M, &a, &elemMat[cind*totDim*totDim], &M, &y[cind*totDim], &one, &b, z, &one));
2266     if (mesh->printFEM > 1) {
2267       ierr = DMPrintCellMatrix(c, name, totDim, totDim, &elemMat[cind*totDim*totDim]);CHKERRQ(ierr);
2268       ierr = DMPrintCellVector(c, "Y",  totDim, &y[cind*totDim]);CHKERRQ(ierr);
2269       ierr = DMPrintCellVector(c, "Z",  totDim, z);CHKERRQ(ierr);
2270     }
2271     ierr = DMPlexVecSetClosure(dm, section, Z, cell, z, ADD_VALUES);CHKERRQ(ierr);
2272   }
2273   ierr = PetscFree6(u,u_t,elemMat,elemMatD,y,z);CHKERRQ(ierr);
2274   if (mesh->printFEM) {
2275     ierr = PetscPrintf(PETSC_COMM_WORLD, "Z:\n");CHKERRQ(ierr);
2276     ierr = VecView(Z, NULL);CHKERRQ(ierr);
2277   }
2278   ierr = PetscFree(a);CHKERRQ(ierr);
2279   ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
2280   ierr = DMDestroy(&plexAux);CHKERRQ(ierr);
2281   ierr = DMDestroy(&plex);CHKERRQ(ierr);
2282   ierr = PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
2283   PetscFunctionReturn(0);
2284 }
2285 
2286 /*@
2287   DMPlexSNESComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user.
2288 
2289   Input Parameters:
2290 + dm - The mesh
2291 . X  - Local input vector
2292 - user - The user context
2293 
2294   Output Parameter:
2295 . Jac  - Jacobian matrix
2296 
2297   Note:
2298   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
2299   like a GPU, or vectorize on a multicore machine.
2300 
2301   Level: developer
2302 
2303 .seealso: FormFunctionLocal()
2304 @*/
2305 PetscErrorCode DMPlexSNESComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP,void *user)
2306 {
2307   DM             plex;
2308   PetscDS        prob;
2309   IS             cellIS;
2310   PetscBool      hasJac, hasPrec;
2311   PetscInt       depth;
2312   PetscErrorCode ierr;
2313 
2314   PetscFunctionBegin;
2315   ierr = DMSNESConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
2316   ierr = DMPlexGetDepth(plex, &depth);CHKERRQ(ierr);
2317   ierr = DMGetStratumIS(plex, "dim", depth, &cellIS);CHKERRQ(ierr);
2318   if (!cellIS) {ierr = DMGetStratumIS(plex, "depth", depth, &cellIS);CHKERRQ(ierr);}
2319   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
2320   ierr = PetscDSHasJacobian(prob, &hasJac);CHKERRQ(ierr);
2321   ierr = PetscDSHasJacobianPreconditioner(prob, &hasPrec);CHKERRQ(ierr);
2322   if (hasJac && hasPrec) {ierr = MatZeroEntries(Jac);CHKERRQ(ierr);}
2323   ierr = MatZeroEntries(JacP);CHKERRQ(ierr);
2324   ierr = DMPlexComputeJacobian_Internal(plex, cellIS, 0.0, 0.0, X, NULL, Jac, JacP, user);CHKERRQ(ierr);
2325   ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
2326   ierr = DMDestroy(&plex);CHKERRQ(ierr);
2327   PetscFunctionReturn(0);
2328 }
2329 
2330 /*@
2331   DMPlexSetSNESLocalFEM - Use DMPlex's internal FEM routines to compute SNES boundary values, residual, and Jacobian.
2332 
2333   Input Parameters:
2334 + dm - The DM object
2335 . boundaryctx - the user context that will be passed to pointwise evaluation of boundary values (see PetscDSAddBoundary())
2336 . residualctx - the user context that will be passed to pointwise evaluation of finite element residual computations (see PetscDSSetResidual())
2337 - jacobianctx - the user context that will be passed to pointwise evaluation of finite element Jacobian construction (see PetscDSSetJacobian())
2338 
2339   Level: developer
2340 @*/
2341 PetscErrorCode DMPlexSetSNESLocalFEM(DM dm, void *boundaryctx, void *residualctx, void *jacobianctx)
2342 {
2343   PetscErrorCode ierr;
2344 
2345   PetscFunctionBegin;
2346   ierr = DMSNESSetBoundaryLocal(dm,DMPlexSNESComputeBoundaryFEM,boundaryctx);CHKERRQ(ierr);
2347   ierr = DMSNESSetFunctionLocal(dm,DMPlexSNESComputeResidualFEM,residualctx);CHKERRQ(ierr);
2348   ierr = DMSNESSetJacobianLocal(dm,DMPlexSNESComputeJacobianFEM,jacobianctx);CHKERRQ(ierr);
2349   PetscFunctionReturn(0);
2350 }
2351 
2352 /*@C
2353   DMSNESCheckDiscretization - Check the discretization error of the exact solution
2354 
2355   Input Parameters:
2356 + snes - the SNES object
2357 . dm   - the DM
2358 . u    - a DM vector
2359 . exactFuncs - pointwise functions of the exact solution for each field
2360 . ctxs - contexts for the functions
2361 . tol  - A tolerance for the check, or -1 to print the results instead
2362 
2363   Level: developer
2364 
2365 .seealso: DNSNESCheckFromOptions(), DMSNESCheckResidual(), DMSNESCheckJacobian()
2366 @*/
2367 PetscErrorCode DMSNESCheckDiscretization(SNES snes, DM dm, Vec u, PetscErrorCode (**exactFuncs)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx), void **ctxs, PetscReal tol)
2368 {
2369   PetscErrorCode (**exacts)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
2370   void            **ectxs;
2371   MPI_Comm          comm;
2372   PetscDS           ds;
2373   PetscReal        *error;
2374   PetscInt          Nf, f;
2375   PetscErrorCode    ierr;
2376 
2377   PetscFunctionBegin;
2378   ierr = PetscObjectGetComm((PetscObject) snes, &comm);CHKERRQ(ierr);
2379   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
2380   ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr);
2381   ierr = PetscMalloc3(Nf, &exacts, Nf, &ectxs, PetscMax(1, Nf), &error);CHKERRQ(ierr);
2382   for (f = 0; f < Nf; ++f) {ierr = PetscDSGetExactSolution(ds, f, &exacts[f], &ectxs[f]);CHKERRQ(ierr);}
2383   ierr = DMProjectFunction(dm, 0.0, exactFuncs ? exactFuncs : exacts, ctxs ? ctxs : ectxs, INSERT_ALL_VALUES, u);CHKERRQ(ierr);
2384   ierr = PetscObjectSetName((PetscObject) u, "Exact Solution");CHKERRQ(ierr);
2385   ierr = PetscObjectSetOptionsPrefix((PetscObject) u, "exact_");CHKERRQ(ierr);
2386   ierr = VecViewFromOptions(u, NULL, "-vec_view");CHKERRQ(ierr);
2387   if (Nf > 1) {
2388     ierr = DMComputeL2FieldDiff(dm, 0.0, exactFuncs ? exactFuncs : exacts, ctxs ? ctxs : ectxs, u, error);CHKERRQ(ierr);
2389     if (tol >= 0.0) {
2390       for (f = 0; f < Nf; ++f) {
2391         if (error[f] > tol) SETERRQ3(comm, PETSC_ERR_ARG_WRONG, "L_2 Error %g for field %D exceeds tolerance %g", (double) error[0], f, (double) tol);
2392       }
2393     } else {
2394       ierr = PetscPrintf(comm, "L_2 Error: [");CHKERRQ(ierr);
2395       for (f = 0; f < Nf; ++f) {
2396         if (f) {ierr = PetscPrintf(comm, ", ");CHKERRQ(ierr);}
2397         ierr = PetscPrintf(comm, "%g", (double)error[f]);CHKERRQ(ierr);
2398       }
2399       ierr = PetscPrintf(comm, "]\n");CHKERRQ(ierr);
2400     }
2401   } else {
2402     ierr = DMComputeL2Diff(dm, 0.0, exactFuncs ? exactFuncs : exacts, ctxs ? ctxs : ectxs , u, &error[0]);CHKERRQ(ierr);
2403     if (tol >= 0.0) {
2404       if (error[0] > tol) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "L_2 Error %g exceeds tolerance %g", (double) error[0], (double) tol);
2405     } else {
2406       ierr = PetscPrintf(comm, "L_2 Error: %g\n", (double)error[0]);CHKERRQ(ierr);
2407     }
2408   }
2409   ierr = PetscFree3(exacts, ectxs, error);CHKERRQ(ierr);
2410   PetscFunctionReturn(0);
2411 }
2412 
2413 /*@C
2414   DMSNESCheckResidual - Check the residual of the exact solution
2415 
2416   Input Parameters:
2417 + snes - the SNES object
2418 . dm   - the DM
2419 . u    - a DM vector
2420 . tol  - A tolerance for the check, or -1 to print the results instead
2421 
2422   Level: developer
2423 
2424 .seealso: DNSNESCheckFromOptions(), DMSNESCheckDiscretization(), DMSNESCheckJacobian()
2425 @*/
2426 PetscErrorCode DMSNESCheckResidual(SNES snes, DM dm, Vec u, PetscReal tol)
2427 {
2428   MPI_Comm       comm;
2429   Vec            r;
2430   PetscReal      res;
2431   PetscErrorCode ierr;
2432 
2433   PetscFunctionBegin;
2434   ierr = PetscObjectGetComm((PetscObject) snes, &comm);CHKERRQ(ierr);
2435   ierr = VecDuplicate(u, &r);CHKERRQ(ierr);
2436   ierr = SNESComputeFunction(snes, u, r);CHKERRQ(ierr);
2437   ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr);
2438   if (tol >= 0.0) {
2439     if (res > tol) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "L_2 Residual %g exceeds tolerance %g", (double) res, (double) tol);
2440   } else {
2441     ierr = PetscPrintf(comm, "L_2 Residual: %g\n", (double)res);CHKERRQ(ierr);
2442     ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr);
2443     ierr = PetscObjectSetName((PetscObject) r, "Initial Residual");CHKERRQ(ierr);
2444     ierr = PetscObjectSetOptionsPrefix((PetscObject)r,"res_");CHKERRQ(ierr);
2445     ierr = VecViewFromOptions(r, NULL, "-vec_view");CHKERRQ(ierr);
2446   }
2447   ierr = VecDestroy(&r);CHKERRQ(ierr);
2448   PetscFunctionReturn(0);
2449 }
2450 
2451 /*@C
2452   DMSNESCheckJacobian - Check the Jacobian of the exact solution against the residual, and using the Taylor Test
2453 
2454   Input Parameters:
2455 + snes - the SNES object
2456 . dm   - the DM
2457 . u    - a DM vector
2458 . tol  - A tolerance for the check, or -1 to print the results instead
2459 
2460   Level: developer
2461 
2462 .seealso: DNSNESCheckFromOptions(), DMSNESCheckDiscretization(), DMSNESCheckResidual()
2463 @*/
2464 PetscErrorCode DMSNESCheckJacobian(SNES snes, DM dm, Vec u, PetscReal tol)
2465 {
2466   MPI_Comm       comm;
2467   PetscDS        ds;
2468   Mat            J, M;
2469   MatNullSpace   nullspace;
2470   PetscBool      hasJac, hasPrec;
2471   PetscErrorCode ierr;
2472 
2473   PetscFunctionBegin;
2474   ierr = PetscObjectGetComm((PetscObject) snes, &comm);CHKERRQ(ierr);
2475   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
2476   /* Create and view matrices */
2477   ierr = DMCreateMatrix(dm, &J);CHKERRQ(ierr);
2478   ierr = PetscDSHasJacobian(ds, &hasJac);CHKERRQ(ierr);
2479   ierr = PetscDSHasJacobianPreconditioner(ds, &hasPrec);CHKERRQ(ierr);
2480   if (hasJac && hasPrec) {
2481     ierr = DMCreateMatrix(dm, &M);CHKERRQ(ierr);
2482     ierr = SNESComputeJacobian(snes, u, J, M);CHKERRQ(ierr);
2483     ierr = PetscObjectSetName((PetscObject) M, "Preconditioning Matrix");CHKERRQ(ierr);
2484     ierr = PetscObjectSetOptionsPrefix((PetscObject) M, "jacpre_");CHKERRQ(ierr);
2485     ierr = MatViewFromOptions(M, NULL, "-mat_view");CHKERRQ(ierr);
2486     ierr = MatDestroy(&M);CHKERRQ(ierr);
2487   } else {
2488     ierr = SNESComputeJacobian(snes, u, J, J);CHKERRQ(ierr);
2489   }
2490   ierr = PetscObjectSetName((PetscObject) J, "Jacobian");CHKERRQ(ierr);
2491   ierr = PetscObjectSetOptionsPrefix((PetscObject) J, "jac_");CHKERRQ(ierr);
2492   ierr = MatViewFromOptions(J, NULL, "-mat_view");CHKERRQ(ierr);
2493   /* Check nullspace */
2494   ierr = MatGetNullSpace(J, &nullspace);CHKERRQ(ierr);
2495   if (nullspace) {
2496     PetscBool isNull;
2497     ierr = MatNullSpaceTest(nullspace, J, &isNull);CHKERRQ(ierr);
2498     if (!isNull) SETERRQ(comm, PETSC_ERR_PLIB, "The null space calculated for the system operator is invalid.");
2499   }
2500   ierr = MatNullSpaceDestroy(&nullspace);CHKERRQ(ierr);
2501   /* Taylor test */
2502   {
2503     PetscRandom rand;
2504     Vec         du, uhat, r, rhat, df;
2505     PetscReal   h, slope, intercept;
2506     PetscReal  *es, *hs, *errors;
2507     PetscReal   hMax = 1.0, hMin = 1e-6, hMult = 0.1;
2508     PetscInt    Nv, v;
2509     PetscBool   isLinear = PETSC_FALSE;
2510 
2511     /* Choose a perturbation direction */
2512     ierr = PetscRandomCreate(comm, &rand);CHKERRQ(ierr);
2513     ierr = VecDuplicate(u, &du);CHKERRQ(ierr);
2514     ierr = VecSetRandom(du, rand); CHKERRQ(ierr);
2515     ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr);
2516     ierr = VecDuplicate(u, &df);CHKERRQ(ierr);
2517     ierr = MatMult(J, du, df);CHKERRQ(ierr);
2518     /* Evaluate residual at u, F(u), save in vector r */
2519     ierr = VecDuplicate(u, &r);CHKERRQ(ierr);
2520     ierr = SNESComputeFunction(snes, u, r);CHKERRQ(ierr);
2521     /* Look at the convergence of our Taylor approximation as we approach u */
2522     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv);
2523     ierr = PetscCalloc3(Nv, &es, Nv, &hs, Nv, &errors);CHKERRQ(ierr);
2524     ierr = VecDuplicate(u, &uhat);CHKERRQ(ierr);
2525     ierr = VecDuplicate(u, &rhat);CHKERRQ(ierr);
2526     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv) {
2527       ierr = VecWAXPY(uhat, h, du, u);CHKERRQ(ierr);
2528       /* F(\hat u) \approx F(u) + J(u) (uhat - u) = F(u) + h * J(u) du */
2529       ierr = SNESComputeFunction(snes, uhat, rhat);CHKERRQ(ierr);
2530       ierr = VecAXPBYPCZ(rhat, -1.0, -h, 1.0, r, df);CHKERRQ(ierr);
2531       ierr = VecNorm(rhat, NORM_2, &errors[Nv]);CHKERRQ(ierr);
2532 
2533       es[Nv] = PetscLog10Real(errors[Nv]);
2534       hs[Nv] = PetscLog10Real(h);
2535     }
2536     ierr = VecDestroy(&uhat);CHKERRQ(ierr);
2537     ierr = VecDestroy(&rhat);CHKERRQ(ierr);
2538     ierr = VecDestroy(&df);CHKERRQ(ierr);
2539     ierr = VecDestroy(&r);CHKERRQ(ierr);
2540     ierr = VecDestroy(&du);CHKERRQ(ierr);
2541     for (v = 0; v < Nv; ++v) {
2542       if ((tol >= 0) && (errors[v] > tol)) break;
2543       else if (errors[v] > PETSC_SMALL)    break;
2544     }
2545     if (v == Nv) isLinear = PETSC_TRUE;
2546     ierr = PetscLinearRegression(Nv, hs, es, &slope, &intercept);CHKERRQ(ierr);
2547     ierr = PetscFree3(es, hs, errors);CHKERRQ(ierr);
2548     /* Slope should be about 2 */
2549     if (tol >= 0) {
2550       if (!isLinear && PetscAbsReal(2 - slope) > tol) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Taylor approximation convergence rate should be 2, not %0.2f", (double) slope);
2551     } else {
2552       if (!isLinear) {ierr = PetscPrintf(comm, "Taylor approximation converging at order %3.2f\n", (double) slope);CHKERRQ(ierr);}
2553       else           {ierr = PetscPrintf(comm, "Function appears to be linear\n");CHKERRQ(ierr);}
2554     }
2555   }
2556   ierr = MatDestroy(&J);CHKERRQ(ierr);
2557   PetscFunctionReturn(0);
2558 }
2559 
2560 PetscErrorCode DMSNESCheck_Internal(SNES snes, DM dm, Vec u, PetscErrorCode (**exactFuncs)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx), void **ctxs)
2561 {
2562   PetscErrorCode ierr;
2563 
2564   PetscFunctionBegin;
2565   ierr = DMSNESCheckDiscretization(snes, dm, u, exactFuncs, ctxs, -1.0);CHKERRQ(ierr);
2566   ierr = DMSNESCheckResidual(snes, dm, u, -1.0);CHKERRQ(ierr);
2567   ierr = DMSNESCheckJacobian(snes, dm, u, -1.0);CHKERRQ(ierr);
2568   PetscFunctionReturn(0);
2569 }
2570 
2571 /*@C
2572   DMSNESCheckFromOptions - Check the residual and Jacobian functions using the exact solution by outputting some diagnostic information
2573 
2574   Input Parameters:
2575 + snes - the SNES object
2576 . u    - representative SNES vector
2577 . exactFuncs - pointwise functions of the exact solution for each field
2578 - ctxs - contexts for the functions
2579 
2580   Level: developer
2581 @*/
2582 PetscErrorCode DMSNESCheckFromOptions(SNES snes, Vec u, PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx), void **ctxs)
2583 {
2584   DM             dm;
2585   Vec            sol;
2586   PetscBool      check;
2587   PetscErrorCode ierr;
2588 
2589   PetscFunctionBegin;
2590   ierr = PetscOptionsHasName(((PetscObject)snes)->options,((PetscObject)snes)->prefix, "-dmsnes_check", &check);CHKERRQ(ierr);
2591   if (!check) PetscFunctionReturn(0);
2592   ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr);
2593   ierr = VecDuplicate(u, &sol);CHKERRQ(ierr);
2594   ierr = SNESSetSolution(snes, sol);CHKERRQ(ierr);
2595   ierr = DMSNESCheck_Internal(snes, dm, sol, exactFuncs, ctxs);CHKERRQ(ierr);
2596   ierr = VecDestroy(&sol);CHKERRQ(ierr);
2597   PetscFunctionReturn(0);
2598 }
2599