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