xref: /petsc/src/snes/utils/dmplexsnes.c (revision dced61a5cfeeeda68dfa4ee3e53b34d3cfb8da9f)
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  - PetscViewer of type ASCII
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 = (PetscViewer) dummy;
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   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4);
781   ierr = SNESGetFunction(snes, &res, 0, 0);CHKERRQ(ierr);
782   ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr);
783   ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr);
784   ierr = PetscSectionGetNumFields(s, &numFields);CHKERRQ(ierr);
785   ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr);
786   ierr = PetscCalloc2(numFields, &lnorms, numFields, &norms);CHKERRQ(ierr);
787   ierr = VecGetArrayRead(res, &r);CHKERRQ(ierr);
788   for (p = pStart; p < pEnd; ++p) {
789     for (f = 0; f < numFields; ++f) {
790       PetscInt fdof, foff, d;
791 
792       ierr = PetscSectionGetFieldDof(s, p, f, &fdof);CHKERRQ(ierr);
793       ierr = PetscSectionGetFieldOffset(s, p, f, &foff);CHKERRQ(ierr);
794       for (d = 0; d < fdof; ++d) lnorms[f] += PetscRealPart(PetscSqr(r[foff+d]));
795     }
796   }
797   ierr = VecRestoreArrayRead(res, &r);CHKERRQ(ierr);
798   ierr = MPI_Allreduce(lnorms, norms, numFields, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject) dm));CHKERRQ(ierr);
799   ierr = PetscViewerASCIIAddTab(viewer, ((PetscObject) snes)->tablevel);CHKERRQ(ierr);
800   ierr = PetscViewerASCIIPrintf(viewer, "%3D SNES Function norm %14.12e [", its, (double) fgnorm);CHKERRQ(ierr);
801   for (f = 0; f < numFields; ++f) {
802     if (f > 0) {ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr);}
803     ierr = PetscViewerASCIIPrintf(viewer, "%14.12e", (double) PetscSqrtReal(norms[f]));CHKERRQ(ierr);
804   }
805   ierr = PetscViewerASCIIPrintf(viewer, "]\n");CHKERRQ(ierr);
806   ierr = PetscViewerASCIISubtractTab(viewer, ((PetscObject) snes)->tablevel);CHKERRQ(ierr);
807   ierr = PetscFree2(lnorms, norms);CHKERRQ(ierr);
808   PetscFunctionReturn(0);
809 }
810 
811 /********************* Residual Computation **************************/
812 
813 #undef __FUNCT__
814 #define __FUNCT__ "DMPlexSNESGetGeometryFEM"
815 /*@
816   DMPlexSNESGetGeometryFEM - Return precomputed geometric data
817 
818   Input Parameter:
819 . dm - The DM
820 
821   Output Parameters:
822 . cellgeom - The values precomputed from cell geometry
823 
824   Level: developer
825 
826 .seealso: DMPlexSNESSetFunctionLocal()
827 @*/
828 PetscErrorCode DMPlexSNESGetGeometryFEM(DM dm, Vec *cellgeom)
829 {
830   DMSNES         dmsnes;
831   PetscObject    obj;
832   PetscErrorCode ierr;
833 
834   PetscFunctionBegin;
835   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
836   ierr = DMGetDMSNES(dm, &dmsnes);CHKERRQ(ierr);
837   ierr = PetscObjectQuery((PetscObject) dmsnes, "DMPlexSNES_cellgeom_fem", &obj);CHKERRQ(ierr);
838   if (!obj) {
839     Vec cellgeom;
840 
841     ierr = DMPlexComputeGeometryFEM(dm, &cellgeom);CHKERRQ(ierr);
842     ierr = PetscObjectCompose((PetscObject) dmsnes, "DMPlexSNES_cellgeom_fem", (PetscObject) cellgeom);CHKERRQ(ierr);
843     ierr = VecDestroy(&cellgeom);CHKERRQ(ierr);
844   }
845   if (cellgeom) {PetscValidPointer(cellgeom, 3); ierr = PetscObjectQuery((PetscObject) dmsnes, "DMPlexSNES_cellgeom_fem", (PetscObject *) cellgeom);CHKERRQ(ierr);}
846   PetscFunctionReturn(0);
847 }
848 
849 #undef __FUNCT__
850 #define __FUNCT__ "DMPlexSNESGetGeometryFVM"
851 /*@
852   DMPlexSNESGetGeometryFVM - Return precomputed geometric data
853 
854   Input Parameter:
855 . dm - The DM
856 
857   Output Parameters:
858 + facegeom - The values precomputed from face geometry
859 . cellgeom - The values precomputed from cell geometry
860 - minRadius - The minimum radius over the mesh of an inscribed sphere in a cell
861 
862   Level: developer
863 
864 .seealso: DMPlexTSSetRHSFunctionLocal()
865 @*/
866 PetscErrorCode DMPlexSNESGetGeometryFVM(DM dm, Vec *facegeom, Vec *cellgeom, PetscReal *minRadius)
867 {
868   DMSNES         dmsnes;
869   PetscObject    obj;
870   PetscErrorCode ierr;
871 
872   PetscFunctionBegin;
873   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
874   ierr = DMGetDMSNES(dm, &dmsnes);CHKERRQ(ierr);
875   ierr = PetscObjectQuery((PetscObject) dmsnes, "DMPlexSNES_facegeom_fvm", &obj);CHKERRQ(ierr);
876   if (!obj) {
877     Vec cellgeom, facegeom;
878 
879     ierr = DMPlexComputeGeometryFVM(dm, &cellgeom, &facegeom);CHKERRQ(ierr);
880     ierr = PetscObjectCompose((PetscObject) dmsnes, "DMPlexSNES_facegeom_fvm", (PetscObject) facegeom);CHKERRQ(ierr);
881     ierr = PetscObjectCompose((PetscObject) dmsnes, "DMPlexSNES_cellgeom_fvm", (PetscObject) cellgeom);CHKERRQ(ierr);
882     ierr = VecDestroy(&facegeom);CHKERRQ(ierr);
883     ierr = VecDestroy(&cellgeom);CHKERRQ(ierr);
884   }
885   if (facegeom) {PetscValidPointer(facegeom, 2); ierr = PetscObjectQuery((PetscObject) dmsnes, "DMPlexSNES_facegeom_fvm", (PetscObject *) facegeom);CHKERRQ(ierr);}
886   if (cellgeom) {PetscValidPointer(cellgeom, 3); ierr = PetscObjectQuery((PetscObject) dmsnes, "DMPlexSNES_cellgeom_fvm", (PetscObject *) cellgeom);CHKERRQ(ierr);}
887   if (minRadius) {ierr = DMPlexGetMinRadius(dm, minRadius);CHKERRQ(ierr);}
888   PetscFunctionReturn(0);
889 }
890 
891 #undef __FUNCT__
892 #define __FUNCT__ "DMPlexSNESGetGradientDM"
893 /*@
894   DMPlexSNESGetGradientDM - Return gradient data layout
895 
896   Input Parameters:
897 + dm - The DM
898 - fv - The PetscFV
899 
900   Output Parameter:
901 . dmGrad - The layout for gradient values
902 
903   Level: developer
904 
905 .seealso: DMPlexSNESGetGeometryFVM()
906 @*/
907 PetscErrorCode DMPlexSNESGetGradientDM(DM dm, PetscFV fv, DM *dmGrad)
908 {
909   DMSNES         dmsnes;
910   PetscObject    obj;
911   PetscBool      computeGradients;
912   PetscErrorCode ierr;
913 
914   PetscFunctionBegin;
915   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
916   PetscValidHeaderSpecific(fv,PETSCFV_CLASSID,2);
917   PetscValidPointer(dmGrad,3);
918   ierr = PetscFVGetComputeGradients(fv, &computeGradients);CHKERRQ(ierr);
919   if (!computeGradients) {*dmGrad = NULL; PetscFunctionReturn(0);}
920   ierr = DMGetDMSNES(dm, &dmsnes);CHKERRQ(ierr);
921   ierr = PetscObjectQuery((PetscObject) dmsnes, "DMPlexSNES_dmgrad_fvm", &obj);CHKERRQ(ierr);
922   if (!obj) {
923     DM  dmGrad;
924     Vec faceGeometry, cellGeometry;
925 
926     ierr = DMPlexSNESGetGeometryFVM(dm, &faceGeometry, &cellGeometry, NULL);CHKERRQ(ierr);
927     ierr = DMPlexComputeGradientFVM(dm, fv, faceGeometry, cellGeometry, &dmGrad);CHKERRQ(ierr);
928     ierr = PetscObjectCompose((PetscObject) dmsnes, "DMPlexSNES_dmgrad_fvm", (PetscObject) dmGrad);CHKERRQ(ierr);
929     ierr = DMDestroy(&dmGrad);CHKERRQ(ierr);
930   }
931   ierr = PetscObjectQuery((PetscObject) dmsnes, "DMPlexSNES_dmgrad_fvm", (PetscObject *) dmGrad);CHKERRQ(ierr);
932   PetscFunctionReturn(0);
933 }
934 
935 #undef __FUNCT__
936 #define __FUNCT__ "DMPlexGetCellFields"
937 /*@C
938   DMPlexGetCellFields - Retrieve the field values values for a chunk of cells
939 
940   Input Parameters:
941 + dm     - The DM
942 . cStart - The first cell to include
943 . cEnd   - The first cell to exclude
944 . locX   - A local vector with the solution fields
945 . locX_t - A local vector with solution field time derivatives, or NULL
946 - locA   - A local vector with auxiliary fields, or NULL
947 
948   Output Parameters:
949 + u   - The field coefficients
950 . u_t - The fields derivative coefficients
951 - a   - The auxiliary field coefficients
952 
953   Level: developer
954 
955 .seealso: DMPlexGetFaceFields()
956 @*/
957 PetscErrorCode DMPlexGetCellFields(DM dm, PetscInt cStart, PetscInt cEnd, Vec locX, Vec locX_t, Vec locA, PetscScalar **u, PetscScalar **u_t, PetscScalar **a)
958 {
959   DM             dmAux;
960   PetscSection   section, sectionAux;
961   PetscDS        prob;
962   PetscInt       numCells = cEnd - cStart, totDim, totDimAux, c;
963   PetscErrorCode ierr;
964 
965   PetscFunctionBegin;
966   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
967   PetscValidHeaderSpecific(locX, VEC_CLASSID, 4);
968   if (locX_t) {PetscValidHeaderSpecific(locX_t, VEC_CLASSID, 5);}
969   if (locA)   {PetscValidHeaderSpecific(locA, VEC_CLASSID, 6);}
970   PetscValidPointer(u, 7);
971   PetscValidPointer(u_t, 8);
972   PetscValidPointer(a, 9);
973   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
974   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
975   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
976   if (locA) {
977     PetscDS probAux;
978 
979     ierr = VecGetDM(locA, &dmAux);CHKERRQ(ierr);
980     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
981     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
982     ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
983   }
984   ierr = DMGetWorkArray(dm, numCells*totDim, PETSC_SCALAR, u);CHKERRQ(ierr);
985   if (locX_t) {ierr = DMGetWorkArray(dm, numCells*totDim, PETSC_SCALAR, u_t);CHKERRQ(ierr);} else {*u_t = NULL;}
986   if (locA)   {ierr = DMGetWorkArray(dm, numCells*totDimAux, PETSC_SCALAR, a);CHKERRQ(ierr);} else {*a = NULL;}
987   for (c = cStart; c < cEnd; ++c) {
988     PetscScalar *x = NULL, *x_t = NULL, *ul = *u, *ul_t = *u_t, *al = *a;
989     PetscInt     i;
990 
991     ierr = DMPlexVecGetClosure(dm, section, locX, c, NULL, &x);CHKERRQ(ierr);
992     for (i = 0; i < totDim; ++i) ul[(c-cStart)*totDim+i] = x[i];
993     ierr = DMPlexVecRestoreClosure(dm, section, locX, c, NULL, &x);CHKERRQ(ierr);
994     if (locX_t) {
995       ierr = DMPlexVecGetClosure(dm, section, locX_t, c, NULL, &x_t);CHKERRQ(ierr);
996       for (i = 0; i < totDim; ++i) ul_t[(c-cStart)*totDim+i] = x_t[i];
997       ierr = DMPlexVecRestoreClosure(dm, section, locX_t, c, NULL, &x_t);CHKERRQ(ierr);
998     }
999     if (locA) {
1000       ierr = DMPlexVecGetClosure(dmAux, sectionAux, locA, c, NULL, &x);CHKERRQ(ierr);
1001       for (i = 0; i < totDimAux; ++i) al[(c-cStart)*totDimAux+i] = x[i];
1002       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, locA, c, NULL, &x);CHKERRQ(ierr);
1003     }
1004   }
1005   PetscFunctionReturn(0);
1006 }
1007 
1008 #undef __FUNCT__
1009 #define __FUNCT__ "DMPlexRestoreCellFields"
1010 /*@C
1011   DMPlexRestoreCellFields - Restore the field values values for a chunk of cells
1012 
1013   Input Parameters:
1014 + dm     - The DM
1015 . cStart - The first cell to include
1016 . cEnd   - The first cell to exclude
1017 . locX   - A local vector with the solution fields
1018 . locX_t - A local vector with solution field time derivatives, or NULL
1019 - locA   - A local vector with auxiliary fields, or NULL
1020 
1021   Output Parameters:
1022 + u   - The field coefficients
1023 . u_t - The fields derivative coefficients
1024 - a   - The auxiliary field coefficients
1025 
1026   Level: developer
1027 
1028 .seealso: DMPlexGetFaceFields()
1029 @*/
1030 PetscErrorCode DMPlexRestoreCellFields(DM dm, PetscInt cStart, PetscInt cEnd, Vec locX, Vec locX_t, Vec locA, PetscScalar **u, PetscScalar **u_t, PetscScalar **a)
1031 {
1032   PetscErrorCode ierr;
1033 
1034   PetscFunctionBegin;
1035   ierr = DMRestoreWorkArray(dm, 0, PETSC_SCALAR, u);CHKERRQ(ierr);
1036   if (*u_t) {ierr = DMRestoreWorkArray(dm, 0, PETSC_SCALAR, u_t);CHKERRQ(ierr);}
1037   if (*a)   {ierr = DMRestoreWorkArray(dm, 0, PETSC_SCALAR, a);CHKERRQ(ierr);}
1038   PetscFunctionReturn(0);
1039 }
1040 
1041 #undef __FUNCT__
1042 #define __FUNCT__ "DMPlexGetFaceFields"
1043 /*@C
1044   DMPlexGetFaceFields - Retrieve the field values values for a chunk of faces
1045 
1046   Input Parameters:
1047 + dm     - The DM
1048 . fStart - The first face to include
1049 . fEnd   - The first face to exclude
1050 . locX   - A local vector with the solution fields
1051 . locX_t - A local vector with solution field time derivatives, or NULL
1052 . faceGeometry - A local vector with face geometry
1053 . cellGeometry - A local vector with cell geometry
1054 - locaGrad - A local vector with field gradients, or NULL
1055 
1056   Output Parameters:
1057 + uL - The field values at the left side of the face
1058 - uR - The field values at the right side of the face
1059 
1060   Level: developer
1061 
1062 .seealso: DMPlexGetCellFields()
1063 @*/
1064 PetscErrorCode DMPlexGetFaceFields(DM dm, PetscInt fStart, PetscInt fEnd, Vec locX, Vec locX_t, Vec faceGeometry, Vec cellGeometry, Vec locGrad, PetscScalar **uL, PetscScalar **uR)
1065 {
1066   DM                 dmFace, dmCell, dmGrad = NULL;
1067   PetscSection       section;
1068   PetscDS            prob;
1069   DMLabel            ghostLabel;
1070   const PetscScalar *facegeom, *cellgeom, *x, *lgrad;
1071   PetscBool         *isFE;
1072   PetscInt           dim, Nf, f, Nc, numFaces = fEnd - fStart, iface, face;
1073   PetscErrorCode     ierr;
1074 
1075   PetscFunctionBegin;
1076   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1077   PetscValidHeaderSpecific(locX, VEC_CLASSID, 4);
1078   if (locX_t) {PetscValidHeaderSpecific(locX_t, VEC_CLASSID, 5);}
1079   PetscValidHeaderSpecific(faceGeometry, VEC_CLASSID, 6);
1080   PetscValidHeaderSpecific(cellGeometry, VEC_CLASSID, 7);
1081   if (locGrad) {PetscValidHeaderSpecific(locGrad, VEC_CLASSID, 8);}
1082   PetscValidPointer(uL, 9);
1083   PetscValidPointer(uR, 10);
1084   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1085   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
1086   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1087   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
1088   ierr = PetscDSGetTotalComponents(prob, &Nc);CHKERRQ(ierr);
1089   ierr = PetscMalloc1(Nf, &isFE);CHKERRQ(ierr);
1090   for (f = 0; f < Nf; ++f) {
1091     PetscObject  obj;
1092     PetscClassId id;
1093 
1094     ierr = DMGetField(dm, f, &obj);CHKERRQ(ierr);
1095     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1096     if (id == PETSCFE_CLASSID)      {isFE[f] = PETSC_TRUE;}
1097     else if (id == PETSCFV_CLASSID) {isFE[f] = PETSC_FALSE;}
1098     else                            {isFE[f] = PETSC_FALSE;}
1099   }
1100   ierr = DMPlexGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr);
1101   ierr = VecGetArrayRead(locX, &x);CHKERRQ(ierr);
1102   ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr);
1103   ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
1104   ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr);
1105   ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
1106   if (locGrad) {
1107     ierr = VecGetDM(locGrad, &dmGrad);CHKERRQ(ierr);
1108     ierr = VecGetArrayRead(locGrad, &lgrad);CHKERRQ(ierr);
1109   }
1110   ierr = DMGetWorkArray(dm, numFaces*Nc, PETSC_SCALAR, uL);CHKERRQ(ierr);
1111   ierr = DMGetWorkArray(dm, numFaces*Nc, PETSC_SCALAR, uR);CHKERRQ(ierr);
1112   /* Right now just eat the extra work for FE (could make a cell loop) */
1113   for (face = fStart, iface = 0; face < fEnd; ++face) {
1114     const PetscInt        *cells;
1115     const PetscFVFaceGeom *fg;
1116     const PetscFVCellGeom *cgL, *cgR;
1117     const PetscScalar     *xL, *xR, *gL, *gR;
1118     PetscScalar           *uLl = *uL, *uRl = *uR;
1119     PetscInt               ghost;
1120 
1121     ierr = DMLabelGetValue(ghostLabel, face, &ghost);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         ierr = DMPlexPointLocalFieldRead(dm, cells[0], f, x, &xL);CHKERRQ(ierr);
1159         ierr = DMPlexPointLocalFieldRead(dm, cells[1], f, x, &xR);CHKERRQ(ierr);
1160         if (dmGrad) {
1161           PetscReal dxL[3], dxR[3];
1162 
1163           ierr = DMPlexPointLocalRead(dmGrad, cells[0], lgrad, &gL);CHKERRQ(ierr);
1164           ierr = DMPlexPointLocalRead(dmGrad, cells[1], lgrad, &gR);CHKERRQ(ierr);
1165           DMPlex_WaxpyD_Internal(dim, -1, cgL->centroid, fg->centroid, dxL);
1166           DMPlex_WaxpyD_Internal(dim, -1, cgR->centroid, fg->centroid, dxR);
1167           for (c = 0; c < numComp; ++c) {
1168             uLl[iface*Nc+off+c] = xL[c] + DMPlex_DotD_Internal(dim, &gL[c*dim], dxL);
1169             uRl[iface*Nc+off+c] = xR[c] + DMPlex_DotD_Internal(dim, &gR[c*dim], dxR);
1170           }
1171         } else {
1172           for (c = 0; c < numComp; ++c) {
1173             uLl[iface*Nc+off+c] = xL[c];
1174             uRl[iface*Nc+off+c] = xR[c];
1175           }
1176         }
1177       }
1178     }
1179     ++iface;
1180   }
1181   ierr = VecRestoreArrayRead(locX, &x);CHKERRQ(ierr);
1182   ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
1183   ierr = VecRestoreArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
1184   if (locGrad) {
1185     ierr = VecRestoreArrayRead(locGrad, &lgrad);CHKERRQ(ierr);
1186   }
1187   ierr = PetscFree(isFE);CHKERRQ(ierr);
1188   PetscFunctionReturn(0);
1189 }
1190 
1191 #undef __FUNCT__
1192 #define __FUNCT__ "DMPlexRestoreFaceFields"
1193 /*@C
1194   DMPlexRestoreFaceFields - Restore the field values values for a chunk of faces
1195 
1196   Input Parameters:
1197 + dm     - The DM
1198 . fStart - The first face to include
1199 . fEnd   - The first face to exclude
1200 . locX   - A local vector with the solution fields
1201 . locX_t - A local vector with solution field time derivatives, or NULL
1202 . faceGeometry - A local vector with face geometry
1203 . cellGeometry - A local vector with cell geometry
1204 - locaGrad - A local vector with field gradients, or NULL
1205 
1206   Output Parameters:
1207 + uL - The field values at the left side of the face
1208 - uR - The field values at the right side of the face
1209 
1210   Level: developer
1211 
1212 .seealso: DMPlexGetFaceFields()
1213 @*/
1214 PetscErrorCode DMPlexRestoreFaceFields(DM dm, PetscInt fStart, PetscInt fEnd, Vec locX, Vec locX_t, Vec faceGeometry, Vec cellGeometry, Vec locGrad, PetscScalar **uL, PetscScalar **uR)
1215 {
1216   PetscErrorCode ierr;
1217 
1218   PetscFunctionBegin;
1219   ierr = DMRestoreWorkArray(dm, 0, PETSC_SCALAR, uL);CHKERRQ(ierr);
1220   ierr = DMRestoreWorkArray(dm, 0, PETSC_SCALAR, uR);CHKERRQ(ierr);
1221   PetscFunctionReturn(0);
1222 }
1223 
1224 #undef __FUNCT__
1225 #define __FUNCT__ "DMPlexGetFaceGeometry"
1226 /*@C
1227   DMPlexGetFaceGeometry - Retrieve the geometric values for a chunk of faces
1228 
1229   Input Parameters:
1230 + dm     - The DM
1231 . fStart - The first face to include
1232 . fEnd   - The first face to exclude
1233 . faceGeometry - A local vector with face geometry
1234 - cellGeometry - A local vector with cell geometry
1235 
1236   Output Parameters:
1237 + fgeom - The extract the face centroid and normal
1238 - vol   - The cell volume
1239 
1240   Level: developer
1241 
1242 .seealso: DMPlexGetCellFields()
1243 @*/
1244 PetscErrorCode DMPlexGetFaceGeometry(DM dm, PetscInt fStart, PetscInt fEnd, Vec faceGeometry, Vec cellGeometry, PetscFVFaceGeom **fgeom, PetscReal **vol)
1245 {
1246   DM                 dmFace, dmCell;
1247   DMLabel            ghostLabel;
1248   const PetscScalar *facegeom, *cellgeom;
1249   PetscInt           dim, numFaces = fEnd - fStart, iface, face;
1250   PetscErrorCode     ierr;
1251 
1252   PetscFunctionBegin;
1253   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1254   PetscValidHeaderSpecific(faceGeometry, VEC_CLASSID, 4);
1255   PetscValidHeaderSpecific(cellGeometry, VEC_CLASSID, 5);
1256   PetscValidPointer(fgeom, 6);
1257   PetscValidPointer(vol, 7);
1258   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1259   ierr = DMPlexGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr);
1260   ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr);
1261   ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
1262   ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr);
1263   ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
1264   ierr = PetscMalloc1(numFaces, fgeom);CHKERRQ(ierr);
1265   ierr = DMGetWorkArray(dm, numFaces*2, PETSC_SCALAR, vol);CHKERRQ(ierr);
1266   for (face = fStart, iface = 0; face < fEnd; ++face) {
1267     const PetscInt        *cells;
1268     const PetscFVFaceGeom *fg;
1269     const PetscFVCellGeom *cgL, *cgR;
1270     PetscFVFaceGeom       *fgeoml = *fgeom;
1271     PetscReal             *voll   = *vol;
1272     PetscInt               ghost, d;
1273 
1274     ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr);
1275     if (ghost >= 0) continue;
1276     ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr);
1277     ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr);
1278     ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cgL);CHKERRQ(ierr);
1279     ierr = DMPlexPointLocalRead(dmCell, cells[1], cellgeom, &cgR);CHKERRQ(ierr);
1280     for (d = 0; d < dim; ++d) {
1281       fgeoml[iface].centroid[d] = fg->centroid[d];
1282       fgeoml[iface].normal[d]   = fg->normal[d];
1283     }
1284     voll[iface*2+0] = cgL->volume;
1285     voll[iface*2+1] = cgR->volume;
1286     ++iface;
1287   }
1288   ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
1289   ierr = VecRestoreArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
1290   PetscFunctionReturn(0);
1291 }
1292 
1293 #undef __FUNCT__
1294 #define __FUNCT__ "DMPlexRestoreFaceGeometry"
1295 /*@C
1296   DMPlexRestoreFaceGeometry - Restore the field values values for a chunk of faces
1297 
1298   Input Parameters:
1299 + dm     - The DM
1300 . fStart - The first face to include
1301 . fEnd   - The first face to exclude
1302 . faceGeometry - A local vector with face geometry
1303 - cellGeometry - A local vector with cell geometry
1304 
1305   Output Parameters:
1306 + fgeom - The extract the face centroid and normal
1307 - vol   - The cell volume
1308 
1309   Level: developer
1310 
1311 .seealso: DMPlexGetFaceFields()
1312 @*/
1313 PetscErrorCode DMPlexRestoreFaceGeometry(DM dm, PetscInt fStart, PetscInt fEnd, Vec faceGeometry, Vec cellGeometry, PetscFVFaceGeom **fgeom, PetscReal **vol)
1314 {
1315   PetscErrorCode ierr;
1316 
1317   PetscFunctionBegin;
1318   ierr = PetscFree(*fgeom);CHKERRQ(ierr);
1319   ierr = DMRestoreWorkArray(dm, 0, PETSC_REAL, vol);CHKERRQ(ierr);
1320   PetscFunctionReturn(0);
1321 }
1322 
1323 #undef __FUNCT__
1324 #define __FUNCT__ "DMPlexReconstructGradients_Internal"
1325 PetscErrorCode DMPlexReconstructGradients_Internal(DM dm, PetscInt fStart, PetscInt fEnd, Vec faceGeometry, Vec cellGeometry, Vec locX, Vec grad)
1326 {
1327   DM                 dmFace, dmCell, dmGrad;
1328   DMLabel            ghostLabel;
1329   PetscDS            prob;
1330   PetscFV            fvm;
1331   PetscLimiter       lim;
1332   const PetscScalar *facegeom, *cellgeom, *x;
1333   PetscScalar       *gr;
1334   PetscReal         *cellPhi;
1335   PetscInt           dim, face, cell, totDim, cStart, cEnd, cEndInterior;
1336   PetscErrorCode     ierr;
1337 
1338   PetscFunctionBegin;
1339   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1340   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
1341   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
1342   ierr = DMPlexGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr);
1343   ierr = PetscDSGetDiscretization(prob, 0, (PetscObject *) &fvm);CHKERRQ(ierr);
1344   ierr = PetscFVGetLimiter(fvm, &lim);CHKERRQ(ierr);
1345   ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr);
1346   ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
1347   ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr);
1348   ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
1349   ierr = VecGetArrayRead(locX, &x);CHKERRQ(ierr);
1350   ierr = VecGetDM(grad, &dmGrad);CHKERRQ(ierr);
1351   ierr = VecZeroEntries(grad);CHKERRQ(ierr);
1352   ierr = VecGetArray(grad, &gr);CHKERRQ(ierr);
1353   /* Reconstruct gradients */
1354   for (face = fStart; face < fEnd; ++face) {
1355     const PetscInt        *cells;
1356     const PetscFVFaceGeom *fg;
1357     const PetscScalar     *cx[2];
1358     PetscScalar           *cgrad[2];
1359     PetscBool              boundary;
1360     PetscInt               ghost, c, pd, d;
1361 
1362     ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr);
1363     if (ghost >= 0) continue;
1364     ierr = DMPlexIsBoundaryPoint(dm, face, &boundary);CHKERRQ(ierr);
1365     if (boundary) continue;
1366     ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr);
1367     ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr);
1368     for (c = 0; c < 2; ++c) {
1369       ierr = DMPlexPointLocalRead(dm, cells[c], x, &cx[c]);CHKERRQ(ierr);
1370       ierr = DMPlexPointGlobalRef(dmGrad, cells[c], gr, &cgrad[c]);CHKERRQ(ierr);
1371     }
1372     for (pd = 0; pd < totDim; ++pd) {
1373       PetscScalar delta = cx[1][pd] - cx[0][pd];
1374 
1375       for (d = 0; d < dim; ++d) {
1376         if (cgrad[0]) cgrad[0][pd*dim+d] += fg->grad[0][d] * delta;
1377         if (cgrad[1]) cgrad[1][pd*dim+d] -= fg->grad[1][d] * delta;
1378       }
1379     }
1380   }
1381   /* Limit interior gradients (using cell-based loop because it generalizes better to vector limiters) */
1382   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1383   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
1384   cEndInterior = cEndInterior < 0 ? cEnd : cEndInterior;
1385   ierr = DMGetWorkArray(dm, totDim, PETSC_REAL, &cellPhi);CHKERRQ(ierr);
1386   for (cell = dmGrad && lim ? cStart : cEnd; cell < cEndInterior; ++cell) {
1387     const PetscInt        *faces;
1388     const PetscScalar     *cx;
1389     const PetscFVCellGeom *cg;
1390     PetscScalar           *cgrad;
1391     PetscInt               coneSize, f, pd, d;
1392 
1393     ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
1394     ierr = DMPlexGetCone(dm, cell, &faces);CHKERRQ(ierr);
1395     ierr = DMPlexPointLocalRead(dm, cell, x, &cx);CHKERRQ(ierr);
1396     ierr = DMPlexPointLocalRead(dmCell, cell, cellgeom, &cg);CHKERRQ(ierr);
1397     ierr = DMPlexPointGlobalRef(dmGrad, cell, gr, &cgrad);CHKERRQ(ierr);
1398     if (!cgrad) continue; /* Unowned overlap cell, we do not compute */
1399     /* Limiter will be minimum value over all neighbors */
1400     for (d = 0; d < totDim; ++d) cellPhi[d] = PETSC_MAX_REAL;
1401     for (f = 0; f < coneSize; ++f) {
1402       const PetscScalar     *ncx;
1403       const PetscFVCellGeom *ncg;
1404       const PetscInt        *fcells;
1405       PetscInt               face = faces[f], ncell, ghost;
1406       PetscReal              v[3];
1407       PetscBool              boundary;
1408 
1409       ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr);
1410       ierr = DMPlexIsBoundaryPoint(dm, face, &boundary);CHKERRQ(ierr);
1411       if ((ghost >= 0) || boundary) continue;
1412       ierr  = DMPlexGetSupport(dm, face, &fcells);CHKERRQ(ierr);
1413       ncell = cell == fcells[0] ? fcells[1] : fcells[0];
1414       ierr  = DMPlexPointLocalRead(dm, ncell, x, &ncx);CHKERRQ(ierr);
1415       ierr  = DMPlexPointLocalRead(dmCell, ncell, cellgeom, &ncg);CHKERRQ(ierr);
1416       DMPlex_WaxpyD_Internal(dim, -1, cg->centroid, ncg->centroid, v);
1417       for (d = 0; d < totDim; ++d) {
1418         /* We use the symmetric slope limited form of Berger, Aftosmis, and Murman 2005 */
1419         PetscReal phi, flim = 0.5 * PetscRealPart(ncx[d] - cx[d]) / DMPlex_DotD_Internal(dim, &cgrad[d*dim], v);
1420 
1421         ierr = PetscLimiterLimit(lim, flim, &phi);CHKERRQ(ierr);
1422         cellPhi[d] = PetscMin(cellPhi[d], phi);
1423       }
1424     }
1425     /* Apply limiter to gradient */
1426     for (pd = 0; pd < totDim; ++pd)
1427       /* Scalar limiter applied to each component separately */
1428       for (d = 0; d < dim; ++d) cgrad[pd*dim+d] *= cellPhi[pd];
1429   }
1430   ierr = DMRestoreWorkArray(dm, totDim, PETSC_REAL, &cellPhi);CHKERRQ(ierr);
1431   ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
1432   ierr = VecRestoreArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
1433   ierr = VecRestoreArrayRead(locX, &x);CHKERRQ(ierr);
1434   ierr = VecRestoreArray(grad, &gr);CHKERRQ(ierr);
1435   PetscFunctionReturn(0);
1436 }
1437 
1438 #undef __FUNCT__
1439 #define __FUNCT__ "DMPlexComputeBdResidual_Internal"
1440 PetscErrorCode DMPlexComputeBdResidual_Internal(DM dm, Vec locX, Vec locX_t, Vec locF, void *user)
1441 {
1442   DM_Plex         *mesh = (DM_Plex *) dm->data;
1443   PetscSection     section;
1444   PetscDS          prob;
1445   DMLabel          depth;
1446   PetscFECellGeom *cgeom;
1447   PetscScalar     *u = NULL, *u_t = NULL, *elemVec = NULL;
1448   PetscInt         dim, Nf, f, totDimBd, numBd, bd;
1449   PetscErrorCode   ierr;
1450 
1451   PetscFunctionBegin;
1452   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1453   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1454   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
1455   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
1456   ierr = PetscDSGetTotalBdDimension(prob, &totDimBd);CHKERRQ(ierr);
1457   ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr);
1458   ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
1459   for (bd = 0; bd < numBd; ++bd) {
1460     const char     *bdLabel;
1461     DMLabel         label;
1462     IS              pointIS;
1463     const PetscInt *points;
1464     const PetscInt *values;
1465     PetscInt        field, numValues, numPoints, p, dep, numFaces;
1466     PetscBool       isEssential;
1467 
1468     ierr = DMPlexGetBoundary(dm, bd, &isEssential, NULL, &bdLabel, &field, NULL, NULL, NULL, &numValues, &values, NULL);CHKERRQ(ierr);
1469     if (isEssential) continue;
1470     if (numValues != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Bug me and I will fix this");
1471     ierr = DMPlexGetLabel(dm, bdLabel, &label);CHKERRQ(ierr);
1472     ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr);
1473     ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr);
1474     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
1475     for (p = 0, numFaces = 0; p < numPoints; ++p) {
1476       ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
1477       if (dep == dim-1) ++numFaces;
1478     }
1479     ierr = PetscMalloc3(numFaces*totDimBd,&u,numFaces,&cgeom,numFaces*totDimBd,&elemVec);CHKERRQ(ierr);
1480     if (locX_t) {ierr = PetscMalloc1(numFaces*totDimBd,&u_t);CHKERRQ(ierr);}
1481     for (p = 0, f = 0; p < numPoints; ++p) {
1482       const PetscInt point = points[p];
1483       PetscScalar   *x     = NULL;
1484       PetscInt       i;
1485 
1486       ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
1487       if (dep != dim-1) continue;
1488       ierr = DMPlexComputeCellGeometryFEM(dm, point, NULL, cgeom[f].v0, cgeom[f].J, cgeom[f].invJ, &cgeom[f].detJ);CHKERRQ(ierr);
1489       ierr = DMPlexComputeCellGeometryFVM(dm, point, NULL, NULL, cgeom[f].n);CHKERRQ(ierr);
1490       if (cgeom[f].detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", cgeom[f].detJ, point);
1491       ierr = DMPlexVecGetClosure(dm, section, locX, point, NULL, &x);CHKERRQ(ierr);
1492       for (i = 0; i < totDimBd; ++i) u[f*totDimBd+i] = x[i];
1493       ierr = DMPlexVecRestoreClosure(dm, section, locX, point, NULL, &x);CHKERRQ(ierr);
1494       if (locX_t) {
1495         ierr = DMPlexVecGetClosure(dm, section, locX_t, point, NULL, &x);CHKERRQ(ierr);
1496         for (i = 0; i < totDimBd; ++i) u_t[f*totDimBd+i] = x[i];
1497         ierr = DMPlexVecRestoreClosure(dm, section, locX_t, point, NULL, &x);CHKERRQ(ierr);
1498       }
1499       ++f;
1500     }
1501     for (f = 0; f < Nf; ++f) {
1502       PetscFE         fe;
1503       PetscQuadrature q;
1504       PetscInt        numQuadPoints, Nb;
1505       /* Conforming batches */
1506       PetscInt        numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1507       /* Remainder */
1508       PetscInt        Nr, offset;
1509 
1510       ierr = PetscDSGetBdDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
1511       ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
1512       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
1513       ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
1514       ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
1515       blockSize = Nb*numQuadPoints;
1516       batchSize = numBlocks * blockSize;
1517       ierr =  PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1518       numChunks = numFaces / (numBatches*batchSize);
1519       Ne        = numChunks*numBatches*batchSize;
1520       Nr        = numFaces % (numBatches*batchSize);
1521       offset    = numFaces - Nr;
1522       ierr = PetscFEIntegrateBdResidual(fe, prob, f, Ne, cgeom, u, u_t, NULL, NULL, elemVec);CHKERRQ(ierr);
1523       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);
1524     }
1525     for (p = 0, f = 0; p < numPoints; ++p) {
1526       const PetscInt point = points[p];
1527 
1528       ierr = DMLabelGetValue(depth, point, &dep);CHKERRQ(ierr);
1529       if (dep != dim-1) continue;
1530       if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "BdResidual", totDimBd, &elemVec[f*totDimBd]);CHKERRQ(ierr);}
1531       ierr = DMPlexVecSetClosure(dm, NULL, locF, point, &elemVec[f*totDimBd], ADD_VALUES);CHKERRQ(ierr);
1532       ++f;
1533     }
1534     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
1535     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
1536     ierr = PetscFree3(u,cgeom,elemVec);CHKERRQ(ierr);
1537     if (locX_t) {ierr = PetscFree(u_t);CHKERRQ(ierr);}
1538   }
1539   PetscFunctionReturn(0);
1540 }
1541 
1542 #undef __FUNCT__
1543 #define __FUNCT__ "DMPlexComputeResidual_Internal"
1544 PetscErrorCode DMPlexComputeResidual_Internal(DM dm, PetscInt cStart, PetscInt cEnd, PetscReal time, Vec locX, Vec locX_t, Vec locF, void *user)
1545 {
1546   DM_Plex          *mesh       = (DM_Plex *) dm->data;
1547   const char       *name       = "Residual";
1548   DM                dmAux      = NULL;
1549   DM                dmGrad     = NULL;
1550   DMLabel           ghostLabel = NULL;
1551   PetscDS           prob       = NULL;
1552   PetscDS           probAux    = NULL;
1553   PetscSection      section    = NULL;
1554   PetscBool         useFEM     = PETSC_FALSE;
1555   PetscBool         useFVM     = PETSC_FALSE;
1556   PetscBool         isImplicit = (locX_t || time == PETSC_MIN_REAL) ? PETSC_TRUE : PETSC_FALSE;
1557   PetscFV           fvm        = NULL;
1558   PetscFECellGeom  *cgeomFEM   = NULL;
1559   PetscFVCellGeom  *cgeomFVM   = NULL;
1560   PetscFVFaceGeom  *fgeomFVM   = NULL;
1561   Vec               locA, cellGeometryFEM = NULL, cellGeometryFVM = NULL, faceGeometryFVM = NULL, grad, locGrad = NULL;
1562   PetscScalar      *u, *u_t, *a, *uL, *uR;
1563   PetscInt          Nf, f, totDim, totDimAux, numChunks, cellChunkSize, faceChunkSize, chunk, fStart, fEnd;
1564   PetscErrorCode    ierr;
1565 
1566   PetscFunctionBegin;
1567   ierr = PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
1568   /* TODO The places where we have to use isFE are probably the member functions for the PetscDisc class */
1569   /* TODO The FVM geometry is over-manipulated. Make the precalc functions return exactly what we need */
1570   /* FEM+FVM */
1571   /* 1: Get sizes from dm and dmAux */
1572   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1573   ierr = DMPlexGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr);
1574   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
1575   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
1576   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
1577   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &locA);CHKERRQ(ierr);
1578   if (locA) {
1579     ierr = VecGetDM(locA, &dmAux);CHKERRQ(ierr);
1580     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
1581     ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
1582   }
1583   /* 2: Get geometric data */
1584   for (f = 0; f < Nf; ++f) {
1585     PetscObject  obj;
1586     PetscClassId id;
1587     PetscBool    fimp;
1588 
1589     ierr = PetscDSGetImplicit(prob, f, &fimp);CHKERRQ(ierr);
1590     if (isImplicit != fimp) continue;
1591     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
1592     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1593     if (id == PETSCFE_CLASSID) {useFEM = PETSC_TRUE;}
1594     if (id == PETSCFV_CLASSID) {useFVM = PETSC_TRUE; fvm = (PetscFV) obj;}
1595   }
1596   if (useFEM) {
1597     ierr = DMPlexSNESGetGeometryFEM(dm, &cellGeometryFEM);CHKERRQ(ierr);
1598     ierr = VecGetArrayRead(cellGeometryFEM, (const PetscScalar **) &cgeomFEM);CHKERRQ(ierr);
1599   }
1600   if (useFVM) {
1601     ierr = DMPlexSNESGetGeometryFVM(dm, &faceGeometryFVM, &cellGeometryFVM, NULL);CHKERRQ(ierr);
1602     ierr = VecGetArrayRead(faceGeometryFVM, (const PetscScalar **) &fgeomFVM);CHKERRQ(ierr);
1603     ierr = VecGetArrayRead(cellGeometryFVM, (const PetscScalar **) &cgeomFVM);CHKERRQ(ierr);
1604     /* Reconstruct and limit cell gradients */
1605     ierr = DMPlexSNESGetGradientDM(dm, fvm, &dmGrad);CHKERRQ(ierr);
1606     if (dmGrad) {
1607       ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
1608       ierr = DMGetGlobalVector(dmGrad, &grad);CHKERRQ(ierr);
1609       ierr = DMPlexReconstructGradients_Internal(dm, fStart, fEnd, faceGeometryFVM, cellGeometryFVM, locX, grad);CHKERRQ(ierr);
1610       /* Communicate gradient values */
1611       ierr = DMGetLocalVector(dmGrad, &locGrad);CHKERRQ(ierr);
1612       ierr = DMGlobalToLocalBegin(dmGrad, grad, INSERT_VALUES, locGrad);CHKERRQ(ierr);
1613       ierr = DMGlobalToLocalEnd(dmGrad, grad, INSERT_VALUES, locGrad);CHKERRQ(ierr);
1614       ierr = DMRestoreGlobalVector(dmGrad, &grad);CHKERRQ(ierr);
1615     }
1616   }
1617   /* Handle boundary values */
1618   ierr = DMPlexInsertBoundaryValues(dm, locX, time, faceGeometryFVM, cellGeometryFVM, locGrad);CHKERRQ(ierr);
1619   /* Loop over chunks */
1620   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
1621   numChunks     = 1;
1622   cellChunkSize = (cEnd - cStart)/numChunks;
1623   faceChunkSize = (fEnd - fStart)/numChunks;
1624   for (chunk = 0; chunk < numChunks; ++chunk) {
1625     PetscScalar     *elemVec, *fluxL, *fluxR;
1626     PetscReal       *vol;
1627     PetscFVFaceGeom *fgeom;
1628     PetscInt         cS = cStart+chunk*cellChunkSize, cE = PetscMin(cS+cellChunkSize, cEnd), numCells = cE - cS, cell;
1629     PetscInt         fS = fStart+chunk*faceChunkSize, fE = PetscMin(fS+faceChunkSize, fEnd), numFaces = fE - fS, face;
1630 
1631     /* Extract field coefficients */
1632     if (useFEM) {
1633       ierr = DMPlexGetCellFields(dm, cS, cE, locX, locX_t, locA, &u, &u_t, &a);CHKERRQ(ierr);
1634       ierr = DMGetWorkArray(dm, numCells*totDim, PETSC_SCALAR, &elemVec);CHKERRQ(ierr);
1635       ierr = PetscMemzero(elemVec, numCells*totDim * sizeof(PetscScalar));CHKERRQ(ierr);
1636     }
1637     if (useFVM) {
1638       ierr = DMPlexGetFaceFields(dm, fS, fE, locX, locX_t, faceGeometryFVM, cellGeometryFVM, locGrad, &uL, &uR);CHKERRQ(ierr);
1639       ierr = DMPlexGetFaceGeometry(dm, fS, fE, faceGeometryFVM, cellGeometryFVM, &fgeom, &vol);CHKERRQ(ierr);
1640       ierr = DMGetWorkArray(dm, numFaces*totDim, PETSC_SCALAR, &fluxL);CHKERRQ(ierr);
1641       ierr = DMGetWorkArray(dm, numFaces*totDim, PETSC_SCALAR, &fluxR);CHKERRQ(ierr);
1642       ierr = PetscMemzero(fluxL, numFaces*totDim * sizeof(PetscScalar));CHKERRQ(ierr);
1643       ierr = PetscMemzero(fluxR, numFaces*totDim * sizeof(PetscScalar));CHKERRQ(ierr);
1644     }
1645     /* TODO We will interlace both our field coefficients (u, u_t, uL, uR, etc.) and our output (elemVec, fL, fR). I think this works */
1646     /* Loop over fields */
1647     for (f = 0; f < Nf; ++f) {
1648       PetscObject  obj;
1649       PetscClassId id;
1650       PetscBool    fimp;
1651       PetscInt     numChunks, numBatches, batchSize, numBlocks, blockSize, Ne, Nr, offset;
1652 
1653       ierr = PetscDSGetImplicit(prob, f, &fimp);CHKERRQ(ierr);
1654       if (isImplicit != fimp) continue;
1655       ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
1656       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1657       if (id == PETSCFE_CLASSID) {
1658         PetscFE         fe = (PetscFE) obj;
1659         PetscQuadrature q;
1660         PetscInt        Nq, Nb;
1661 
1662         ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
1663 
1664         ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
1665         ierr = PetscQuadratureGetData(q, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);
1666         ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
1667         blockSize = Nb*Nq;
1668         batchSize = numBlocks * blockSize;
1669         ierr =  PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1670         numChunks = numCells / (numBatches*batchSize);
1671         Ne        = numChunks*numBatches*batchSize;
1672         Nr        = numCells % (numBatches*batchSize);
1673         offset    = numCells - Nr;
1674         /* Integrate FE residual to get elemVec (need fields at quadrature points) */
1675         /*   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) */
1676         ierr = PetscFEIntegrateResidual(fe, prob, f, Ne, cgeomFEM, u, u_t, probAux, a, elemVec);CHKERRQ(ierr);
1677         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);
1678       } else if (id == PETSCFV_CLASSID) {
1679         PetscFV fv = (PetscFV) obj;
1680 
1681         Ne = numFaces;
1682         Nr = 0;
1683         /* Riemann solve over faces (need fields at face centroids) */
1684         /*   We need to evaluate FE fields at those coordinates */
1685         ierr = PetscFVIntegrateRHSFunction(fv, prob, f, Ne, fgeom, vol, uL, uR, fluxL, fluxR);CHKERRQ(ierr);
1686       } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
1687     }
1688     /* Loop over domain */
1689     if (useFEM) {
1690       /* Add elemVec to locX */
1691       for (cell = cS; cell < cE; ++cell) {
1692         if (mesh->printFEM > 1) {ierr = DMPrintCellVector(cell, name, totDim, &elemVec[cell*totDim]);CHKERRQ(ierr);}
1693         ierr = DMPlexVecSetClosure(dm, section, locF, cell, &elemVec[cell*totDim], ADD_VALUES);CHKERRQ(ierr);
1694       }
1695     }
1696     if (useFVM) {
1697       PetscScalar *fa;
1698       PetscInt     iface;
1699 
1700       ierr = VecGetArray(locF, &fa);CHKERRQ(ierr);
1701       for (f = 0; f < Nf; ++f) {
1702         PetscFV      fv;
1703         PetscObject  obj;
1704         PetscClassId id;
1705         PetscInt     foff, pdim;
1706 
1707         ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
1708         ierr = PetscDSGetFieldOffset(prob, f, &foff);CHKERRQ(ierr);
1709         ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1710         if (id != PETSCFV_CLASSID) continue;
1711         fv   = (PetscFV) obj;
1712         ierr = PetscFVGetNumComponents(fv, &pdim);CHKERRQ(ierr);
1713         /* Accumulate fluxes to cells */
1714         for (face = fS, iface = 0; face < fE; ++face) {
1715           const PetscInt *cells;
1716           PetscScalar    *fL, *fR;
1717           PetscInt        ghost, d;
1718 
1719           ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr);
1720           if (ghost >= 0) continue;
1721           ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr);
1722           ierr = DMPlexPointGlobalFieldRef(dm, cells[0], f, fa, &fL);CHKERRQ(ierr);
1723           ierr = DMPlexPointGlobalFieldRef(dm, cells[1], f, fa, &fR);CHKERRQ(ierr);
1724           for (d = 0; d < pdim; ++d) {
1725             if (fL) fL[d] -= fluxL[iface*totDim+foff+d];
1726             if (fR) fR[d] += fluxR[iface*totDim+foff+d];
1727           }
1728           ++iface;
1729         }
1730       }
1731       ierr = VecRestoreArray(locF, &fa);CHKERRQ(ierr);
1732     }
1733     /* Handle time derivative */
1734     if (locX_t) {
1735       PetscScalar *x_t, *fa;
1736 
1737       ierr = VecGetArray(locF, &fa);CHKERRQ(ierr);
1738       ierr = VecGetArray(locX_t, &x_t);CHKERRQ(ierr);
1739       for (f = 0; f < Nf; ++f) {
1740         PetscFV      fv;
1741         PetscObject  obj;
1742         PetscClassId id;
1743         PetscInt     pdim, d;
1744 
1745         ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
1746         ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1747         if (id != PETSCFV_CLASSID) continue;
1748         fv   = (PetscFV) obj;
1749         ierr = PetscFVGetNumComponents(fv, &pdim);CHKERRQ(ierr);
1750         for (cell = cS; cell < cE; ++cell) {
1751           PetscScalar *u_t, *r;
1752 
1753           ierr = DMPlexPointLocalFieldRead(dm, cell, f, x_t, &u_t);CHKERRQ(ierr);
1754           ierr = DMPlexPointLocalFieldRef(dm, cell, f, fa, &r);CHKERRQ(ierr);
1755           for (d = 0; d < pdim; ++d) r[d] += u_t[d];
1756         }
1757       }
1758       ierr = VecRestoreArray(locX_t, &x_t);CHKERRQ(ierr);
1759       ierr = VecRestoreArray(locF, &fa);CHKERRQ(ierr);
1760     }
1761     if (useFEM) {
1762       ierr = DMPlexRestoreCellFields(dm, cS, cE, locX, locX_t, locA, &u, &u_t, &a);CHKERRQ(ierr);
1763       ierr = DMRestoreWorkArray(dm, numCells*totDim, PETSC_SCALAR, &elemVec);CHKERRQ(ierr);
1764     }
1765     if (useFVM) {
1766       ierr = DMPlexRestoreFaceFields(dm, fS, fE, locX, locX_t, faceGeometryFVM, cellGeometryFVM, locGrad, &uL, &uR);CHKERRQ(ierr);
1767       ierr = DMPlexRestoreFaceGeometry(dm, fS, fE, faceGeometryFVM, cellGeometryFVM, &fgeom, &vol);CHKERRQ(ierr);
1768       ierr = DMRestoreWorkArray(dm, numFaces*totDim, PETSC_SCALAR, &fluxL);CHKERRQ(ierr);
1769       ierr = DMRestoreWorkArray(dm, numFaces*totDim, PETSC_SCALAR, &fluxR);CHKERRQ(ierr);
1770       if (dmGrad) {ierr = DMRestoreLocalVector(dmGrad, &locGrad);CHKERRQ(ierr);}
1771     }
1772   }
1773 
1774   if (useFEM) {ierr = DMPlexComputeBdResidual_Internal(dm, locX, locX_t, locF, user);CHKERRQ(ierr);}
1775 
1776   /* FEM */
1777   /* 1: Get sizes from dm and dmAux */
1778   /* 2: Get geometric data */
1779   /* 3: Handle boundary values */
1780   /* 4: Loop over domain */
1781   /*   Extract coefficients */
1782   /* Loop over fields */
1783   /*   Set tiling for FE*/
1784   /*   Integrate FE residual to get elemVec */
1785   /*     Loop over subdomain */
1786   /*       Loop over quad points */
1787   /*         Transform coords to real space */
1788   /*         Evaluate field and aux fields at point */
1789   /*         Evaluate residual at point */
1790   /*         Transform residual to real space */
1791   /*       Add residual to elemVec */
1792   /* Loop over domain */
1793   /*   Add elemVec to locX */
1794 
1795   /* FVM */
1796   /* Get geometric data */
1797   /* If using gradients */
1798   /*   Compute gradient data */
1799   /*   Loop over domain faces */
1800   /*     Count computational faces */
1801   /*     Reconstruct cell gradient */
1802   /*   Loop over domain cells */
1803   /*     Limit cell gradients */
1804   /* Handle boundary values */
1805   /* Loop over domain faces */
1806   /*   Read out field, centroid, normal, volume for each side of face */
1807   /* Riemann solve over faces */
1808   /* Loop over domain faces */
1809   /*   Accumulate fluxes to cells */
1810   /* TODO Change printFEM to printDisc here */
1811   if (mesh->printFEM) {ierr = DMPrintLocalVec(dm, name, mesh->printTol, locF);CHKERRQ(ierr);}
1812   ierr = PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
1813   PetscFunctionReturn(0);
1814 }
1815 
1816 #undef __FUNCT__
1817 #define __FUNCT__ "DMPlexComputeResidualFEM_Check_Internal"
1818 static PetscErrorCode DMPlexComputeResidualFEM_Check_Internal(DM dm, Vec X, Vec X_t, Vec F, void *user)
1819 {
1820   DM                dmCh, dmAux;
1821   Vec               A, cellgeom;
1822   PetscDS           prob, probCh, probAux = NULL;
1823   PetscQuadrature   q;
1824   PetscSection      section, sectionAux;
1825   PetscFECellGeom  *cgeom;
1826   PetscScalar      *elemVec, *elemVecCh, *u, *u_t, *a = NULL;
1827   PetscInt          dim, Nf, f, numCells, cStart, cEnd, c;
1828   PetscInt          totDim, totDimAux, diffCell = 0;
1829   PetscErrorCode    ierr;
1830 
1831   PetscFunctionBegin;
1832   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1833   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1834   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
1835   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
1836   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
1837   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1838   numCells = cEnd - cStart;
1839   ierr = PetscObjectQuery((PetscObject) dm, "dmCh", (PetscObject *) &dmCh);CHKERRQ(ierr);
1840   ierr = DMGetDS(dmCh, &probCh);CHKERRQ(ierr);
1841   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
1842   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
1843   if (dmAux) {
1844     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
1845     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
1846     ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
1847   }
1848   ierr = DMPlexInsertBoundaryValues(dm, X, 0.0, NULL, NULL, NULL);CHKERRQ(ierr);
1849   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
1850   ierr = PetscMalloc3(numCells*totDim,&u,X_t ? numCells*totDim : 0,&u_t,numCells*totDim,&elemVec);CHKERRQ(ierr);
1851   ierr = PetscMalloc1(numCells*totDim,&elemVecCh);CHKERRQ(ierr);
1852   if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);}
1853   ierr = DMPlexSNESGetGeometryFEM(dm, &cellgeom);CHKERRQ(ierr);
1854   ierr = VecGetArray(cellgeom, (PetscScalar **) &cgeom);CHKERRQ(ierr);
1855   for (c = cStart; c < cEnd; ++c) {
1856     PetscScalar *x = NULL, *x_t = NULL;
1857     PetscInt     i;
1858 
1859     ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
1860     for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i];
1861     ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
1862     if (X_t) {
1863       ierr = DMPlexVecGetClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr);
1864       for (i = 0; i < totDim; ++i) u_t[c*totDim+i] = x_t[i];
1865       ierr = DMPlexVecRestoreClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr);
1866     }
1867     if (dmAux) {
1868       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
1869       for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i];
1870       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
1871     }
1872   }
1873   for (f = 0; f < Nf; ++f) {
1874     PetscFE  fe, feCh;
1875     PetscInt numQuadPoints, Nb;
1876     /* Conforming batches */
1877     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1878     /* Remainder */
1879     PetscInt Nr, offset;
1880 
1881     ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
1882     ierr = PetscDSGetDiscretization(probCh, f, (PetscObject *) &feCh);CHKERRQ(ierr);
1883     ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
1884     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
1885     ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
1886     ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
1887     blockSize = Nb*numQuadPoints;
1888     batchSize = numBlocks * blockSize;
1889     ierr =  PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1890     numChunks = numCells / (numBatches*batchSize);
1891     Ne        = numChunks*numBatches*batchSize;
1892     Nr        = numCells % (numBatches*batchSize);
1893     offset    = numCells - Nr;
1894     ierr = PetscFEIntegrateResidual(fe, prob, f, Ne, cgeom, u, u_t, probAux, a, elemVec);CHKERRQ(ierr);
1895     ierr = PetscFEIntegrateResidual(feCh, prob, f, Ne, cgeom, u, u_t, probAux, a, elemVecCh);CHKERRQ(ierr);
1896     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);
1897     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);
1898   }
1899   for (c = cStart; c < cEnd; ++c) {
1900     PetscBool diff = PETSC_FALSE;
1901     PetscInt  d;
1902 
1903     for (d = 0; d < totDim; ++d) if (PetscAbsScalar(elemVec[c*totDim+d] - elemVecCh[c*totDim+d]) > 1.0e-7) {diff = PETSC_TRUE;break;}
1904     if (diff) {
1905       ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Different cell %d\n", c);CHKERRQ(ierr);
1906       ierr = DMPrintCellVector(c, "Residual", totDim, &elemVec[c*totDim]);CHKERRQ(ierr);
1907       ierr = DMPrintCellVector(c, "Check Residual", totDim, &elemVecCh[c*totDim]);CHKERRQ(ierr);
1908       ++diffCell;
1909     }
1910     if (diffCell > 9) break;
1911     ierr = DMPlexVecSetClosure(dm, section, F, c, &elemVec[c*totDim], ADD_VALUES);CHKERRQ(ierr);
1912   }
1913   ierr = VecRestoreArray(cellgeom, (PetscScalar **) &cgeom);CHKERRQ(ierr);
1914   ierr = PetscFree3(u,u_t,elemVec);CHKERRQ(ierr);
1915   ierr = PetscFree(elemVecCh);CHKERRQ(ierr);
1916   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
1917   PetscFunctionReturn(0);
1918 }
1919 
1920 #undef __FUNCT__
1921 #define __FUNCT__ "DMPlexSNESComputeResidualFEM"
1922 /*@
1923   DMPlexSNESComputeResidualFEM - Form the local residual F from the local input X using pointwise functions specified by the user
1924 
1925   Input Parameters:
1926 + dm - The mesh
1927 . X  - Local solution
1928 - user - The user context
1929 
1930   Output Parameter:
1931 . F  - Local output vector
1932 
1933   Level: developer
1934 
1935 .seealso: DMPlexComputeJacobianActionFEM()
1936 @*/
1937 PetscErrorCode DMPlexSNESComputeResidualFEM(DM dm, Vec X, Vec F, void *user)
1938 {
1939   PetscObject    check;
1940   PetscInt       cStart, cEnd, cEndInterior;
1941   PetscErrorCode ierr;
1942 
1943   PetscFunctionBegin;
1944   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1945   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
1946   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
1947   /* The dmCh is used to check two mathematically equivalent discretizations for computational equivalence */
1948   ierr = PetscObjectQuery((PetscObject) dm, "dmCh", &check);CHKERRQ(ierr);
1949   if (check) {ierr = DMPlexComputeResidualFEM_Check_Internal(dm, X, NULL, F, user);CHKERRQ(ierr);}
1950   else       {ierr = DMPlexComputeResidual_Internal(dm, cStart, cEnd, PETSC_MIN_REAL, X, NULL, F, user);CHKERRQ(ierr);}
1951   PetscFunctionReturn(0);
1952 }
1953 
1954 #undef __FUNCT__
1955 #define __FUNCT__ "DMPlexComputeJacobian_Internal"
1956 PetscErrorCode DMPlexComputeJacobian_Internal(DM dm, PetscInt cStart, PetscInt cEnd, PetscReal t, PetscReal X_tShift, Vec X, Vec X_t, Mat Jac, Mat JacP,void *user)
1957 {
1958   DM_Plex          *mesh  = (DM_Plex *) dm->data;
1959   const char       *name  = "Jacobian";
1960   DM                dmAux;
1961   DMLabel           depth;
1962   Vec               A, cellgeom;
1963   PetscDS           prob, probAux = NULL;
1964   PetscQuadrature   quad;
1965   PetscSection      section, globalSection, sectionAux;
1966   PetscFECellGeom  *cgeom;
1967   PetscScalar      *elemMat, *u, *u_t, *a = NULL;
1968   PetscInt          dim, Nf, f, fieldI, fieldJ, numCells, c;
1969   PetscInt          totDim, totDimBd, totDimAux, numBd, bd;
1970   PetscBool         isShell;
1971   PetscErrorCode    ierr;
1972 
1973   PetscFunctionBegin;
1974   ierr = PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
1975   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1976   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1977   ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);
1978   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
1979   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
1980   ierr = PetscDSGetTotalBdDimension(prob, &totDimBd);CHKERRQ(ierr);
1981   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
1982   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1983   numCells = cEnd - cStart;
1984   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
1985   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
1986   if (dmAux) {
1987     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
1988     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
1989     ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
1990   }
1991   ierr = DMPlexInsertBoundaryValues(dm, X, 0.0, NULL, NULL, NULL);CHKERRQ(ierr);
1992   ierr = MatZeroEntries(JacP);CHKERRQ(ierr);
1993   ierr = PetscMalloc3(numCells*totDim,&u,X_t ? numCells*totDim : 0,&u_t,numCells*totDim*totDim,&elemMat);CHKERRQ(ierr);
1994   if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);}
1995   ierr = DMPlexSNESGetGeometryFEM(dm, &cellgeom);CHKERRQ(ierr);
1996   ierr = VecGetArray(cellgeom, (PetscScalar **) &cgeom);CHKERRQ(ierr);
1997   for (c = cStart; c < cEnd; ++c) {
1998     PetscScalar *x = NULL,  *x_t = NULL;
1999     PetscInt     i;
2000 
2001     ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
2002     for (i = 0; i < totDim; ++i) u[(c-cStart)*totDim+i] = x[i];
2003     ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
2004     if (X_t) {
2005       ierr = DMPlexVecGetClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr);
2006       for (i = 0; i < totDim; ++i) u_t[(c-cStart)*totDim+i] = x_t[i];
2007       ierr = DMPlexVecRestoreClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr);
2008     }
2009     if (dmAux) {
2010       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
2011       for (i = 0; i < totDimAux; ++i) a[(c-cStart)*totDimAux+i] = x[i];
2012       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
2013     }
2014   }
2015   ierr = PetscMemzero(elemMat, numCells*totDim*totDim * sizeof(PetscScalar));CHKERRQ(ierr);
2016   for (fieldI = 0; fieldI < Nf; ++fieldI) {
2017     PetscFE  fe;
2018     PetscInt numQuadPoints, Nb;
2019     /* Conforming batches */
2020     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
2021     /* Remainder */
2022     PetscInt Nr, offset;
2023 
2024     ierr = PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);CHKERRQ(ierr);
2025     ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
2026     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
2027     ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
2028     ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
2029     blockSize = Nb*numQuadPoints;
2030     batchSize = numBlocks * blockSize;
2031     ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
2032     numChunks = numCells / (numBatches*batchSize);
2033     Ne        = numChunks*numBatches*batchSize;
2034     Nr        = numCells % (numBatches*batchSize);
2035     offset    = numCells - Nr;
2036     for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
2037       ierr = PetscFEIntegrateJacobian(fe, prob, fieldI, fieldJ, Ne, cgeom, u, u_t, probAux, a, elemMat);CHKERRQ(ierr);
2038       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);
2039     }
2040   }
2041   for (c = cStart; c < cEnd; ++c) {
2042     if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, name, totDim, totDim, &elemMat[(c-cStart)*totDim*totDim]);CHKERRQ(ierr);}
2043     ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, c, &elemMat[(c-cStart)*totDim*totDim], ADD_VALUES);CHKERRQ(ierr);
2044   }
2045   ierr = VecGetArray(cellgeom, (PetscScalar **) &cgeom);CHKERRQ(ierr);
2046   ierr = PetscFree3(u,u_t,elemMat);CHKERRQ(ierr);
2047   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
2048   ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr);
2049   ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
2050   ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr);
2051   ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
2052   for (bd = 0; bd < numBd; ++bd) {
2053     const char     *bdLabel;
2054     DMLabel         label;
2055     IS              pointIS;
2056     const PetscInt *points;
2057     const PetscInt *values;
2058     PetscInt        field, numValues, numPoints, p, dep, numFaces;
2059     PetscBool       isEssential;
2060 
2061     ierr = DMPlexGetBoundary(dm, bd, &isEssential, NULL, &bdLabel, &field, NULL, NULL, NULL, &numValues, &values, NULL);CHKERRQ(ierr);
2062     if (isEssential) continue;
2063     if (numValues != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Bug me and I will fix this");
2064     ierr = DMPlexGetLabel(dm, bdLabel, &label);CHKERRQ(ierr);
2065     ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr);
2066     ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr);
2067     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
2068     for (p = 0, numFaces = 0; p < numPoints; ++p) {
2069       ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
2070       if (dep == dim-1) ++numFaces;
2071     }
2072     ierr = PetscMalloc3(numFaces*totDimBd,&u,numFaces,&cgeom,numFaces*totDimBd*totDimBd,&elemMat);CHKERRQ(ierr);
2073     if (X_t) {ierr = PetscMalloc1(numFaces*totDimBd,&u_t);CHKERRQ(ierr);}
2074     for (p = 0, f = 0; p < numPoints; ++p) {
2075       const PetscInt point = points[p];
2076       PetscScalar   *x     = NULL;
2077       PetscInt       i;
2078 
2079       ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
2080       if (dep != dim-1) continue;
2081       ierr = DMPlexComputeCellGeometryFEM(dm, point, NULL, cgeom[f].v0, cgeom[f].J, cgeom[f].invJ, &cgeom[f].detJ);CHKERRQ(ierr);
2082       ierr = DMPlexComputeCellGeometryFVM(dm, point, NULL, NULL, cgeom[f].n);CHKERRQ(ierr);
2083       if (cgeom[f].detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", cgeom[f].detJ, point);
2084       ierr = DMPlexVecGetClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
2085       for (i = 0; i < totDimBd; ++i) u[f*totDimBd+i] = x[i];
2086       ierr = DMPlexVecRestoreClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
2087       if (X_t) {
2088         ierr = DMPlexVecGetClosure(dm, section, X_t, point, NULL, &x);CHKERRQ(ierr);
2089         for (i = 0; i < totDimBd; ++i) u_t[f*totDimBd+i] = x[i];
2090         ierr = DMPlexVecRestoreClosure(dm, section, X_t, point, NULL, &x);CHKERRQ(ierr);
2091       }
2092       ++f;
2093     }
2094     ierr = PetscMemzero(elemMat, numFaces*totDimBd*totDimBd * sizeof(PetscScalar));CHKERRQ(ierr);
2095     for (fieldI = 0; fieldI < Nf; ++fieldI) {
2096       PetscFE  fe;
2097       PetscInt numQuadPoints, Nb;
2098       /* Conforming batches */
2099       PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
2100       /* Remainder */
2101       PetscInt Nr, offset;
2102 
2103       ierr = PetscDSGetBdDiscretization(prob, fieldI, (PetscObject *) &fe);CHKERRQ(ierr);
2104       ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
2105       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
2106       ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
2107       ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
2108       blockSize = Nb*numQuadPoints;
2109       batchSize = numBlocks * blockSize;
2110       ierr =  PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
2111       numChunks = numFaces / (numBatches*batchSize);
2112       Ne        = numChunks*numBatches*batchSize;
2113       Nr        = numFaces % (numBatches*batchSize);
2114       offset    = numFaces - Nr;
2115       for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
2116         ierr = PetscFEIntegrateBdJacobian(fe, prob, fieldI, fieldJ, Ne, cgeom, u, u_t, NULL, NULL, elemMat);CHKERRQ(ierr);
2117         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);
2118       }
2119     }
2120     for (p = 0, f = 0; p < numPoints; ++p) {
2121       const PetscInt point = points[p];
2122 
2123       ierr = DMLabelGetValue(depth, point, &dep);CHKERRQ(ierr);
2124       if (dep != dim-1) continue;
2125       if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(point, "BdJacobian", totDimBd, totDimBd, &elemMat[f*totDimBd*totDimBd]);CHKERRQ(ierr);}
2126       ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, point, &elemMat[f*totDimBd*totDimBd], ADD_VALUES);CHKERRQ(ierr);
2127       ++f;
2128     }
2129     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
2130     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
2131     ierr = PetscFree3(u,cgeom,elemMat);CHKERRQ(ierr);
2132     if (X_t) {ierr = PetscFree(u_t);CHKERRQ(ierr);}
2133   }
2134   ierr = MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2135   ierr = MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2136   if (mesh->printFEM) {
2137     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
2138     ierr = MatChop(JacP, 1.0e-10);CHKERRQ(ierr);
2139     ierr = MatView(JacP, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
2140   }
2141   ierr = PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
2142   ierr = PetscObjectTypeCompare((PetscObject) Jac, MATSHELL, &isShell);CHKERRQ(ierr);
2143   if (isShell) {
2144     JacActionCtx *jctx;
2145 
2146     ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr);
2147     ierr = VecCopy(X, jctx->u);CHKERRQ(ierr);
2148   }
2149   PetscFunctionReturn(0);
2150 }
2151 
2152 #undef __FUNCT__
2153 #define __FUNCT__ "DMPlexSNESComputeJacobianFEM"
2154 /*@
2155   DMPlexSNESComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user.
2156 
2157   Input Parameters:
2158 + dm - The mesh
2159 . X  - Local input vector
2160 - user - The user context
2161 
2162   Output Parameter:
2163 . Jac  - Jacobian matrix
2164 
2165   Note:
2166   The first member of the user context must be an FEMContext.
2167 
2168   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
2169   like a GPU, or vectorize on a multicore machine.
2170 
2171   Level: developer
2172 
2173 .seealso: FormFunctionLocal()
2174 @*/
2175 PetscErrorCode DMPlexSNESComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP,void *user)
2176 {
2177   PetscInt       cStart, cEnd, cEndInterior;
2178   PetscErrorCode ierr;
2179 
2180   PetscFunctionBegin;
2181   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
2182   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
2183   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
2184   ierr = DMPlexComputeJacobian_Internal(dm, cStart, cEnd, 0.0, 0.0, X, NULL, Jac, JacP, user);CHKERRQ(ierr);
2185   PetscFunctionReturn(0);
2186 }
2187 
2188 #undef __FUNCT__
2189 #define __FUNCT__ "DMSNESCheckFromOptions_Internal"
2190 PetscErrorCode DMSNESCheckFromOptions_Internal(SNES snes, DM dm, Vec u, Vec sol, PetscErrorCode (**exactFuncs)(PetscInt, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx), void **ctxs)
2191 {
2192   Mat            J, M;
2193   Vec            r, b;
2194   MatNullSpace   nullSpace;
2195   PetscReal     *error, res = 0.0;
2196   PetscInt       numFields;
2197   PetscErrorCode ierr;
2198 
2199   PetscFunctionBegin;
2200   ierr = VecDuplicate(u, &r);CHKERRQ(ierr);
2201   ierr = DMCreateMatrix(dm, &J);CHKERRQ(ierr);
2202   M    = J;
2203   /* TODO Null space for J */
2204   /* Check discretization error */
2205   ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr);
2206   ierr = PetscMalloc1(PetscMax(1, numFields), &error);CHKERRQ(ierr);
2207   if (numFields > 1) {
2208     PetscInt f;
2209 
2210     ierr = DMPlexComputeL2FieldDiff(dm, exactFuncs, ctxs, u, error);CHKERRQ(ierr);
2211     ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: [");CHKERRQ(ierr);
2212     for (f = 0; f < numFields; ++f) {
2213       if (f) {ierr = PetscPrintf(PETSC_COMM_WORLD, ", ");CHKERRQ(ierr);}
2214       if (error[f] >= 1.0e-11) {ierr = PetscPrintf(PETSC_COMM_WORLD, "%g", error[f]);CHKERRQ(ierr);}
2215       else                     {ierr = PetscPrintf(PETSC_COMM_WORLD, "< 1.0e-11");CHKERRQ(ierr);}
2216     }
2217     ierr = PetscPrintf(PETSC_COMM_WORLD, "]\n");CHKERRQ(ierr);
2218   } else {
2219     ierr = DMPlexComputeL2Diff(dm, exactFuncs, ctxs, u, &error[0]);CHKERRQ(ierr);
2220     if (error[0] >= 1.0e-11) {ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: %g\n", error[0]);CHKERRQ(ierr);}
2221     else                     {ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: < 1.0e-11\n");CHKERRQ(ierr);}
2222   }
2223   ierr = PetscFree(error);CHKERRQ(ierr);
2224   /* Check residual */
2225   ierr = SNESComputeFunction(snes, u, r);CHKERRQ(ierr);
2226   ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr);
2227   ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Residual: %g\n", res);CHKERRQ(ierr);
2228   ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr);
2229   ierr = PetscObjectSetName((PetscObject) r, "Initial Residual");CHKERRQ(ierr);
2230   ierr = PetscObjectSetOptionsPrefix((PetscObject)r,"res_");CHKERRQ(ierr);
2231   ierr = VecViewFromOptions(r, NULL, "-vec_view");CHKERRQ(ierr);
2232   /* Check Jacobian */
2233   ierr = SNESComputeJacobian(snes, u, M, M);CHKERRQ(ierr);
2234   ierr = MatGetNullSpace(J, &nullSpace);CHKERRQ(ierr);
2235   if (nullSpace) {
2236     PetscBool isNull;
2237     ierr = MatNullSpaceTest(nullSpace, J, &isNull);CHKERRQ(ierr);
2238     if (!isNull) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "The null space calculated for the system operator is invalid.");
2239   }
2240   ierr = VecDuplicate(u, &b);CHKERRQ(ierr);
2241   ierr = VecSet(r, 0.0);CHKERRQ(ierr);
2242   ierr = SNESComputeFunction(snes, r, b);CHKERRQ(ierr);
2243   ierr = MatMult(M, u, r);CHKERRQ(ierr);
2244   ierr = VecAXPY(r, 1.0, b);CHKERRQ(ierr);
2245   ierr = VecDestroy(&b);CHKERRQ(ierr);
2246   ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr);
2247   ierr = PetscPrintf(PETSC_COMM_WORLD, "Linear L_2 Residual: %g\n", res);CHKERRQ(ierr);
2248   ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr);
2249   ierr = PetscObjectSetName((PetscObject) r, "Au - b = Au + F(0)");CHKERRQ(ierr);
2250   ierr = PetscObjectSetOptionsPrefix((PetscObject)r,"linear_res_");CHKERRQ(ierr);
2251   ierr = VecViewFromOptions(r, NULL, "-vec_view");CHKERRQ(ierr);
2252   ierr = VecDestroy(&r);CHKERRQ(ierr);
2253   ierr = MatNullSpaceDestroy(&nullSpace);CHKERRQ(ierr);
2254   ierr = MatDestroy(&J);CHKERRQ(ierr);
2255   PetscFunctionReturn(0);
2256 }
2257 
2258 #undef __FUNCT__
2259 #define __FUNCT__ "DMSNESCheckFromOptions"
2260 PetscErrorCode DMSNESCheckFromOptions(SNES snes, Vec u, PetscErrorCode (**exactFuncs)(PetscInt dim, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx), void **ctxs)
2261 {
2262   DM             dm;
2263   Vec            sol;
2264   PetscBool      check;
2265   PetscErrorCode ierr;
2266 
2267   PetscFunctionBegin;
2268   ierr = PetscOptionsHasName(snes->hdr.prefix, "-dmsnes_check", &check);CHKERRQ(ierr);
2269   if (!check) PetscFunctionReturn(0);
2270   ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr);
2271   ierr = VecDuplicate(u, &sol);CHKERRQ(ierr);
2272   ierr = SNESSetSolution(snes, sol);CHKERRQ(ierr);
2273   ierr = DMSNESCheckFromOptions_Internal(snes, dm, u, sol, exactFuncs, ctxs);CHKERRQ(ierr);
2274   ierr = VecDestroy(&sol);CHKERRQ(ierr);
2275   PetscFunctionReturn(0);
2276 }
2277