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