xref: /petsc/src/snes/utils/dmplexsnes.c (revision 3e64cd2f68bf05bac619d1207a7eafd1bfb17961)
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   DM             plex;
908   PetscErrorCode ierr;
909 
910   PetscFunctionBegin;
911   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
912   ierr = DMSNESConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
913   ierr = DMPlexGetDataFVM(plex, NULL, cellgeom, facegeom, NULL);CHKERRQ(ierr);
914   if (minRadius) {ierr = DMPlexGetMinRadius(plex, minRadius);CHKERRQ(ierr);}
915   ierr = DMDestroy(&plex);CHKERRQ(ierr);
916   PetscFunctionReturn(0);
917 }
918 
919 #undef __FUNCT__
920 #define __FUNCT__ "DMPlexSNESGetGradientDM"
921 /*@
922   DMPlexSNESGetGradientDM - Return gradient data layout
923 
924   Input Parameters:
925 + dm - The DM
926 - fv - The PetscFV
927 
928   Output Parameter:
929 . dmGrad - The layout for gradient values
930 
931   Level: developer
932 
933 .seealso: DMPlexSNESGetGeometryFVM()
934 @*/
935 PetscErrorCode DMPlexSNESGetGradientDM(DM dm, PetscFV fv, DM *dmGrad)
936 {
937   DM             plex;
938   PetscBool      computeGradients;
939   PetscErrorCode ierr;
940 
941   PetscFunctionBegin;
942   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
943   PetscValidHeaderSpecific(fv,PETSCFV_CLASSID,2);
944   PetscValidPointer(dmGrad,3);
945   ierr = PetscFVGetComputeGradients(fv, &computeGradients);CHKERRQ(ierr);
946   if (!computeGradients) {*dmGrad = NULL; PetscFunctionReturn(0);}
947   ierr = DMSNESConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
948   ierr = DMPlexGetDataFVM(plex, fv, NULL, NULL, dmGrad);CHKERRQ(ierr);
949   ierr = DMDestroy(&plex);CHKERRQ(ierr);
950   PetscFunctionReturn(0);
951 }
952 
953 #undef __FUNCT__
954 #define __FUNCT__ "DMPlexGetCellFields"
955 /*@C
956   DMPlexGetCellFields - Retrieve the field values values for a chunk of cells
957 
958   Input Parameters:
959 + dm     - The DM
960 . cStart - The first cell to include
961 . cEnd   - The first cell to exclude
962 . locX   - A local vector with the solution fields
963 . locX_t - A local vector with solution field time derivatives, or NULL
964 - locA   - A local vector with auxiliary fields, or NULL
965 
966   Output Parameters:
967 + u   - The field coefficients
968 . u_t - The fields derivative coefficients
969 - a   - The auxiliary field coefficients
970 
971   Level: developer
972 
973 .seealso: DMPlexGetFaceFields()
974 @*/
975 PetscErrorCode DMPlexGetCellFields(DM dm, PetscInt cStart, PetscInt cEnd, Vec locX, Vec locX_t, Vec locA, PetscScalar **u, PetscScalar **u_t, PetscScalar **a)
976 {
977   DM             dmAux;
978   PetscSection   section, sectionAux;
979   PetscDS        prob;
980   PetscInt       numCells = cEnd - cStart, totDim, totDimAux, c;
981   PetscErrorCode ierr;
982 
983   PetscFunctionBegin;
984   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
985   PetscValidHeaderSpecific(locX, VEC_CLASSID, 4);
986   if (locX_t) {PetscValidHeaderSpecific(locX_t, VEC_CLASSID, 5);}
987   if (locA)   {PetscValidHeaderSpecific(locA, VEC_CLASSID, 6);}
988   PetscValidPointer(u, 7);
989   PetscValidPointer(u_t, 8);
990   PetscValidPointer(a, 9);
991   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
992   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
993   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
994   if (locA) {
995     PetscDS probAux;
996 
997     ierr = VecGetDM(locA, &dmAux);CHKERRQ(ierr);
998     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
999     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
1000     ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
1001   }
1002   ierr = DMGetWorkArray(dm, numCells*totDim, PETSC_SCALAR, u);CHKERRQ(ierr);
1003   if (locX_t) {ierr = DMGetWorkArray(dm, numCells*totDim, PETSC_SCALAR, u_t);CHKERRQ(ierr);} else {*u_t = NULL;}
1004   if (locA)   {ierr = DMGetWorkArray(dm, numCells*totDimAux, PETSC_SCALAR, a);CHKERRQ(ierr);} else {*a = NULL;}
1005   for (c = cStart; c < cEnd; ++c) {
1006     PetscScalar *x = NULL, *x_t = NULL, *ul = *u, *ul_t = *u_t, *al = *a;
1007     PetscInt     i;
1008 
1009     ierr = DMPlexVecGetClosure(dm, section, locX, c, NULL, &x);CHKERRQ(ierr);
1010     for (i = 0; i < totDim; ++i) ul[(c-cStart)*totDim+i] = x[i];
1011     ierr = DMPlexVecRestoreClosure(dm, section, locX, c, NULL, &x);CHKERRQ(ierr);
1012     if (locX_t) {
1013       ierr = DMPlexVecGetClosure(dm, section, locX_t, c, NULL, &x_t);CHKERRQ(ierr);
1014       for (i = 0; i < totDim; ++i) ul_t[(c-cStart)*totDim+i] = x_t[i];
1015       ierr = DMPlexVecRestoreClosure(dm, section, locX_t, c, NULL, &x_t);CHKERRQ(ierr);
1016     }
1017     if (locA) {
1018       DM dmAuxPlex;
1019 
1020       ierr = DMSNESConvertPlex(dmAux, &dmAuxPlex, PETSC_FALSE);CHKERRQ(ierr);
1021       ierr = DMPlexVecGetClosure(dmAuxPlex, sectionAux, locA, c, NULL, &x);CHKERRQ(ierr);
1022       for (i = 0; i < totDimAux; ++i) al[(c-cStart)*totDimAux+i] = x[i];
1023       ierr = DMPlexVecRestoreClosure(dmAuxPlex, sectionAux, locA, c, NULL, &x);CHKERRQ(ierr);
1024       ierr = DMDestroy(&dmAuxPlex);CHKERRQ(ierr);
1025     }
1026   }
1027   PetscFunctionReturn(0);
1028 }
1029 
1030 #undef __FUNCT__
1031 #define __FUNCT__ "DMPlexRestoreCellFields"
1032 /*@C
1033   DMPlexRestoreCellFields - Restore the field values values for a chunk of cells
1034 
1035   Input Parameters:
1036 + dm     - The DM
1037 . cStart - The first cell to include
1038 . cEnd   - The first cell to exclude
1039 . locX   - A local vector with the solution fields
1040 . locX_t - A local vector with solution field time derivatives, or NULL
1041 - locA   - A local vector with auxiliary fields, or NULL
1042 
1043   Output Parameters:
1044 + u   - The field coefficients
1045 . u_t - The fields derivative coefficients
1046 - a   - The auxiliary field coefficients
1047 
1048   Level: developer
1049 
1050 .seealso: DMPlexGetFaceFields()
1051 @*/
1052 PetscErrorCode DMPlexRestoreCellFields(DM dm, PetscInt cStart, PetscInt cEnd, Vec locX, Vec locX_t, Vec locA, PetscScalar **u, PetscScalar **u_t, PetscScalar **a)
1053 {
1054   PetscErrorCode ierr;
1055 
1056   PetscFunctionBegin;
1057   ierr = DMRestoreWorkArray(dm, 0, PETSC_SCALAR, u);CHKERRQ(ierr);
1058   if (*u_t) {ierr = DMRestoreWorkArray(dm, 0, PETSC_SCALAR, u_t);CHKERRQ(ierr);}
1059   if (*a)   {ierr = DMRestoreWorkArray(dm, 0, PETSC_SCALAR, a);CHKERRQ(ierr);}
1060   PetscFunctionReturn(0);
1061 }
1062 
1063 #undef __FUNCT__
1064 #define __FUNCT__ "DMPlexGetFaceFields"
1065 /*@C
1066   DMPlexGetFaceFields - Retrieve the field values values for a chunk of faces
1067 
1068   Input Parameters:
1069 + dm     - The DM
1070 . fStart - The first face to include
1071 . fEnd   - The first face to exclude
1072 . locX   - A local vector with the solution fields
1073 . locX_t - A local vector with solution field time derivatives, or NULL
1074 . faceGeometry - A local vector with face geometry
1075 . cellGeometry - A local vector with cell geometry
1076 - locaGrad - A local vector with field gradients, or NULL
1077 
1078   Output Parameters:
1079 + uL - The field values at the left side of the face
1080 - uR - The field values at the right side of the face
1081 
1082   Level: developer
1083 
1084 .seealso: DMPlexGetCellFields()
1085 @*/
1086 PetscErrorCode DMPlexGetFaceFields(DM dm, PetscInt fStart, PetscInt fEnd, Vec locX, Vec locX_t, Vec faceGeometry, Vec cellGeometry, Vec locGrad, PetscScalar **uL, PetscScalar **uR)
1087 {
1088   DM                 dmFace, dmCell, dmGrad = NULL;
1089   PetscSection       section;
1090   PetscDS            prob;
1091   DMLabel            ghostLabel;
1092   const PetscScalar *facegeom, *cellgeom, *x, *lgrad;
1093   PetscBool         *isFE;
1094   PetscInt           dim, Nf, f, Nc, numFaces = fEnd - fStart, iface, face;
1095   PetscErrorCode     ierr;
1096 
1097   PetscFunctionBegin;
1098   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1099   PetscValidHeaderSpecific(locX, VEC_CLASSID, 4);
1100   if (locX_t) {PetscValidHeaderSpecific(locX_t, VEC_CLASSID, 5);}
1101   PetscValidHeaderSpecific(faceGeometry, VEC_CLASSID, 6);
1102   PetscValidHeaderSpecific(cellGeometry, VEC_CLASSID, 7);
1103   if (locGrad) {PetscValidHeaderSpecific(locGrad, VEC_CLASSID, 8);}
1104   PetscValidPointer(uL, 9);
1105   PetscValidPointer(uR, 10);
1106   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1107   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
1108   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1109   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
1110   ierr = PetscDSGetTotalComponents(prob, &Nc);CHKERRQ(ierr);
1111   ierr = PetscMalloc1(Nf, &isFE);CHKERRQ(ierr);
1112   for (f = 0; f < Nf; ++f) {
1113     PetscObject  obj;
1114     PetscClassId id;
1115 
1116     ierr = DMGetField(dm, f, &obj);CHKERRQ(ierr);
1117     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1118     if (id == PETSCFE_CLASSID)      {isFE[f] = PETSC_TRUE;}
1119     else if (id == PETSCFV_CLASSID) {isFE[f] = PETSC_FALSE;}
1120     else                            {isFE[f] = PETSC_FALSE;}
1121   }
1122   ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr);
1123   ierr = VecGetArrayRead(locX, &x);CHKERRQ(ierr);
1124   ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr);
1125   ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
1126   ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr);
1127   ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
1128   if (locGrad) {
1129     ierr = VecGetDM(locGrad, &dmGrad);CHKERRQ(ierr);
1130     ierr = VecGetArrayRead(locGrad, &lgrad);CHKERRQ(ierr);
1131   }
1132   ierr = DMGetWorkArray(dm, numFaces*Nc, PETSC_SCALAR, uL);CHKERRQ(ierr);
1133   ierr = DMGetWorkArray(dm, numFaces*Nc, PETSC_SCALAR, uR);CHKERRQ(ierr);
1134   /* Right now just eat the extra work for FE (could make a cell loop) */
1135   for (face = fStart, iface = 0; face < fEnd; ++face) {
1136     const PetscInt        *cells;
1137     PetscFVFaceGeom       *fg;
1138     PetscFVCellGeom       *cgL, *cgR;
1139     PetscScalar           *xL, *xR, *gL, *gR;
1140     PetscScalar           *uLl = *uL, *uRl = *uR;
1141     PetscInt               ghost, nsupp, nchild;
1142 
1143     ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr);
1144     ierr = DMPlexGetSupportSize(dm, face, &nsupp);CHKERRQ(ierr);
1145     ierr = DMPlexGetTreeChildren(dm, face, &nchild, NULL);CHKERRQ(ierr);
1146     if (ghost >= 0 || nchild > 0) continue;
1147     ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr);
1148     ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr);
1149     ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cgL);CHKERRQ(ierr);
1150     ierr = DMPlexPointLocalRead(dmCell, cells[1], cellgeom, &cgR);CHKERRQ(ierr);
1151     for (f = 0; f < Nf; ++f) {
1152       PetscInt off;
1153 
1154       ierr = PetscDSGetComponentOffset(prob, f, &off);CHKERRQ(ierr);
1155       if (isFE[f]) {
1156         const PetscInt *cone;
1157         PetscInt        comp, coneSizeL, coneSizeR, faceLocL, faceLocR, ldof, rdof, d;
1158 
1159         xL = xR = NULL;
1160         ierr = DMPlexVecGetClosure(dm, section, locX, cells[0], &ldof, (PetscScalar **) &xL);CHKERRQ(ierr);
1161         ierr = DMPlexVecGetClosure(dm, section, locX, cells[1], &rdof, (PetscScalar **) &xR);CHKERRQ(ierr);
1162         ierr = DMPlexGetCone(dm, cells[0], &cone);CHKERRQ(ierr);
1163         ierr = DMPlexGetConeSize(dm, cells[0], &coneSizeL);CHKERRQ(ierr);
1164         for (faceLocL = 0; faceLocL < coneSizeL; ++faceLocL) if (cone[faceLocL] == face) break;
1165         ierr = DMPlexGetCone(dm, cells[1], &cone);CHKERRQ(ierr);
1166         ierr = DMPlexGetConeSize(dm, cells[1], &coneSizeR);CHKERRQ(ierr);
1167         for (faceLocR = 0; faceLocR < coneSizeR; ++faceLocR) if (cone[faceLocR] == face) break;
1168         if (faceLocL == coneSizeL && faceLocR == coneSizeR) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not find face %d in cone of cell %d or cell %d", face, cells[0], cells[1]);
1169         /* Check that FEM field has values in the right cell (sometimes its an FV ghost cell) */
1170         /* TODO: this is a hack that might not be right for nonconforming */
1171         if (faceLocL < coneSizeL) {
1172           ierr = EvaluateFaceFields(prob, f, faceLocL, xL, &uLl[iface*Nc+off]);CHKERRQ(ierr);
1173           if (rdof == ldof && faceLocR < coneSizeR) {ierr = EvaluateFaceFields(prob, f, faceLocR, xR, &uRl[iface*Nc+off]);CHKERRQ(ierr);}
1174           else              {ierr = PetscSectionGetFieldComponents(section, f, &comp);CHKERRQ(ierr); for(d = 0; d < comp; ++d) uRl[iface*Nc+off+d] = uLl[iface*Nc+off+d];}
1175         }
1176         else {
1177           ierr = EvaluateFaceFields(prob, f, faceLocR, xR, &uRl[iface*Nc+off]);CHKERRQ(ierr);
1178           ierr = PetscSectionGetFieldComponents(section, f, &comp);CHKERRQ(ierr);
1179           for(d = 0; d < comp; ++d) uLl[iface*Nc+off+d] = uRl[iface*Nc+off+d];
1180         }
1181         ierr = DMPlexVecRestoreClosure(dm, section, locX, cells[0], &ldof, (PetscScalar **) &xL);CHKERRQ(ierr);
1182         ierr = DMPlexVecRestoreClosure(dm, section, locX, cells[1], &rdof, (PetscScalar **) &xR);CHKERRQ(ierr);
1183       } else {
1184         PetscFV  fv;
1185         PetscInt numComp, c;
1186 
1187         ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fv);CHKERRQ(ierr);
1188         ierr = PetscFVGetNumComponents(fv, &numComp);CHKERRQ(ierr);
1189         if (nsupp > 2) {
1190           for (f = 0; f < Nf; ++f) {
1191             PetscInt off;
1192 
1193             ierr = PetscDSGetComponentOffset(prob, f, &off);CHKERRQ(ierr);
1194             ierr = PetscFVGetNumComponents(fv, &numComp);CHKERRQ(ierr);
1195             for (c = 0; c < numComp; ++c) {
1196               uLl[iface*Nc+off+c] = 0.;
1197               uRl[iface*Nc+off+c] = 0.;
1198             }
1199           }
1200           continue;
1201         }
1202         ierr = DMPlexPointLocalFieldRead(dm, cells[0], f, x, &xL);CHKERRQ(ierr);
1203         ierr = DMPlexPointLocalFieldRead(dm, cells[1], f, x, &xR);CHKERRQ(ierr);
1204         if (dmGrad) {
1205           PetscReal dxL[3], dxR[3];
1206 
1207           ierr = DMPlexPointLocalRead(dmGrad, cells[0], lgrad, &gL);CHKERRQ(ierr);
1208           ierr = DMPlexPointLocalRead(dmGrad, cells[1], lgrad, &gR);CHKERRQ(ierr);
1209           DMPlex_WaxpyD_Internal(dim, -1, cgL->centroid, fg->centroid, dxL);
1210           DMPlex_WaxpyD_Internal(dim, -1, cgR->centroid, fg->centroid, dxR);
1211           for (c = 0; c < numComp; ++c) {
1212             uLl[iface*Nc+off+c] = xL[c] + DMPlex_DotD_Internal(dim, &gL[c*dim], dxL);
1213             uRl[iface*Nc+off+c] = xR[c] + DMPlex_DotD_Internal(dim, &gR[c*dim], dxR);
1214           }
1215         } else {
1216           for (c = 0; c < numComp; ++c) {
1217             uLl[iface*Nc+off+c] = xL[c];
1218             uRl[iface*Nc+off+c] = xR[c];
1219           }
1220         }
1221       }
1222     }
1223     ++iface;
1224   }
1225   ierr = VecRestoreArrayRead(locX, &x);CHKERRQ(ierr);
1226   ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
1227   ierr = VecRestoreArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
1228   if (locGrad) {
1229     ierr = VecRestoreArrayRead(locGrad, &lgrad);CHKERRQ(ierr);
1230   }
1231   ierr = PetscFree(isFE);CHKERRQ(ierr);
1232   PetscFunctionReturn(0);
1233 }
1234 
1235 #undef __FUNCT__
1236 #define __FUNCT__ "DMPlexRestoreFaceFields"
1237 /*@C
1238   DMPlexRestoreFaceFields - Restore the field values values for a chunk of faces
1239 
1240   Input Parameters:
1241 + dm     - The DM
1242 . fStart - The first face to include
1243 . fEnd   - The first face to exclude
1244 . locX   - A local vector with the solution fields
1245 . locX_t - A local vector with solution field time derivatives, or NULL
1246 . faceGeometry - A local vector with face geometry
1247 . cellGeometry - A local vector with cell geometry
1248 - locaGrad - A local vector with field gradients, or NULL
1249 
1250   Output Parameters:
1251 + uL - The field values at the left side of the face
1252 - uR - The field values at the right side of the face
1253 
1254   Level: developer
1255 
1256 .seealso: DMPlexGetFaceFields()
1257 @*/
1258 PetscErrorCode DMPlexRestoreFaceFields(DM dm, PetscInt fStart, PetscInt fEnd, Vec locX, Vec locX_t, Vec faceGeometry, Vec cellGeometry, Vec locGrad, PetscScalar **uL, PetscScalar **uR)
1259 {
1260   PetscErrorCode ierr;
1261 
1262   PetscFunctionBegin;
1263   ierr = DMRestoreWorkArray(dm, 0, PETSC_SCALAR, uL);CHKERRQ(ierr);
1264   ierr = DMRestoreWorkArray(dm, 0, PETSC_SCALAR, uR);CHKERRQ(ierr);
1265   PetscFunctionReturn(0);
1266 }
1267 
1268 #undef __FUNCT__
1269 #define __FUNCT__ "DMPlexGetFaceGeometry"
1270 /*@C
1271   DMPlexGetFaceGeometry - Retrieve the geometric values for a chunk of faces
1272 
1273   Input Parameters:
1274 + dm     - The DM
1275 . fStart - The first face to include
1276 . fEnd   - The first face to exclude
1277 . faceGeometry - A local vector with face geometry
1278 - cellGeometry - A local vector with cell geometry
1279 
1280   Output Parameters:
1281 + fgeom - The extract the face centroid and normal
1282 - vol   - The cell volume
1283 
1284   Level: developer
1285 
1286 .seealso: DMPlexGetCellFields()
1287 @*/
1288 PetscErrorCode DMPlexGetFaceGeometry(DM dm, PetscInt fStart, PetscInt fEnd, Vec faceGeometry, Vec cellGeometry, PetscFVFaceGeom **fgeom, PetscReal **vol)
1289 {
1290   DM                 dmFace, dmCell;
1291   DMLabel            ghostLabel;
1292   const PetscScalar *facegeom, *cellgeom;
1293   PetscInt           dim, numFaces = fEnd - fStart, iface, face;
1294   PetscErrorCode     ierr;
1295 
1296   PetscFunctionBegin;
1297   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1298   PetscValidHeaderSpecific(faceGeometry, VEC_CLASSID, 4);
1299   PetscValidHeaderSpecific(cellGeometry, VEC_CLASSID, 5);
1300   PetscValidPointer(fgeom, 6);
1301   PetscValidPointer(vol, 7);
1302   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1303   ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr);
1304   ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr);
1305   ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
1306   ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr);
1307   ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
1308   ierr = PetscMalloc1(numFaces, fgeom);CHKERRQ(ierr);
1309   ierr = DMGetWorkArray(dm, numFaces*2, PETSC_SCALAR, vol);CHKERRQ(ierr);
1310   for (face = fStart, iface = 0; face < fEnd; ++face) {
1311     const PetscInt        *cells;
1312     PetscFVFaceGeom       *fg;
1313     PetscFVCellGeom       *cgL, *cgR;
1314     PetscFVFaceGeom       *fgeoml = *fgeom;
1315     PetscReal             *voll   = *vol;
1316     PetscInt               ghost, d;
1317 
1318     ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr);
1319     if (ghost >= 0) continue;
1320     ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr);
1321     ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr);
1322     ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cgL);CHKERRQ(ierr);
1323     ierr = DMPlexPointLocalRead(dmCell, cells[1], cellgeom, &cgR);CHKERRQ(ierr);
1324     for (d = 0; d < dim; ++d) {
1325       fgeoml[iface].centroid[d] = fg->centroid[d];
1326       fgeoml[iface].normal[d]   = fg->normal[d];
1327     }
1328     voll[iface*2+0] = cgL->volume;
1329     voll[iface*2+1] = cgR->volume;
1330     ++iface;
1331   }
1332   ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
1333   ierr = VecRestoreArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
1334   PetscFunctionReturn(0);
1335 }
1336 
1337 #undef __FUNCT__
1338 #define __FUNCT__ "DMPlexRestoreFaceGeometry"
1339 /*@C
1340   DMPlexRestoreFaceGeometry - Restore the field values values for a chunk of faces
1341 
1342   Input Parameters:
1343 + dm     - The DM
1344 . fStart - The first face to include
1345 . fEnd   - The first face to exclude
1346 . faceGeometry - A local vector with face geometry
1347 - cellGeometry - A local vector with cell geometry
1348 
1349   Output Parameters:
1350 + fgeom - The extract the face centroid and normal
1351 - vol   - The cell volume
1352 
1353   Level: developer
1354 
1355 .seealso: DMPlexGetFaceFields()
1356 @*/
1357 PetscErrorCode DMPlexRestoreFaceGeometry(DM dm, PetscInt fStart, PetscInt fEnd, Vec faceGeometry, Vec cellGeometry, PetscFVFaceGeom **fgeom, PetscReal **vol)
1358 {
1359   PetscErrorCode ierr;
1360 
1361   PetscFunctionBegin;
1362   ierr = PetscFree(*fgeom);CHKERRQ(ierr);
1363   ierr = DMRestoreWorkArray(dm, 0, PETSC_REAL, vol);CHKERRQ(ierr);
1364   PetscFunctionReturn(0);
1365 }
1366 
1367 #undef __FUNCT__
1368 #define __FUNCT__ "DMPlexComputeBdResidual_Internal"
1369 PetscErrorCode DMPlexComputeBdResidual_Internal(DM dm, Vec locX, Vec locX_t, PetscReal t, Vec locF, void *user)
1370 {
1371   DM_Plex         *mesh = (DM_Plex *) dm->data;
1372   PetscSection     section;
1373   PetscDS          prob;
1374   DMLabel          depth;
1375   PetscFECellGeom *cgeom;
1376   PetscScalar     *u = NULL, *u_t = NULL, *elemVec = NULL;
1377   PetscInt         dim, Nf, f, totDimBd, numBd, bd;
1378   PetscErrorCode   ierr;
1379 
1380   PetscFunctionBegin;
1381   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1382   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1383   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
1384   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
1385   ierr = PetscDSGetTotalBdDimension(prob, &totDimBd);CHKERRQ(ierr);
1386   ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr);
1387   ierr = DMGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
1388   for (bd = 0; bd < numBd; ++bd) {
1389     const char     *bdLabel;
1390     DMLabel         label;
1391     IS              pointIS;
1392     const PetscInt *points;
1393     const PetscInt *values;
1394     PetscInt        field, numValues, v, numPoints, p, dep, numFaces;
1395     PetscBool       isEssential;
1396     PetscObject     obj;
1397     PetscClassId    id;
1398 
1399     ierr = DMGetBoundary(dm, bd, &isEssential, NULL, &bdLabel, &field, NULL, NULL, NULL, &numValues, &values, NULL);CHKERRQ(ierr);
1400     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
1401     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1402     if ((id != PETSCFE_CLASSID) || isEssential) continue;
1403     ierr = DMGetLabel(dm, bdLabel, &label);CHKERRQ(ierr);
1404     for (v = 0; v < numValues; ++v) {
1405       ierr = DMLabelGetStratumSize(label, values[v], &numPoints);CHKERRQ(ierr);
1406       ierr = DMLabelGetStratumIS(label, values[v], &pointIS);CHKERRQ(ierr);
1407       if (!pointIS) continue; /* No points with that id on this process */
1408       ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
1409       for (p = 0, numFaces = 0; p < numPoints; ++p) {
1410         ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
1411         if (dep == dim-1) ++numFaces;
1412       }
1413       ierr = PetscMalloc3(numFaces*totDimBd,&u,numFaces,&cgeom,numFaces*totDimBd,&elemVec);CHKERRQ(ierr);
1414       if (locX_t) {ierr = PetscMalloc1(numFaces*totDimBd,&u_t);CHKERRQ(ierr);}
1415       for (p = 0, f = 0; p < numPoints; ++p) {
1416         const PetscInt point = points[p];
1417         PetscScalar   *x     = NULL;
1418         PetscInt       i;
1419 
1420         ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
1421         if (dep != dim-1) continue;
1422         ierr = DMPlexComputeCellGeometryFEM(dm, point, NULL, cgeom[f].v0, cgeom[f].J, cgeom[f].invJ, &cgeom[f].detJ);CHKERRQ(ierr);
1423         ierr = DMPlexComputeCellGeometryFVM(dm, point, NULL, NULL, cgeom[f].n);CHKERRQ(ierr);
1424         if (cgeom[f].detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", cgeom[f].detJ, point);
1425         /* TODO: Matt, this is wrong if feBd does not match fe: i.e., if the order differs. */
1426         ierr = DMPlexVecGetClosure(dm, section, locX, point, NULL, &x);CHKERRQ(ierr);
1427         for (i = 0; i < totDimBd; ++i) u[f*totDimBd+i] = x[i];
1428         ierr = DMPlexVecRestoreClosure(dm, section, locX, point, NULL, &x);CHKERRQ(ierr);
1429         if (locX_t) {
1430           ierr = DMPlexVecGetClosure(dm, section, locX_t, point, NULL, &x);CHKERRQ(ierr);
1431           for (i = 0; i < totDimBd; ++i) u_t[f*totDimBd+i] = x[i];
1432           ierr = DMPlexVecRestoreClosure(dm, section, locX_t, point, NULL, &x);CHKERRQ(ierr);
1433         }
1434         ++f;
1435       }
1436       for (f = 0; f < Nf; ++f) {
1437         PetscFE         fe;
1438         PetscQuadrature q;
1439         PetscInt        numQuadPoints, Nb;
1440         /* Conforming batches */
1441         PetscInt        numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1442         /* Remainder */
1443         PetscInt        Nr, offset;
1444 
1445         ierr = PetscDSGetBdDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
1446         ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
1447         ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
1448         ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
1449         ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
1450         blockSize = Nb*numQuadPoints;
1451         batchSize = numBlocks * blockSize;
1452         ierr =  PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1453         numChunks = numFaces / (numBatches*batchSize);
1454         Ne        = numChunks*numBatches*batchSize;
1455         Nr        = numFaces % (numBatches*batchSize);
1456         offset    = numFaces - Nr;
1457         ierr = PetscFEIntegrateBdResidual(fe, prob, f, Ne, cgeom, u, u_t, NULL, NULL, t, elemVec);CHKERRQ(ierr);
1458         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);
1459       }
1460       for (p = 0, f = 0; p < numPoints; ++p) {
1461         const PetscInt point = points[p];
1462 
1463         ierr = DMLabelGetValue(depth, point, &dep);CHKERRQ(ierr);
1464         if (dep != dim-1) continue;
1465         if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "BdResidual", totDimBd, &elemVec[f*totDimBd]);CHKERRQ(ierr);}
1466         ierr = DMPlexVecSetClosure(dm, NULL, locF, point, &elemVec[f*totDimBd], ADD_ALL_VALUES);CHKERRQ(ierr);
1467         ++f;
1468       }
1469       ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
1470       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
1471       ierr = PetscFree3(u,cgeom,elemVec);CHKERRQ(ierr);
1472       if (locX_t) {ierr = PetscFree(u_t);CHKERRQ(ierr);}
1473     }
1474   }
1475   PetscFunctionReturn(0);
1476 }
1477 
1478 #undef __FUNCT__
1479 #define __FUNCT__ "DMPlexComputeResidual_Internal"
1480 PetscErrorCode DMPlexComputeResidual_Internal(DM dm, PetscInt cStart, PetscInt cEnd, PetscReal time, Vec locX, Vec locX_t, PetscReal t, Vec locF, void *user)
1481 {
1482   DM_Plex          *mesh       = (DM_Plex *) dm->data;
1483   const char       *name       = "Residual";
1484   DM                dmAux      = NULL;
1485   DM                dmGrad     = NULL;
1486   DMLabel           ghostLabel = NULL;
1487   PetscDS           prob       = NULL;
1488   PetscDS           probAux    = NULL;
1489   PetscSection      section    = NULL;
1490   PetscBool         useFEM     = PETSC_FALSE;
1491   PetscBool         useFVM     = PETSC_FALSE;
1492   PetscBool         isImplicit = (locX_t || time == PETSC_MIN_REAL) ? PETSC_TRUE : PETSC_FALSE;
1493   PetscFV           fvm        = NULL;
1494   PetscFECellGeom  *cgeomFEM   = NULL;
1495   PetscScalar      *cgeomScal;
1496   PetscFVCellGeom  *cgeomFVM   = NULL;
1497   PetscFVFaceGeom  *fgeomFVM   = NULL;
1498   Vec               locA, cellGeometryFEM = NULL, cellGeometryFVM = NULL, faceGeometryFVM = NULL, grad, locGrad = NULL;
1499   PetscScalar      *u, *u_t, *a, *uL, *uR;
1500   PetscInt          Nf, f, totDim, totDimAux, numChunks, cellChunkSize, faceChunkSize, chunk, fStart, fEnd;
1501   PetscErrorCode    ierr;
1502 
1503   PetscFunctionBegin;
1504   ierr = PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
1505   /* TODO The places where we have to use isFE are probably the member functions for the PetscDisc class */
1506   /* TODO The FVM geometry is over-manipulated. Make the precalc functions return exactly what we need */
1507   /* FEM+FVM */
1508   /* 1: Get sizes from dm and dmAux */
1509   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1510   ierr = DMGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr);
1511   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
1512   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
1513   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
1514   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &locA);CHKERRQ(ierr);
1515   if (locA) {
1516     ierr = VecGetDM(locA, &dmAux);CHKERRQ(ierr);
1517     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
1518     ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
1519   }
1520   /* 2: Get geometric data */
1521   for (f = 0; f < Nf; ++f) {
1522     PetscObject  obj;
1523     PetscClassId id;
1524     PetscBool    fimp;
1525 
1526     ierr = PetscDSGetImplicit(prob, f, &fimp);CHKERRQ(ierr);
1527     if (isImplicit != fimp) continue;
1528     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
1529     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1530     if (id == PETSCFE_CLASSID) {useFEM = PETSC_TRUE;}
1531     if (id == PETSCFV_CLASSID) {useFVM = PETSC_TRUE; fvm = (PetscFV) obj;}
1532   }
1533   if (useFEM) {
1534     ierr = DMPlexSNESGetGeometryFEM(dm, &cellGeometryFEM);CHKERRQ(ierr);
1535     ierr = VecGetArray(cellGeometryFEM, &cgeomScal);CHKERRQ(ierr);
1536     if (sizeof(PetscFECellGeom) % sizeof(PetscScalar)) {
1537       DM dmCell;
1538       PetscInt c;
1539 
1540       ierr = VecGetDM(cellGeometryFEM,&dmCell);CHKERRQ(ierr);
1541       ierr = PetscMalloc1(cEnd-cStart,&cgeomFEM);CHKERRQ(ierr);
1542       for (c = 0; c < cEnd - cStart; c++) {
1543         PetscScalar *thisgeom;
1544 
1545         ierr = DMPlexPointLocalRef(dmCell, c + cStart, cgeomScal, &thisgeom);CHKERRQ(ierr);
1546         cgeomFEM[c] = *((PetscFECellGeom *) thisgeom);
1547       }
1548     }
1549     else {
1550       cgeomFEM = (PetscFECellGeom *) cgeomScal;
1551     }
1552   }
1553   if (useFVM) {
1554     ierr = DMPlexSNESGetGeometryFVM(dm, &faceGeometryFVM, &cellGeometryFVM, NULL);CHKERRQ(ierr);
1555     ierr = VecGetArrayRead(faceGeometryFVM, (const PetscScalar **) &fgeomFVM);CHKERRQ(ierr);
1556     ierr = VecGetArrayRead(cellGeometryFVM, (const PetscScalar **) &cgeomFVM);CHKERRQ(ierr);
1557     /* Reconstruct and limit cell gradients */
1558     ierr = DMPlexSNESGetGradientDM(dm, fvm, &dmGrad);CHKERRQ(ierr);
1559     if (dmGrad) {
1560       ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
1561       ierr = DMGetGlobalVector(dmGrad, &grad);CHKERRQ(ierr);
1562       ierr = DMPlexReconstructGradients_Internal(dm, fvm, fStart, fEnd, faceGeometryFVM, cellGeometryFVM, locX, grad);CHKERRQ(ierr);
1563       /* Communicate gradient values */
1564       ierr = DMGetLocalVector(dmGrad, &locGrad);CHKERRQ(ierr);
1565       ierr = DMGlobalToLocalBegin(dmGrad, grad, INSERT_VALUES, locGrad);CHKERRQ(ierr);
1566       ierr = DMGlobalToLocalEnd(dmGrad, grad, INSERT_VALUES, locGrad);CHKERRQ(ierr);
1567       ierr = DMRestoreGlobalVector(dmGrad, &grad);CHKERRQ(ierr);
1568     }
1569     /* Handle non-essential (e.g. outflow) boundary values */
1570     ierr = DMPlexInsertBoundaryValues(dm, PETSC_FALSE, locX, time, faceGeometryFVM, cellGeometryFVM, locGrad);CHKERRQ(ierr);
1571   }
1572   /* Loop over chunks */
1573   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
1574   numChunks     = 1;
1575   cellChunkSize = (cEnd - cStart)/numChunks;
1576   faceChunkSize = (fEnd - fStart)/numChunks;
1577   for (chunk = 0; chunk < numChunks; ++chunk) {
1578     PetscScalar     *elemVec, *fluxL, *fluxR;
1579     PetscReal       *vol;
1580     PetscFVFaceGeom *fgeom;
1581     PetscInt         cS = cStart+chunk*cellChunkSize, cE = PetscMin(cS+cellChunkSize, cEnd), numCells = cE - cS, cell;
1582     PetscInt         fS = fStart+chunk*faceChunkSize, fE = PetscMin(fS+faceChunkSize, fEnd), numFaces = fE - fS, face;
1583 
1584     /* Extract field coefficients */
1585     if (useFEM) {
1586       ierr = DMPlexGetCellFields(dm, cS, cE, locX, locX_t, locA, &u, &u_t, &a);CHKERRQ(ierr);
1587       ierr = DMGetWorkArray(dm, numCells*totDim, PETSC_SCALAR, &elemVec);CHKERRQ(ierr);
1588       ierr = PetscMemzero(elemVec, numCells*totDim * sizeof(PetscScalar));CHKERRQ(ierr);
1589     }
1590     if (useFVM) {
1591       ierr = DMPlexGetFaceFields(dm, fS, fE, locX, locX_t, faceGeometryFVM, cellGeometryFVM, locGrad, &uL, &uR);CHKERRQ(ierr);
1592       ierr = DMPlexGetFaceGeometry(dm, fS, fE, faceGeometryFVM, cellGeometryFVM, &fgeom, &vol);CHKERRQ(ierr);
1593       ierr = DMGetWorkArray(dm, numFaces*totDim, PETSC_SCALAR, &fluxL);CHKERRQ(ierr);
1594       ierr = DMGetWorkArray(dm, numFaces*totDim, PETSC_SCALAR, &fluxR);CHKERRQ(ierr);
1595       ierr = PetscMemzero(fluxL, numFaces*totDim * sizeof(PetscScalar));CHKERRQ(ierr);
1596       ierr = PetscMemzero(fluxR, numFaces*totDim * sizeof(PetscScalar));CHKERRQ(ierr);
1597     }
1598     /* TODO We will interlace both our field coefficients (u, u_t, uL, uR, etc.) and our output (elemVec, fL, fR). I think this works */
1599     /* Loop over fields */
1600     for (f = 0; f < Nf; ++f) {
1601       PetscObject  obj;
1602       PetscClassId id;
1603       PetscBool    fimp;
1604       PetscInt     numChunks, numBatches, batchSize, numBlocks, blockSize, Ne, Nr, offset;
1605 
1606       ierr = PetscDSGetImplicit(prob, f, &fimp);CHKERRQ(ierr);
1607       if (isImplicit != fimp) continue;
1608       ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
1609       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1610       if (id == PETSCFE_CLASSID) {
1611         PetscFE         fe = (PetscFE) obj;
1612         PetscQuadrature q;
1613         PetscInt        Nq, Nb;
1614 
1615         ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
1616 
1617         ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
1618         ierr = PetscQuadratureGetData(q, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);
1619         ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
1620         blockSize = Nb*Nq;
1621         batchSize = numBlocks * blockSize;
1622         ierr =  PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1623         numChunks = numCells / (numBatches*batchSize);
1624         Ne        = numChunks*numBatches*batchSize;
1625         Nr        = numCells % (numBatches*batchSize);
1626         offset    = numCells - Nr;
1627         /* Integrate FE residual to get elemVec (need fields at quadrature points) */
1628         /*   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) */
1629         ierr = PetscFEIntegrateResidual(fe, prob, f, Ne, cgeomFEM, u, u_t, probAux, a, t, elemVec);CHKERRQ(ierr);
1630         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);
1631       } else if (id == PETSCFV_CLASSID) {
1632         PetscFV fv = (PetscFV) obj;
1633 
1634         Ne = numFaces;
1635         /* Riemann solve over faces (need fields at face centroids) */
1636         /*   We need to evaluate FE fields at those coordinates */
1637         ierr = PetscFVIntegrateRHSFunction(fv, prob, f, Ne, fgeom, vol, uL, uR, fluxL, fluxR);CHKERRQ(ierr);
1638       } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
1639     }
1640     if (sizeof(PetscFECellGeom) % sizeof(PetscScalar)) {
1641       ierr = PetscFree(cgeomFEM);CHKERRQ(ierr);
1642     }
1643     else {
1644       cgeomFEM = NULL;
1645     }
1646     if (cellGeometryFEM) {ierr = VecRestoreArray(cellGeometryFEM, &cgeomScal);CHKERRQ(ierr);}
1647     /* Loop over domain */
1648     if (useFEM) {
1649       /* Add elemVec to locX */
1650       for (cell = cS; cell < cE; ++cell) {
1651         if (mesh->printFEM > 1) {ierr = DMPrintCellVector(cell, name, totDim, &elemVec[cell*totDim]);CHKERRQ(ierr);}
1652         ierr = DMPlexVecSetClosure(dm, section, locF, cell, &elemVec[cell*totDim], ADD_ALL_VALUES);CHKERRQ(ierr);
1653       }
1654     }
1655     if (useFVM) {
1656       PetscScalar *fa;
1657       PetscInt     iface;
1658 
1659       ierr = VecGetArray(locF, &fa);CHKERRQ(ierr);
1660       for (f = 0; f < Nf; ++f) {
1661         PetscFV      fv;
1662         PetscObject  obj;
1663         PetscClassId id;
1664         PetscInt     foff, pdim;
1665 
1666         ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
1667         ierr = PetscDSGetFieldOffset(prob, f, &foff);CHKERRQ(ierr);
1668         ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1669         if (id != PETSCFV_CLASSID) continue;
1670         fv   = (PetscFV) obj;
1671         ierr = PetscFVGetNumComponents(fv, &pdim);CHKERRQ(ierr);
1672         /* Accumulate fluxes to cells */
1673         for (face = fS, iface = 0; face < fE; ++face) {
1674           const PetscInt *cells;
1675           PetscScalar    *fL, *fR;
1676           PetscInt        ghost, d, nsupp;
1677 
1678           ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr);
1679           ierr = DMPlexGetSupportSize(dm, face, &nsupp);CHKERRQ(ierr);
1680           if (ghost >= 0) continue;
1681           if (nsupp > 2) { /* noop */
1682             ++iface;
1683             continue;
1684           }
1685           ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr);
1686           ierr = DMPlexPointGlobalFieldRef(dm, cells[0], f, fa, &fL);CHKERRQ(ierr);
1687           ierr = DMPlexPointGlobalFieldRef(dm, cells[1], f, fa, &fR);CHKERRQ(ierr);
1688           for (d = 0; d < pdim; ++d) {
1689             if (fL) fL[d] -= fluxL[iface*totDim+foff+d];
1690             if (fR) fR[d] += fluxR[iface*totDim+foff+d];
1691           }
1692           ++iface;
1693         }
1694       }
1695       ierr = VecRestoreArray(locF, &fa);CHKERRQ(ierr);
1696     }
1697     /* Handle time derivative */
1698     if (locX_t) {
1699       PetscScalar *x_t, *fa;
1700 
1701       ierr = VecGetArray(locF, &fa);CHKERRQ(ierr);
1702       ierr = VecGetArray(locX_t, &x_t);CHKERRQ(ierr);
1703       for (f = 0; f < Nf; ++f) {
1704         PetscFV      fv;
1705         PetscObject  obj;
1706         PetscClassId id;
1707         PetscInt     pdim, d;
1708 
1709         ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
1710         ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1711         if (id != PETSCFV_CLASSID) continue;
1712         fv   = (PetscFV) obj;
1713         ierr = PetscFVGetNumComponents(fv, &pdim);CHKERRQ(ierr);
1714         for (cell = cS; cell < cE; ++cell) {
1715           PetscScalar *u_t, *r;
1716 
1717           ierr = DMPlexPointLocalFieldRead(dm, cell, f, x_t, &u_t);CHKERRQ(ierr);
1718           ierr = DMPlexPointLocalFieldRef(dm, cell, f, fa, &r);CHKERRQ(ierr);
1719           for (d = 0; d < pdim; ++d) r[d] += u_t[d];
1720         }
1721       }
1722       ierr = VecRestoreArray(locX_t, &x_t);CHKERRQ(ierr);
1723       ierr = VecRestoreArray(locF, &fa);CHKERRQ(ierr);
1724     }
1725     if (useFEM) {
1726       ierr = DMPlexRestoreCellFields(dm, cS, cE, locX, locX_t, locA, &u, &u_t, &a);CHKERRQ(ierr);
1727       ierr = DMRestoreWorkArray(dm, numCells*totDim, PETSC_SCALAR, &elemVec);CHKERRQ(ierr);
1728     }
1729     if (useFVM) {
1730       ierr = DMPlexRestoreFaceFields(dm, fS, fE, locX, locX_t, faceGeometryFVM, cellGeometryFVM, locGrad, &uL, &uR);CHKERRQ(ierr);
1731       ierr = DMPlexRestoreFaceGeometry(dm, fS, fE, faceGeometryFVM, cellGeometryFVM, &fgeom, &vol);CHKERRQ(ierr);
1732       ierr = DMRestoreWorkArray(dm, numFaces*totDim, PETSC_SCALAR, &fluxL);CHKERRQ(ierr);
1733       ierr = DMRestoreWorkArray(dm, numFaces*totDim, PETSC_SCALAR, &fluxR);CHKERRQ(ierr);
1734       if (dmGrad) {ierr = DMRestoreLocalVector(dmGrad, &locGrad);CHKERRQ(ierr);}
1735     }
1736   }
1737 
1738   if (useFEM) {ierr = DMPlexComputeBdResidual_Internal(dm, locX, locX_t, t, locF, user);CHKERRQ(ierr);}
1739 
1740   /* FEM */
1741   /* 1: Get sizes from dm and dmAux */
1742   /* 2: Get geometric data */
1743   /* 3: Handle boundary values */
1744   /* 4: Loop over domain */
1745   /*   Extract coefficients */
1746   /* Loop over fields */
1747   /*   Set tiling for FE*/
1748   /*   Integrate FE residual to get elemVec */
1749   /*     Loop over subdomain */
1750   /*       Loop over quad points */
1751   /*         Transform coords to real space */
1752   /*         Evaluate field and aux fields at point */
1753   /*         Evaluate residual at point */
1754   /*         Transform residual to real space */
1755   /*       Add residual to elemVec */
1756   /* Loop over domain */
1757   /*   Add elemVec to locX */
1758 
1759   /* FVM */
1760   /* Get geometric data */
1761   /* If using gradients */
1762   /*   Compute gradient data */
1763   /*   Loop over domain faces */
1764   /*     Count computational faces */
1765   /*     Reconstruct cell gradient */
1766   /*   Loop over domain cells */
1767   /*     Limit cell gradients */
1768   /* Handle boundary values */
1769   /* Loop over domain faces */
1770   /*   Read out field, centroid, normal, volume for each side of face */
1771   /* Riemann solve over faces */
1772   /* Loop over domain faces */
1773   /*   Accumulate fluxes to cells */
1774   /* TODO Change printFEM to printDisc here */
1775   if (mesh->printFEM) {
1776     Vec         locFbc;
1777     PetscInt    pStart, pEnd, p, maxDof;
1778     PetscScalar *zeroes;
1779 
1780     ierr = VecDuplicate(locF,&locFbc);CHKERRQ(ierr);
1781     ierr = VecCopy(locF,locFbc);CHKERRQ(ierr);
1782     ierr = PetscSectionGetChart(section,&pStart,&pEnd);CHKERRQ(ierr);
1783     ierr = PetscSectionGetMaxDof(section,&maxDof);CHKERRQ(ierr);
1784     ierr = PetscCalloc1(maxDof,&zeroes);CHKERRQ(ierr);
1785     for (p = pStart; p < pEnd; p++) {
1786       ierr = VecSetValuesSection(locFbc,section,p,zeroes,INSERT_BC_VALUES);CHKERRQ(ierr);
1787     }
1788     ierr = PetscFree(zeroes);CHKERRQ(ierr);
1789     ierr = DMPrintLocalVec(dm, name, mesh->printTol, locFbc);CHKERRQ(ierr);
1790     ierr = VecDestroy(&locFbc);CHKERRQ(ierr);
1791   }
1792   ierr = PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
1793   PetscFunctionReturn(0);
1794 }
1795 
1796 #undef __FUNCT__
1797 #define __FUNCT__ "DMPlexComputeResidualFEM_Check_Internal"
1798 static PetscErrorCode DMPlexComputeResidualFEM_Check_Internal(DM dm, Vec X, Vec X_t, PetscReal t, Vec F, void *user)
1799 {
1800   DM                dmCh, dmAux;
1801   Vec               A, cellgeom;
1802   PetscDS           prob, probCh, probAux = NULL;
1803   PetscQuadrature   q;
1804   PetscSection      section, sectionAux;
1805   PetscFECellGeom  *cgeom = NULL;
1806   PetscScalar      *cgeomScal;
1807   PetscScalar      *elemVec, *elemVecCh, *u, *u_t, *a = NULL;
1808   PetscInt          dim, Nf, f, numCells, cStart, cEnd, c;
1809   PetscInt          totDim, totDimAux, diffCell = 0;
1810   PetscErrorCode    ierr;
1811 
1812   PetscFunctionBegin;
1813   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1814   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1815   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
1816   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
1817   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
1818   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1819   numCells = cEnd - cStart;
1820   ierr = PetscObjectQuery((PetscObject) dm, "dmCh", (PetscObject *) &dmCh);CHKERRQ(ierr);
1821   ierr = DMGetDS(dmCh, &probCh);CHKERRQ(ierr);
1822   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
1823   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
1824   if (dmAux) {
1825     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
1826     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
1827     ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
1828   }
1829   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
1830   ierr = PetscMalloc3(numCells*totDim,&u,X_t ? numCells*totDim : 0,&u_t,numCells*totDim,&elemVec);CHKERRQ(ierr);
1831   ierr = PetscMalloc1(numCells*totDim,&elemVecCh);CHKERRQ(ierr);
1832   if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);}
1833   ierr = DMPlexSNESGetGeometryFEM(dm, &cellgeom);CHKERRQ(ierr);
1834   ierr = VecGetArray(cellgeom, &cgeomScal);CHKERRQ(ierr);
1835   if (sizeof(PetscFECellGeom) % sizeof(PetscScalar)) {
1836     DM dmCell;
1837 
1838     ierr = VecGetDM(cellgeom,&dmCell);CHKERRQ(ierr);
1839     ierr = PetscMalloc1(cEnd-cStart,&cgeom);CHKERRQ(ierr);
1840     for (c = 0; c < cEnd - cStart; c++) {
1841       PetscScalar *thisgeom;
1842 
1843       ierr = DMPlexPointLocalRef(dmCell, c + cStart, cgeomScal, &thisgeom);CHKERRQ(ierr);
1844       cgeom[c] = *((PetscFECellGeom *) thisgeom);
1845     }
1846   }
1847   else {
1848     cgeom = (PetscFECellGeom *) cgeomScal;
1849   }
1850   for (c = cStart; c < cEnd; ++c) {
1851     PetscScalar *x = NULL, *x_t = NULL;
1852     PetscInt     i;
1853 
1854     ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
1855     for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i];
1856     ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
1857     if (X_t) {
1858       ierr = DMPlexVecGetClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr);
1859       for (i = 0; i < totDim; ++i) u_t[c*totDim+i] = x_t[i];
1860       ierr = DMPlexVecRestoreClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr);
1861     }
1862     if (dmAux) {
1863       DM dmAuxPlex;
1864 
1865       ierr = DMSNESConvertPlex(dmAux,&dmAuxPlex, PETSC_FALSE);CHKERRQ(ierr);
1866       ierr = DMPlexVecGetClosure(dmAuxPlex, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
1867       for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i];
1868       ierr = DMPlexVecRestoreClosure(dmAuxPlex, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
1869       ierr = DMDestroy(&dmAuxPlex);CHKERRQ(ierr);
1870     }
1871   }
1872   for (f = 0; f < Nf; ++f) {
1873     PetscFE  fe, feCh;
1874     PetscInt numQuadPoints, Nb;
1875     /* Conforming batches */
1876     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1877     /* Remainder */
1878     PetscInt Nr, offset;
1879 
1880     ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
1881     ierr = PetscDSGetDiscretization(probCh, f, (PetscObject *) &feCh);CHKERRQ(ierr);
1882     ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
1883     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
1884     ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
1885     ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
1886     blockSize = Nb*numQuadPoints;
1887     batchSize = numBlocks * blockSize;
1888     ierr =  PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1889     numChunks = numCells / (numBatches*batchSize);
1890     Ne        = numChunks*numBatches*batchSize;
1891     Nr        = numCells % (numBatches*batchSize);
1892     offset    = numCells - Nr;
1893     ierr = PetscFEIntegrateResidual(fe, prob, f, Ne, cgeom, u, u_t, probAux, a, t, elemVec);CHKERRQ(ierr);
1894     ierr = PetscFEIntegrateResidual(feCh, prob, f, Ne, cgeom, u, u_t, probAux, a, t, elemVecCh);CHKERRQ(ierr);
1895     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);
1896     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);
1897   }
1898   for (c = cStart; c < cEnd; ++c) {
1899     PetscBool diff = PETSC_FALSE;
1900     PetscInt  d;
1901 
1902     for (d = 0; d < totDim; ++d) if (PetscAbsScalar(elemVec[c*totDim+d] - elemVecCh[c*totDim+d]) > 1.0e-7) {diff = PETSC_TRUE;break;}
1903     if (diff) {
1904       ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Different cell %d\n", c);CHKERRQ(ierr);
1905       ierr = DMPrintCellVector(c, "Residual", totDim, &elemVec[c*totDim]);CHKERRQ(ierr);
1906       ierr = DMPrintCellVector(c, "Check Residual", totDim, &elemVecCh[c*totDim]);CHKERRQ(ierr);
1907       ++diffCell;
1908     }
1909     if (diffCell > 9) break;
1910     ierr = DMPlexVecSetClosure(dm, section, F, c, &elemVec[c*totDim], ADD_ALL_VALUES);CHKERRQ(ierr);
1911   }
1912   if (sizeof(PetscFECellGeom) % sizeof(PetscScalar)) {
1913     ierr = PetscFree(cgeom);CHKERRQ(ierr);
1914   }
1915   else {
1916     cgeom = NULL;
1917   }
1918   ierr = VecRestoreArray(cellgeom, &cgeomScal);CHKERRQ(ierr);
1919   ierr = PetscFree3(u,u_t,elemVec);CHKERRQ(ierr);
1920   ierr = PetscFree(elemVecCh);CHKERRQ(ierr);
1921   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
1922   PetscFunctionReturn(0);
1923 }
1924 
1925 #undef __FUNCT__
1926 #define __FUNCT__ "DMPlexSNESComputeResidualFEM"
1927 /*@
1928   DMPlexSNESComputeResidualFEM - Form the local residual F from the local input X using pointwise functions specified by the user
1929 
1930   Input Parameters:
1931 + dm - The mesh
1932 . X  - Local solution
1933 - user - The user context
1934 
1935   Output Parameter:
1936 . F  - Local output vector
1937 
1938   Level: developer
1939 
1940 .seealso: DMPlexComputeJacobianActionFEM()
1941 @*/
1942 PetscErrorCode DMPlexSNESComputeResidualFEM(DM dm, Vec X, Vec F, void *user)
1943 {
1944   PetscObject    check;
1945   PetscInt       cStart, cEnd, cEndInterior;
1946   DM             plex;
1947   PetscErrorCode ierr;
1948 
1949   PetscFunctionBegin;
1950   ierr = DMSNESConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
1951   ierr = DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd);CHKERRQ(ierr);
1952   ierr = DMPlexGetHybridBounds(plex, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
1953   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
1954   /* The dmCh is used to check two mathematically equivalent discretizations for computational equivalence */
1955   ierr = PetscObjectQuery((PetscObject) plex, "dmCh", &check);CHKERRQ(ierr);
1956   if (check) {ierr = DMPlexComputeResidualFEM_Check_Internal(plex, X, NULL, 0.0, F, user);CHKERRQ(ierr);}
1957   else       {ierr = DMPlexComputeResidual_Internal(plex, cStart, cEnd, PETSC_MIN_REAL, X, NULL, 0.0, F, user);CHKERRQ(ierr);}
1958   ierr = DMDestroy(&plex);CHKERRQ(ierr);
1959   PetscFunctionReturn(0);
1960 }
1961 
1962 #undef __FUNCT__
1963 #define __FUNCT__ "DMPlexSNESComputeBoundaryFEM"
1964 /*@
1965   DMPlexSNESComputeBoundaryFEM - Form the boundary values for the local input X
1966 
1967   Input Parameters:
1968 + dm - The mesh
1969 - user - The user context
1970 
1971   Output Parameter:
1972 . X  - Local solution
1973 
1974   Level: developer
1975 
1976 .seealso: DMPlexComputeJacobianActionFEM()
1977 @*/
1978 PetscErrorCode DMPlexSNESComputeBoundaryFEM(DM dm, Vec X, void *user)
1979 {
1980   DM             plex;
1981   PetscErrorCode ierr;
1982 
1983   PetscFunctionBegin;
1984   ierr = DMSNESConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
1985   ierr = DMPlexInsertBoundaryValues(plex, PETSC_TRUE, X, PETSC_MIN_REAL, NULL, NULL, NULL);CHKERRQ(ierr);
1986   ierr = DMDestroy(&plex);CHKERRQ(ierr);
1987   PetscFunctionReturn(0);
1988 }
1989 
1990 #undef __FUNCT__
1991 #define __FUNCT__ "DMPlexComputeJacobian_Internal"
1992 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)
1993 {
1994   DM_Plex          *mesh  = (DM_Plex *) dm->data;
1995   const char       *name  = "Jacobian";
1996   DM                dmAux, plex;
1997   DMLabel           depth;
1998   Vec               A, cellgeom;
1999   PetscDS           prob, probAux = NULL;
2000   PetscQuadrature   quad;
2001   PetscSection      section, globalSection, sectionAux;
2002   PetscFECellGeom  *cgeom = NULL;
2003   PetscScalar      *cgeomScal;
2004   PetscScalar      *elemMat, *elemMatP, *elemMatD, *u, *u_t, *a = NULL;
2005   PetscInt          dim, Nf, f, fieldI, fieldJ, numCells, c;
2006   PetscInt          totDim, totDimBd, totDimAux, numBd, bd;
2007   PetscBool         isShell, hasJac, hasPrec, hasDyn, hasFV = PETSC_FALSE;
2008   PetscErrorCode    ierr;
2009 
2010   PetscFunctionBegin;
2011   ierr = PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
2012   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
2013   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
2014   ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);
2015   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
2016   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
2017   ierr = PetscDSGetTotalBdDimension(prob, &totDimBd);CHKERRQ(ierr);
2018   ierr = PetscDSHasJacobian(prob, &hasJac);CHKERRQ(ierr);
2019   ierr = PetscDSHasJacobianPreconditioner(prob, &hasPrec);CHKERRQ(ierr);
2020   ierr = PetscDSHasDynamicJacobian(prob, &hasDyn);CHKERRQ(ierr);
2021   hasDyn = hasDyn && (X_tShift != 0.0) ? PETSC_TRUE : PETSC_FALSE;
2022   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
2023   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
2024   numCells = cEnd - cStart;
2025   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
2026   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
2027   if (dmAux) {
2028     ierr = DMConvert(dmAux, DMPLEX, &plex);CHKERRQ(ierr);
2029     ierr = DMGetDefaultSection(plex, &sectionAux);CHKERRQ(ierr);
2030     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
2031     ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
2032   }
2033   ierr = MatZeroEntries(JacP);CHKERRQ(ierr);
2034   ierr = PetscMalloc5(numCells*totDim,&u,X_t ? numCells*totDim : 0,&u_t,hasJac ? numCells*totDim*totDim : 0,&elemMat,hasPrec ? numCells*totDim*totDim : 0, &elemMatP,hasDyn ? numCells*totDim*totDim : 0, &elemMatD);CHKERRQ(ierr);
2035   if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);}
2036   ierr = DMPlexSNESGetGeometryFEM(dm, &cellgeom);CHKERRQ(ierr);
2037   ierr = VecGetArray(cellgeom, &cgeomScal);CHKERRQ(ierr);
2038   if (sizeof(PetscFECellGeom) % sizeof(PetscScalar)) {
2039     DM dmCell;
2040 
2041     ierr = VecGetDM(cellgeom,&dmCell);CHKERRQ(ierr);
2042     ierr = PetscMalloc1(cEnd-cStart,&cgeom);CHKERRQ(ierr);
2043     for (c = 0; c < cEnd - cStart; c++) {
2044       PetscScalar *thisgeom;
2045 
2046       ierr = DMPlexPointLocalRef(dmCell, c + cStart, cgeomScal, &thisgeom);CHKERRQ(ierr);
2047       cgeom[c] = *((PetscFECellGeom *) thisgeom);
2048     }
2049   } else {
2050     cgeom = (PetscFECellGeom *) cgeomScal;
2051   }
2052   for (c = cStart; c < cEnd; ++c) {
2053     PetscScalar *x = NULL,  *x_t = NULL;
2054     PetscInt     i;
2055 
2056     ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
2057     for (i = 0; i < totDim; ++i) u[(c-cStart)*totDim+i] = x[i];
2058     ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
2059     if (X_t) {
2060       ierr = DMPlexVecGetClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr);
2061       for (i = 0; i < totDim; ++i) u_t[(c-cStart)*totDim+i] = x_t[i];
2062       ierr = DMPlexVecRestoreClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr);
2063     }
2064     if (dmAux) {
2065       ierr = DMPlexVecGetClosure(plex, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
2066       for (i = 0; i < totDimAux; ++i) a[(c-cStart)*totDimAux+i] = x[i];
2067       ierr = DMPlexVecRestoreClosure(plex, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
2068     }
2069   }
2070   if (hasJac)  {ierr = PetscMemzero(elemMat,  numCells*totDim*totDim * sizeof(PetscScalar));CHKERRQ(ierr);}
2071   if (hasPrec) {ierr = PetscMemzero(elemMatP, numCells*totDim*totDim * sizeof(PetscScalar));CHKERRQ(ierr);}
2072   if (hasDyn)  {ierr = PetscMemzero(elemMatD, numCells*totDim*totDim * sizeof(PetscScalar));CHKERRQ(ierr);}
2073   for (fieldI = 0; fieldI < Nf; ++fieldI) {
2074     PetscClassId id;
2075     PetscFE  fe;
2076     PetscInt numQuadPoints, Nb;
2077     /* Conforming batches */
2078     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
2079     /* Remainder */
2080     PetscInt Nr, offset;
2081 
2082     ierr = PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);CHKERRQ(ierr);
2083     ierr = PetscObjectGetClassId((PetscObject) fe, &id);CHKERRQ(ierr);
2084     if (id == PETSCFV_CLASSID) {hasFV = PETSC_TRUE; continue;}
2085     ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
2086     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
2087     ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
2088     ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
2089     blockSize = Nb*numQuadPoints;
2090     batchSize = numBlocks * blockSize;
2091     ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
2092     numChunks = numCells / (numBatches*batchSize);
2093     Ne        = numChunks*numBatches*batchSize;
2094     Nr        = numCells % (numBatches*batchSize);
2095     offset    = numCells - Nr;
2096     for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
2097       if (hasJac) {
2098         ierr = PetscFEIntegrateJacobian(fe, prob, PETSCFE_JACOBIAN, fieldI, fieldJ, Ne, cgeom, u, u_t, probAux, a, t, X_tShift, elemMat);CHKERRQ(ierr);
2099         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);
2100       }
2101       if (hasPrec) {
2102         ierr = PetscFEIntegrateJacobian(fe, prob, PETSCFE_JACOBIAN_PRE, fieldI, fieldJ, Ne, cgeom, u, u_t, probAux, a, t, X_tShift, elemMatP);CHKERRQ(ierr);
2103         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);
2104       }
2105       if (hasDyn) {
2106         ierr = PetscFEIntegrateJacobian(fe, prob, PETSCFE_JACOBIAN_DYN, fieldI, fieldJ, Ne, cgeom, u, u_t, probAux, a, t, X_tShift, elemMatD);CHKERRQ(ierr);
2107         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);
2108       }
2109     }
2110   }
2111   if (hasDyn) {
2112     for (c = 0; c < (cEnd - cStart)*totDim*totDim; ++c) elemMat[c] += X_tShift*elemMatD[c];
2113   }
2114   if (hasFV) {
2115     PetscClassId id;
2116     PetscFV      fv;
2117     PetscInt     offsetI, NcI, NbI = 1, fc, f;
2118 
2119     for (fieldI = 0; fieldI < Nf; ++fieldI) {
2120       ierr = PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fv);CHKERRQ(ierr);
2121       ierr = PetscDSGetFieldOffset(prob, fieldI, &offsetI);CHKERRQ(ierr);
2122       ierr = PetscObjectGetClassId((PetscObject) fv, &id);CHKERRQ(ierr);
2123       if (id != PETSCFV_CLASSID) continue;
2124       /* Put in the identity */
2125       ierr = PetscFVGetNumComponents(fv, &NcI);CHKERRQ(ierr);
2126       for (c = cStart; c < cEnd; ++c) {
2127         const PetscInt eOffset = (c-cStart)*totDim*totDim;
2128         for (fc = 0; fc < NcI; ++fc) {
2129           for (f = 0; f < NbI; ++f) {
2130             const PetscInt i = offsetI + f*NcI+fc;
2131             if (hasPrec) {
2132               if (hasJac) {elemMat[eOffset+i*totDim+i] = 1.0;}
2133               elemMatP[eOffset+i*totDim+i] = 1.0;
2134             } else {elemMat[eOffset+i*totDim+i] = 1.0;}
2135           }
2136         }
2137       }
2138     }
2139     /* No allocated space for FV stuff, so ignore the zero entries */
2140     ierr = MatSetOption(JacP, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE);CHKERRQ(ierr);
2141   }
2142   for (c = cStart; c < cEnd; ++c) {
2143     if (hasPrec) {
2144       if (hasJac) {
2145         if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, name, totDim, totDim, &elemMat[(c-cStart)*totDim*totDim]);CHKERRQ(ierr);}
2146         ierr = DMPlexMatSetClosure(dm, section, globalSection, Jac, c, &elemMat[(c-cStart)*totDim*totDim], ADD_VALUES);CHKERRQ(ierr);
2147       }
2148       if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, name, totDim, totDim, &elemMatP[(c-cStart)*totDim*totDim]);CHKERRQ(ierr);}
2149       ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, c, &elemMatP[(c-cStart)*totDim*totDim], ADD_VALUES);CHKERRQ(ierr);
2150     } else {
2151       if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, name, totDim, totDim, &elemMat[(c-cStart)*totDim*totDim]);CHKERRQ(ierr);}
2152       ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, c, &elemMat[(c-cStart)*totDim*totDim], ADD_VALUES);CHKERRQ(ierr);
2153     }
2154   }
2155   if (hasFV) {ierr = MatSetOption(JacP, MAT_IGNORE_ZERO_ENTRIES, PETSC_FALSE);CHKERRQ(ierr);}
2156   if (sizeof(PetscFECellGeom) % sizeof(PetscScalar)) {
2157     ierr = PetscFree(cgeom);CHKERRQ(ierr);
2158   } else {
2159     cgeom = NULL;
2160   }
2161   ierr = VecRestoreArray(cellgeom, &cgeomScal);CHKERRQ(ierr);
2162   ierr = PetscFree5(u,u_t,elemMat,elemMatP,elemMatD);CHKERRQ(ierr);
2163   if (dmAux) {
2164     ierr = PetscFree(a);CHKERRQ(ierr);
2165     ierr = DMDestroy(&plex);CHKERRQ(ierr);
2166   }
2167   ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr);
2168   ierr = DMGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
2169   ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr);
2170   ierr = DMGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
2171   for (bd = 0; bd < numBd; ++bd) {
2172     const char     *bdLabel;
2173     DMLabel         label;
2174     IS              pointIS;
2175     const PetscInt *points;
2176     const PetscInt *values;
2177     PetscInt        field, numValues, v, numPoints, p, dep, numFaces;
2178     PetscBool       isEssential;
2179     PetscObject     obj;
2180     PetscClassId    id;
2181 
2182     ierr = DMGetBoundary(dm, bd, &isEssential, NULL, &bdLabel, &field, NULL, NULL, NULL, &numValues, &values, NULL);CHKERRQ(ierr);
2183     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
2184     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
2185     if ((id != PETSCFE_CLASSID) || isEssential) continue;
2186     ierr = DMGetLabel(dm, bdLabel, &label);CHKERRQ(ierr);
2187     for (v = 0; v < numValues; ++v) {
2188       ierr = DMLabelGetStratumSize(label, values[v], &numPoints);CHKERRQ(ierr);
2189       ierr = DMLabelGetStratumIS(label, values[v], &pointIS);CHKERRQ(ierr);
2190       if (!pointIS) continue; /* No points with that id on this process */
2191       ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
2192       for (p = 0, numFaces = 0; p < numPoints; ++p) {
2193         ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
2194         if (dep == dim-1) ++numFaces;
2195       }
2196       ierr = PetscMalloc3(numFaces*totDimBd,&u,numFaces,&cgeom,numFaces*totDimBd*totDimBd,&elemMat);CHKERRQ(ierr);
2197       if (X_t) {ierr = PetscMalloc1(numFaces*totDimBd,&u_t);CHKERRQ(ierr);}
2198       for (p = 0, f = 0; p < numPoints; ++p) {
2199         const PetscInt point = points[p];
2200         PetscScalar   *x     = NULL;
2201         PetscInt       i;
2202 
2203         ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
2204         if (dep != dim-1) continue;
2205         ierr = DMPlexComputeCellGeometryFEM(dm, point, NULL, cgeom[f].v0, cgeom[f].J, cgeom[f].invJ, &cgeom[f].detJ);CHKERRQ(ierr);
2206         ierr = DMPlexComputeCellGeometryFVM(dm, point, NULL, NULL, cgeom[f].n);CHKERRQ(ierr);
2207         if (cgeom[f].detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", cgeom[f].detJ, point);
2208         ierr = DMPlexVecGetClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
2209         for (i = 0; i < totDimBd; ++i) u[f*totDimBd+i] = x[i];
2210         ierr = DMPlexVecRestoreClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
2211         if (X_t) {
2212           ierr = DMPlexVecGetClosure(dm, section, X_t, point, NULL, &x);CHKERRQ(ierr);
2213           for (i = 0; i < totDimBd; ++i) u_t[f*totDimBd+i] = x[i];
2214           ierr = DMPlexVecRestoreClosure(dm, section, X_t, point, NULL, &x);CHKERRQ(ierr);
2215         }
2216         ++f;
2217       }
2218       ierr = PetscMemzero(elemMat, numFaces*totDimBd*totDimBd * sizeof(PetscScalar));CHKERRQ(ierr);
2219       for (fieldI = 0; fieldI < Nf; ++fieldI) {
2220         PetscFE  fe;
2221         PetscInt numQuadPoints, Nb;
2222         /* Conforming batches */
2223         PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
2224         /* Remainder */
2225         PetscInt Nr, offset;
2226 
2227         ierr = PetscDSGetBdDiscretization(prob, fieldI, (PetscObject *) &fe);CHKERRQ(ierr);
2228         ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
2229         ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
2230         ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
2231         ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
2232         blockSize = Nb*numQuadPoints;
2233         batchSize = numBlocks * blockSize;
2234         ierr =  PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
2235         numChunks = numFaces / (numBatches*batchSize);
2236         Ne        = numChunks*numBatches*batchSize;
2237         Nr        = numFaces % (numBatches*batchSize);
2238         offset    = numFaces - Nr;
2239         for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
2240           ierr = PetscFEIntegrateBdJacobian(fe, prob, fieldI, fieldJ, Ne, cgeom, u, u_t, NULL, NULL, t, X_tShift, elemMat);CHKERRQ(ierr);
2241           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);
2242         }
2243       }
2244       for (p = 0, f = 0; p < numPoints; ++p) {
2245         const PetscInt point = points[p];
2246 
2247         ierr = DMLabelGetValue(depth, point, &dep);CHKERRQ(ierr);
2248         if (dep != dim-1) continue;
2249         if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(point, "BdJacobian", totDimBd, totDimBd, &elemMat[f*totDimBd*totDimBd]);CHKERRQ(ierr);}
2250         ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, point, &elemMat[f*totDimBd*totDimBd], ADD_VALUES);CHKERRQ(ierr);
2251         ++f;
2252       }
2253       ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
2254       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
2255       ierr = PetscFree3(u,cgeom,elemMat);CHKERRQ(ierr);
2256       if (X_t) {ierr = PetscFree(u_t);CHKERRQ(ierr);}
2257     }
2258   }
2259   if (hasJac && hasPrec) {
2260     ierr = MatAssemblyBegin(Jac, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2261     ierr = MatAssemblyEnd(Jac, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2262   }
2263   ierr = MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2264   ierr = MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2265   if (mesh->printFEM) {
2266     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
2267     ierr = MatChop(JacP, 1.0e-10);CHKERRQ(ierr);
2268     ierr = MatView(JacP, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
2269   }
2270   ierr = PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
2271   ierr = PetscObjectTypeCompare((PetscObject) Jac, MATSHELL, &isShell);CHKERRQ(ierr);
2272   if (isShell) {
2273     JacActionCtx *jctx;
2274 
2275     ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr);
2276     ierr = VecCopy(X, jctx->u);CHKERRQ(ierr);
2277   }
2278   PetscFunctionReturn(0);
2279 }
2280 
2281 
2282 #undef __FUNCT__
2283 #define __FUNCT__ "DMPlexComputeJacobianAction_Internal"
2284 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)
2285 {
2286   DM_Plex          *mesh  = (DM_Plex *) dm->data;
2287   const char       *name  = "Jacobian";
2288   DM                dmAux, plex;
2289   Vec               A, cellgeom;
2290   PetscDS           prob, probAux = NULL;
2291   PetscQuadrature   quad;
2292   PetscSection      section, globalSection, sectionAux;
2293   PetscFECellGeom  *cgeom = NULL;
2294   PetscScalar      *cgeomScal;
2295   PetscScalar      *elemMat, *elemMatD, *u, *u_t, *a = NULL, *y, *z;
2296   PetscInt          dim, Nf, fieldI, fieldJ, numCells, c;
2297   PetscInt          totDim, totDimBd, totDimAux;
2298   PetscBool         hasDyn;
2299   PetscErrorCode    ierr;
2300 
2301   PetscFunctionBegin;
2302   ierr = PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
2303   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
2304   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
2305   ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);
2306   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
2307   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
2308   ierr = PetscDSGetTotalBdDimension(prob, &totDimBd);CHKERRQ(ierr);
2309   ierr = PetscDSHasDynamicJacobian(prob, &hasDyn);CHKERRQ(ierr);
2310   hasDyn = hasDyn && (X_tShift != 0.0) ? PETSC_TRUE : PETSC_FALSE;
2311   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
2312   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
2313   numCells = cEnd - cStart;
2314   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
2315   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
2316   if (dmAux) {
2317     ierr = DMConvert(dmAux, DMPLEX, &plex);CHKERRQ(ierr);
2318     ierr = DMGetDefaultSection(plex, &sectionAux);CHKERRQ(ierr);
2319     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
2320     ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
2321   }
2322   ierr = VecSet(Z, 0.0);CHKERRQ(ierr);
2323   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);
2324   if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);}
2325   ierr = DMPlexSNESGetGeometryFEM(dm, &cellgeom);CHKERRQ(ierr);
2326   ierr = VecGetArray(cellgeom, &cgeomScal);CHKERRQ(ierr);
2327   if (sizeof(PetscFECellGeom) % sizeof(PetscScalar)) {
2328     DM dmCell;
2329 
2330     ierr = VecGetDM(cellgeom,&dmCell);CHKERRQ(ierr);
2331     ierr = PetscMalloc1(cEnd-cStart,&cgeom);CHKERRQ(ierr);
2332     for (c = 0; c < cEnd - cStart; c++) {
2333       PetscScalar *thisgeom;
2334 
2335       ierr = DMPlexPointLocalRef(dmCell, c + cStart, cgeomScal, &thisgeom);CHKERRQ(ierr);
2336       cgeom[c] = *((PetscFECellGeom *) thisgeom);
2337     }
2338   } else {
2339     cgeom = (PetscFECellGeom *) cgeomScal;
2340   }
2341   for (c = cStart; c < cEnd; ++c) {
2342     PetscScalar *x = NULL,  *x_t = NULL;
2343     PetscInt     i;
2344 
2345     ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
2346     for (i = 0; i < totDim; ++i) u[(c-cStart)*totDim+i] = x[i];
2347     ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
2348     if (X_t) {
2349       ierr = DMPlexVecGetClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr);
2350       for (i = 0; i < totDim; ++i) u_t[(c-cStart)*totDim+i] = x_t[i];
2351       ierr = DMPlexVecRestoreClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr);
2352     }
2353     if (dmAux) {
2354       ierr = DMPlexVecGetClosure(plex, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
2355       for (i = 0; i < totDimAux; ++i) a[(c-cStart)*totDimAux+i] = x[i];
2356       ierr = DMPlexVecRestoreClosure(plex, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
2357     }
2358     ierr = DMPlexVecGetClosure(dm, section, Y, c, NULL, &x);CHKERRQ(ierr);
2359     for (i = 0; i < totDim; ++i) y[(c-cStart)*totDim+i] = x[i];
2360     ierr = DMPlexVecRestoreClosure(dm, section, Y, c, NULL, &x);CHKERRQ(ierr);
2361   }
2362   ierr = PetscMemzero(elemMat, numCells*totDim*totDim * sizeof(PetscScalar));CHKERRQ(ierr);
2363   if (hasDyn)  {ierr = PetscMemzero(elemMatD, numCells*totDim*totDim * sizeof(PetscScalar));CHKERRQ(ierr);}
2364   for (fieldI = 0; fieldI < Nf; ++fieldI) {
2365     PetscFE  fe;
2366     PetscInt numQuadPoints, Nb;
2367     /* Conforming batches */
2368     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
2369     /* Remainder */
2370     PetscInt Nr, offset;
2371 
2372     ierr = PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);CHKERRQ(ierr);
2373     ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
2374     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
2375     ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
2376     ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
2377     blockSize = Nb*numQuadPoints;
2378     batchSize = numBlocks * blockSize;
2379     ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
2380     numChunks = numCells / (numBatches*batchSize);
2381     Ne        = numChunks*numBatches*batchSize;
2382     Nr        = numCells % (numBatches*batchSize);
2383     offset    = numCells - Nr;
2384     for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
2385       ierr = PetscFEIntegrateJacobian(fe, prob, PETSCFE_JACOBIAN, fieldI, fieldJ, Ne, cgeom, u, u_t, probAux, a, t, X_tShift, elemMat);CHKERRQ(ierr);
2386       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);
2387       if (hasDyn) {
2388         ierr = PetscFEIntegrateJacobian(fe, prob, PETSCFE_JACOBIAN_DYN, fieldI, fieldJ, Ne, cgeom, u, u_t, probAux, a, t, X_tShift, elemMatD);CHKERRQ(ierr);
2389         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);
2390       }
2391     }
2392   }
2393   if (hasDyn) {
2394     for (c = 0; c < (cEnd - cStart)*totDim*totDim; ++c) elemMat[c] += X_tShift*elemMatD[c];
2395   }
2396   for (c = cStart; c < cEnd; ++c) {
2397     const PetscBLASInt M = totDim, one = 1;
2398     const PetscScalar  a = 1.0, b = 0.0;
2399 
2400     PetscStackCallBLAS("BLASgemv", BLASgemv_("N", &M, &M, &a, &elemMat[(c-cStart)*totDim*totDim], &M, &y[(c-cStart)*totDim], &one, &b, z, &one));
2401     if (mesh->printFEM > 1) {
2402       ierr = DMPrintCellMatrix(c, name, totDim, totDim, &elemMat[(c-cStart)*totDim*totDim]);CHKERRQ(ierr);
2403       ierr = DMPrintCellVector(c, "Y",  totDim, &y[(c-cStart)*totDim]);CHKERRQ(ierr);
2404       ierr = DMPrintCellVector(c, "Z",  totDim, z);CHKERRQ(ierr);
2405     }
2406     ierr = DMPlexVecSetClosure(dm, section, Z, c, z, ADD_VALUES);CHKERRQ(ierr);
2407   }
2408   if (sizeof(PetscFECellGeom) % sizeof(PetscScalar)) {ierr = PetscFree(cgeom);CHKERRQ(ierr);}
2409   else                                               {cgeom = NULL;}
2410   ierr = VecRestoreArray(cellgeom, &cgeomScal);CHKERRQ(ierr);
2411   ierr = PetscFree6(u,u_t,elemMat,elemMatD,y,z);CHKERRQ(ierr);
2412   if (dmAux) {
2413     ierr = PetscFree(a);CHKERRQ(ierr);
2414     ierr = DMDestroy(&plex);CHKERRQ(ierr);
2415   }
2416   if (mesh->printFEM) {
2417     ierr = PetscPrintf(PETSC_COMM_WORLD, "Z:\n");CHKERRQ(ierr);
2418     ierr = VecView(Z, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
2419   }
2420   ierr = PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
2421   PetscFunctionReturn(0);
2422 }
2423 
2424 #undef __FUNCT__
2425 #define __FUNCT__ "DMPlexSNESComputeJacobianFEM"
2426 /*@
2427   DMPlexSNESComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user.
2428 
2429   Input Parameters:
2430 + dm - The mesh
2431 . X  - Local input vector
2432 - user - The user context
2433 
2434   Output Parameter:
2435 . Jac  - Jacobian matrix
2436 
2437   Note:
2438   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
2439   like a GPU, or vectorize on a multicore machine.
2440 
2441   Level: developer
2442 
2443 .seealso: FormFunctionLocal()
2444 @*/
2445 PetscErrorCode DMPlexSNESComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP,void *user)
2446 {
2447   PetscInt       cStart, cEnd, cEndInterior;
2448   DM             plex;
2449   PetscErrorCode ierr;
2450 
2451   PetscFunctionBegin;
2452   ierr = DMSNESConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
2453   ierr = DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd);CHKERRQ(ierr);
2454   ierr = DMPlexGetHybridBounds(plex, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
2455   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
2456   ierr = DMPlexComputeJacobian_Internal(plex, cStart, cEnd, 0.0, 0.0, X, NULL, Jac, JacP, user);CHKERRQ(ierr);
2457   ierr = DMDestroy(&plex);CHKERRQ(ierr);
2458   PetscFunctionReturn(0);
2459 }
2460 
2461 #undef __FUNCT__
2462 #define __FUNCT__ "DMPlexSNESComputeJacobianActionFEM"
2463 /*@
2464   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.
2465 
2466   Input Parameters:
2467 + dm - The mesh
2468 . X  - Local solution vector
2469 . Y  - Local input vector
2470 - user - The user context
2471 
2472   Output Parameter:
2473 . Z - Local output vector
2474 
2475   Note:
2476   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
2477   like a GPU, or vectorize on a multicore machine.
2478 
2479   Level: developer
2480 
2481 .seealso: FormFunctionLocal()
2482 @*/
2483 PetscErrorCode DMPlexSNESComputeJacobianActionFEM(DM dm, Vec X, Vec Y, Vec Z, void *user)
2484 {
2485   PetscInt       cStart, cEnd, cEndInterior;
2486   DM             plex;
2487   PetscErrorCode ierr;
2488 
2489   PetscFunctionBegin;
2490   ierr = DMSNESConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
2491   ierr = DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd);CHKERRQ(ierr);
2492   ierr = DMPlexGetHybridBounds(plex, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
2493   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
2494   ierr = DMPlexComputeJacobianAction_Internal(plex, cStart, cEnd, 0.0, 0.0, X, NULL, Y, Z, user);CHKERRQ(ierr);
2495   ierr = DMDestroy(&plex);CHKERRQ(ierr);
2496   PetscFunctionReturn(0);
2497 }
2498 
2499 #undef __FUNCT__
2500 #define __FUNCT__ "DMPlexSetSNESLocalFEM"
2501 /*@
2502   DMPlexSetSNESLocalFEM - Use DMPlex's internal FEM routines to compute SNES boundary values, residual, and Jacobian.
2503 
2504   Input Parameters:
2505 + dm - The DM object
2506 . boundaryctx - the user context that will be passed to pointwise evaluation of boundary values (see DMAddBoundary())
2507 . residualctx - the user context that will be passed to pointwise evaluation of finite element residual computations (see PetscDSSetResidual())
2508 - jacobianctx - the user context that will be passed to pointwise evaluation of finite element Jacobian construction (see PetscDSSetJacobian())
2509 
2510   Level: developer
2511 @*/
2512 PetscErrorCode DMPlexSetSNESLocalFEM(DM dm, void *boundaryctx, void *residualctx, void *jacobianctx)
2513 {
2514   PetscErrorCode ierr;
2515 
2516   PetscFunctionBegin;
2517   ierr = DMSNESSetBoundaryLocal(dm,DMPlexSNESComputeBoundaryFEM,boundaryctx);CHKERRQ(ierr);
2518   ierr = DMSNESSetFunctionLocal(dm,DMPlexSNESComputeResidualFEM,residualctx);CHKERRQ(ierr);
2519   ierr = DMSNESSetJacobianLocal(dm,DMPlexSNESComputeJacobianFEM,jacobianctx);CHKERRQ(ierr);
2520   PetscFunctionReturn(0);
2521 }
2522 
2523 #undef __FUNCT__
2524 #define __FUNCT__ "DMSNESCheckFromOptions_Internal"
2525 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)
2526 {
2527   Mat            J, M;
2528   Vec            r, b;
2529   MatNullSpace   nullSpace;
2530   PetscReal     *error, res = 0.0;
2531   PetscInt       numFields;
2532   PetscErrorCode ierr;
2533 
2534   PetscFunctionBegin;
2535   ierr = VecDuplicate(u, &r);CHKERRQ(ierr);
2536   ierr = DMCreateMatrix(dm, &J);CHKERRQ(ierr);
2537   M    = J;
2538   /* TODO Null space for J */
2539   /* Check discretization error */
2540   ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr);
2541   ierr = PetscMalloc1(PetscMax(1, numFields), &error);CHKERRQ(ierr);
2542   if (numFields > 1) {
2543     PetscInt f;
2544 
2545     ierr = DMComputeL2FieldDiff(dm, 0.0, exactFuncs, ctxs, u, error);CHKERRQ(ierr);
2546     ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: [");CHKERRQ(ierr);
2547     for (f = 0; f < numFields; ++f) {
2548       if (f) {ierr = PetscPrintf(PETSC_COMM_WORLD, ", ");CHKERRQ(ierr);}
2549       if (error[f] >= 1.0e-11) {ierr = PetscPrintf(PETSC_COMM_WORLD, "%g", error[f]);CHKERRQ(ierr);}
2550       else                     {ierr = PetscPrintf(PETSC_COMM_WORLD, "< 1.0e-11");CHKERRQ(ierr);}
2551     }
2552     ierr = PetscPrintf(PETSC_COMM_WORLD, "]\n");CHKERRQ(ierr);
2553   } else {
2554     ierr = DMComputeL2Diff(dm, 0.0, exactFuncs, ctxs, u, &error[0]);CHKERRQ(ierr);
2555     if (error[0] >= 1.0e-11) {ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: %g\n", error[0]);CHKERRQ(ierr);}
2556     else                     {ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: < 1.0e-11\n");CHKERRQ(ierr);}
2557   }
2558   ierr = PetscFree(error);CHKERRQ(ierr);
2559   /* Check residual */
2560   ierr = SNESComputeFunction(snes, u, r);CHKERRQ(ierr);
2561   ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr);
2562   ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Residual: %g\n", res);CHKERRQ(ierr);
2563   ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr);
2564   ierr = PetscObjectSetName((PetscObject) r, "Initial Residual");CHKERRQ(ierr);
2565   ierr = PetscObjectSetOptionsPrefix((PetscObject)r,"res_");CHKERRQ(ierr);
2566   ierr = VecViewFromOptions(r, NULL, "-vec_view");CHKERRQ(ierr);
2567   /* Check Jacobian */
2568   ierr = SNESComputeJacobian(snes, u, M, M);CHKERRQ(ierr);
2569   ierr = MatGetNullSpace(J, &nullSpace);CHKERRQ(ierr);
2570   if (nullSpace) {
2571     PetscBool isNull;
2572     ierr = MatNullSpaceTest(nullSpace, J, &isNull);CHKERRQ(ierr);
2573     if (!isNull) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "The null space calculated for the system operator is invalid.");
2574   }
2575   ierr = VecDuplicate(u, &b);CHKERRQ(ierr);
2576   ierr = VecSet(r, 0.0);CHKERRQ(ierr);
2577   ierr = SNESComputeFunction(snes, r, b);CHKERRQ(ierr);
2578   ierr = MatMult(M, u, r);CHKERRQ(ierr);
2579   ierr = VecAXPY(r, 1.0, b);CHKERRQ(ierr);
2580   ierr = VecDestroy(&b);CHKERRQ(ierr);
2581   ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr);
2582   ierr = PetscPrintf(PETSC_COMM_WORLD, "Linear L_2 Residual: %g\n", res);CHKERRQ(ierr);
2583   ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr);
2584   ierr = PetscObjectSetName((PetscObject) r, "Au - b = Au + F(0)");CHKERRQ(ierr);
2585   ierr = PetscObjectSetOptionsPrefix((PetscObject)r,"linear_res_");CHKERRQ(ierr);
2586   ierr = VecViewFromOptions(r, NULL, "-vec_view");CHKERRQ(ierr);
2587   ierr = VecDestroy(&r);CHKERRQ(ierr);
2588   ierr = MatNullSpaceDestroy(&nullSpace);CHKERRQ(ierr);
2589   ierr = MatDestroy(&J);CHKERRQ(ierr);
2590   PetscFunctionReturn(0);
2591 }
2592 
2593 #undef __FUNCT__
2594 #define __FUNCT__ "DMSNESCheckFromOptions"
2595 PetscErrorCode DMSNESCheckFromOptions(SNES snes, Vec u, PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx), void **ctxs)
2596 {
2597   DM             dm;
2598   Vec            sol;
2599   PetscBool      check;
2600   PetscErrorCode ierr;
2601 
2602   PetscFunctionBegin;
2603   ierr = PetscOptionsHasName(((PetscObject)snes)->options,((PetscObject)snes)->prefix, "-dmsnes_check", &check);CHKERRQ(ierr);
2604   if (!check) PetscFunctionReturn(0);
2605   ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr);
2606   ierr = VecDuplicate(u, &sol);CHKERRQ(ierr);
2607   ierr = SNESSetSolution(snes, sol);CHKERRQ(ierr);
2608   ierr = DMSNESCheckFromOptions_Internal(snes, dm, u, sol, exactFuncs, ctxs);CHKERRQ(ierr);
2609   ierr = VecDestroy(&sol);CHKERRQ(ierr);
2610   PetscFunctionReturn(0);
2611 }
2612