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