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