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