xref: /petsc/src/dm/impls/plex/plexfem.c (revision 249cb5c497384507c8233e20189d5c9ded68a4fc)
1 #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2 #include <petscsf.h>
3 
4 #include <petsc/private/petscfeimpl.h>
5 #include <petsc/private/petscfvimpl.h>
6 
7 /*@
8   DMPlexGetScale - Get the scale for the specified fundamental unit
9 
10   Not collective
11 
12   Input Arguments:
13 + dm   - the DM
14 - unit - The SI unit
15 
16   Output Argument:
17 . scale - The value used to scale all quantities with this unit
18 
19   Level: advanced
20 
21 .seealso: DMPlexSetScale(), PetscUnit
22 @*/
23 PetscErrorCode DMPlexGetScale(DM dm, PetscUnit unit, PetscReal *scale)
24 {
25   DM_Plex *mesh = (DM_Plex*) dm->data;
26 
27   PetscFunctionBegin;
28   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
29   PetscValidPointer(scale, 3);
30   *scale = mesh->scale[unit];
31   PetscFunctionReturn(0);
32 }
33 
34 /*@
35   DMPlexSetScale - Set the scale for the specified fundamental unit
36 
37   Not collective
38 
39   Input Arguments:
40 + dm   - the DM
41 . unit - The SI unit
42 - scale - The value used to scale all quantities with this unit
43 
44   Level: advanced
45 
46 .seealso: DMPlexGetScale(), PetscUnit
47 @*/
48 PetscErrorCode DMPlexSetScale(DM dm, PetscUnit unit, PetscReal scale)
49 {
50   DM_Plex *mesh = (DM_Plex*) dm->data;
51 
52   PetscFunctionBegin;
53   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
54   mesh->scale[unit] = scale;
55   PetscFunctionReturn(0);
56 }
57 
58 static PetscErrorCode DMPlexProjectRigidBody_Private(PetscInt dim, PetscReal t, const PetscReal X[], PetscInt Nf, PetscScalar *mode, void *ctx)
59 {
60   const PetscInt eps[3][3][3] = {{{0, 0, 0}, {0, 0, 1}, {0, -1, 0}}, {{0, 0, -1}, {0, 0, 0}, {1, 0, 0}}, {{0, 1, 0}, {-1, 0, 0}, {0, 0, 0}}};
61   PetscInt *ctxInt  = (PetscInt *) ctx;
62   PetscInt  dim2    = ctxInt[0];
63   PetscInt  d       = ctxInt[1];
64   PetscInt  i, j, k = dim > 2 ? d - dim : d;
65 
66   PetscFunctionBegin;
67   if (dim != dim2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Input dimension %d does not match context dimension %d", dim, dim2);
68   for (i = 0; i < dim; i++) mode[i] = 0.;
69   if (d < dim) {
70     mode[d] = 1.; /* Translation along axis d */
71   } else {
72     for (i = 0; i < dim; i++) {
73       for (j = 0; j < dim; j++) {
74         mode[j] += eps[i][j][k]*X[i]; /* Rotation about axis d */
75       }
76     }
77   }
78   PetscFunctionReturn(0);
79 }
80 
81 /*@C
82   DMPlexCreateRigidBody - For the default global section, create rigid body modes by function space interpolation
83 
84   Collective on DM
85 
86   Input Arguments:
87 . dm - the DM
88 
89   Output Argument:
90 . sp - the null space
91 
92   Note: This is necessary to provide a suitable coarse space for algebraic multigrid
93 
94   Level: advanced
95 
96 .seealso: MatNullSpaceCreate(), PCGAMG
97 @*/
98 PetscErrorCode DMPlexCreateRigidBody(DM dm, MatNullSpace *sp)
99 {
100   MPI_Comm       comm;
101   Vec            mode[6];
102   PetscSection   section, globalSection;
103   PetscInt       dim, dimEmbed, n, m, d, i, j;
104   PetscErrorCode ierr;
105 
106   PetscFunctionBegin;
107   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
108   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
109   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
110   if (dim == 1) {
111     ierr = MatNullSpaceCreate(comm, PETSC_TRUE, 0, NULL, sp);CHKERRQ(ierr);
112     PetscFunctionReturn(0);
113   }
114   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
115   ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);
116   ierr = PetscSectionGetConstrainedStorageSize(globalSection, &n);CHKERRQ(ierr);
117   m    = (dim*(dim+1))/2;
118   ierr = VecCreate(comm, &mode[0]);CHKERRQ(ierr);
119   ierr = VecSetSizes(mode[0], n, PETSC_DETERMINE);CHKERRQ(ierr);
120   ierr = VecSetUp(mode[0]);CHKERRQ(ierr);
121   for (i = 1; i < m; ++i) {ierr = VecDuplicate(mode[0], &mode[i]);CHKERRQ(ierr);}
122   for (d = 0; d < m; d++) {
123     PetscInt         ctx[2];
124     PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal *, PetscInt, PetscScalar *, void *) = DMPlexProjectRigidBody_Private;
125     void            *voidctx = (void *) (&ctx[0]);
126 
127     ctx[0] = dimEmbed;
128     ctx[1] = d;
129     ierr = DMProjectFunction(dm, 0.0, &func, &voidctx, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
130   }
131   for (i = 0; i < dim; ++i) {ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);}
132   /* Orthonormalize system */
133   for (i = dim; i < m; ++i) {
134     PetscScalar dots[6];
135 
136     ierr = VecMDot(mode[i], i, mode, dots);CHKERRQ(ierr);
137     for (j = 0; j < i; ++j) dots[j] *= -1.0;
138     ierr = VecMAXPY(mode[i], i, dots, mode);CHKERRQ(ierr);
139     ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);
140   }
141   ierr = MatNullSpaceCreate(comm, PETSC_FALSE, m, mode, sp);CHKERRQ(ierr);
142   for (i = 0; i< m; ++i) {ierr = VecDestroy(&mode[i]);CHKERRQ(ierr);}
143   PetscFunctionReturn(0);
144 }
145 
146 /*@C
147   DMPlexCreateRigidBodies - For the default global section, create rigid body modes by function space interpolation
148 
149   Collective on DM
150 
151   Input Arguments:
152 + dm    - the DM
153 . nb    - The number of bodies
154 . label - The DMLabel marking each domain
155 . nids  - The number of ids per body
156 - ids   - An array of the label ids in sequence for each domain
157 
158   Output Argument:
159 . sp - the null space
160 
161   Note: This is necessary to provide a suitable coarse space for algebraic multigrid
162 
163   Level: advanced
164 
165 .seealso: MatNullSpaceCreate()
166 @*/
167 PetscErrorCode DMPlexCreateRigidBodies(DM dm, PetscInt nb, DMLabel label, const PetscInt nids[], const PetscInt ids[], MatNullSpace *sp)
168 {
169   MPI_Comm       comm;
170   PetscSection   section, globalSection;
171   Vec           *mode;
172   PetscScalar   *dots;
173   PetscInt       dim, dimEmbed, n, m, b, d, i, j, off;
174   PetscErrorCode ierr;
175 
176   PetscFunctionBegin;
177   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
178   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
179   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
180   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
181   ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);
182   ierr = PetscSectionGetConstrainedStorageSize(globalSection, &n);CHKERRQ(ierr);
183   m    = nb * (dim*(dim+1))/2;
184   ierr = PetscMalloc2(m, &mode, m, &dots);CHKERRQ(ierr);
185   ierr = VecCreate(comm, &mode[0]);CHKERRQ(ierr);
186   ierr = VecSetSizes(mode[0], n, PETSC_DETERMINE);CHKERRQ(ierr);
187   ierr = VecSetUp(mode[0]);CHKERRQ(ierr);
188   for (i = 1; i < m; ++i) {ierr = VecDuplicate(mode[0], &mode[i]);CHKERRQ(ierr);}
189   for (b = 0, off = 0; b < nb; ++b) {
190     for (d = 0; d < m/nb; ++d) {
191       PetscInt         ctx[2];
192       PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal *, PetscInt, PetscScalar *, void *) = DMPlexProjectRigidBody_Private;
193       void            *voidctx = (void *) (&ctx[0]);
194 
195       ctx[0] = dimEmbed;
196       ctx[1] = d;
197       ierr = DMProjectFunctionLabel(dm, 0.0, label, nids[b], &ids[off], 0, NULL, &func, &voidctx, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
198       off   += nids[b];
199     }
200   }
201   for (i = 0; i < dim; ++i) {ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);}
202   /* Orthonormalize system */
203   for (i = 0; i < m; ++i) {
204     ierr = VecMDot(mode[i], i, mode, dots);CHKERRQ(ierr);
205     for (j = 0; j < i; ++j) dots[j] *= -1.0;
206     ierr = VecMAXPY(mode[i], i, dots, mode);CHKERRQ(ierr);
207     ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);
208   }
209   ierr = MatNullSpaceCreate(comm, PETSC_FALSE, m, mode, sp);CHKERRQ(ierr);
210   for (i = 0; i< m; ++i) {ierr = VecDestroy(&mode[i]);CHKERRQ(ierr);}
211   ierr = PetscFree2(mode, dots);CHKERRQ(ierr);
212   PetscFunctionReturn(0);
213 }
214 
215 /*@
216   DMPlexSetMaxProjectionHeight - In DMPlexProjectXXXLocal() functions, the projected values of a basis function's dofs
217   are computed by associating the basis function with one of the mesh points in its transitively-closed support, and
218   evaluating the dual space basis of that point.  A basis function is associated with the point in its
219   transitively-closed support whose mesh height is highest (w.r.t. DAG height), but not greater than the maximum
220   projection height, which is set with this function.  By default, the maximum projection height is zero, which means
221   that only mesh cells are used to project basis functions.  A height of one, for example, evaluates a cell-interior
222   basis functions using its cells dual space basis, but all other basis functions with the dual space basis of a face.
223 
224   Input Parameters:
225 + dm - the DMPlex object
226 - height - the maximum projection height >= 0
227 
228   Level: advanced
229 
230 .seealso: DMPlexGetMaxProjectionHeight(), DMProjectFunctionLocal(), DMProjectFunctionLabelLocal()
231 @*/
232 PetscErrorCode DMPlexSetMaxProjectionHeight(DM dm, PetscInt height)
233 {
234   DM_Plex *plex = (DM_Plex *) dm->data;
235 
236   PetscFunctionBegin;
237   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
238   plex->maxProjectionHeight = height;
239   PetscFunctionReturn(0);
240 }
241 
242 /*@
243   DMPlexGetMaxProjectionHeight - Get the maximum height (w.r.t. DAG) of mesh points used to evaluate dual bases in
244   DMPlexProjectXXXLocal() functions.
245 
246   Input Parameters:
247 . dm - the DMPlex object
248 
249   Output Parameters:
250 . height - the maximum projection height
251 
252   Level: intermediate
253 
254 .seealso: DMPlexSetMaxProjectionHeight(), DMProjectFunctionLocal(), DMProjectFunctionLabelLocal()
255 @*/
256 PetscErrorCode DMPlexGetMaxProjectionHeight(DM dm, PetscInt *height)
257 {
258   DM_Plex *plex = (DM_Plex *) dm->data;
259 
260   PetscFunctionBegin;
261   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
262   *height = plex->maxProjectionHeight;
263   PetscFunctionReturn(0);
264 }
265 
266 /*@C
267   DMPlexInsertBoundaryValuesEssential - Insert boundary values into a local vector
268 
269   Input Parameters:
270 + dm     - The DM, with a PetscDS that matches the problem being constrained
271 . time   - The time
272 . field  - The field to constrain
273 . Nc     - The number of constrained field components, or 0 for all components
274 . comps  - An array of constrained component numbers, or NULL for all components
275 . label  - The DMLabel defining constrained points
276 . numids - The number of DMLabel ids for constrained points
277 . ids    - An array of ids for constrained points
278 . func   - A pointwise function giving boundary values
279 - ctx    - An optional user context for bcFunc
280 
281   Output Parameter:
282 . locX   - A local vector to receives the boundary values
283 
284   Level: developer
285 
286 .seealso: DMPlexInsertBoundaryValuesEssentialField(), DMAddBoundary()
287 @*/
288 PetscErrorCode DMPlexInsertBoundaryValuesEssential(DM dm, PetscReal time, PetscInt field, PetscInt Nc, const PetscInt comps[], DMLabel label, PetscInt numids, const PetscInt ids[], PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void *ctx, Vec locX)
289 {
290   PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
291   void            **ctxs;
292   PetscInt          numFields;
293   PetscErrorCode    ierr;
294 
295   PetscFunctionBegin;
296   ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr);
297   ierr = PetscCalloc2(numFields,&funcs,numFields,&ctxs);CHKERRQ(ierr);
298   funcs[field] = func;
299   ctxs[field]  = ctx;
300   ierr = DMProjectFunctionLabelLocal(dm, time, label, numids, ids, Nc, comps, funcs, ctxs, INSERT_BC_VALUES, locX);CHKERRQ(ierr);
301   ierr = PetscFree2(funcs,ctxs);CHKERRQ(ierr);
302   PetscFunctionReturn(0);
303 }
304 
305 /*@C
306   DMPlexInsertBoundaryValuesEssentialField - Insert boundary values into a local vector
307 
308   Input Parameters:
309 + dm     - The DM, with a PetscDS that matches the problem being constrained
310 . time   - The time
311 . locU   - A local vector with the input solution values
312 . field  - The field to constrain
313 . Nc     - The number of constrained field components, or 0 for all components
314 . comps  - An array of constrained component numbers, or NULL for all components
315 . label  - The DMLabel defining constrained points
316 . numids - The number of DMLabel ids for constrained points
317 . ids    - An array of ids for constrained points
318 . func   - A pointwise function giving boundary values
319 - ctx    - An optional user context for bcFunc
320 
321   Output Parameter:
322 . locX   - A local vector to receives the boundary values
323 
324   Level: developer
325 
326 .seealso: DMPlexInsertBoundaryValuesEssential(), DMAddBoundary()
327 @*/
328 PetscErrorCode DMPlexInsertBoundaryValuesEssentialField(DM dm, PetscReal time, Vec locU, PetscInt field, PetscInt Nc, const PetscInt comps[], DMLabel label, PetscInt numids, const PetscInt ids[],
329                                                         void (*func)(PetscInt, PetscInt, PetscInt,
330                                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
331                                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
332                                                                      PetscReal, const PetscReal[], PetscInt, const PetscScalar[],
333                                                                      PetscScalar[]),
334                                                         void *ctx, Vec locX)
335 {
336   void (**funcs)(PetscInt, PetscInt, PetscInt,
337                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
338                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
339                  PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);
340   void            **ctxs;
341   PetscInt          numFields;
342   PetscErrorCode    ierr;
343 
344   PetscFunctionBegin;
345   ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr);
346   ierr = PetscCalloc2(numFields,&funcs,numFields,&ctxs);CHKERRQ(ierr);
347   funcs[field] = func;
348   ctxs[field]  = ctx;
349   ierr = DMProjectFieldLabelLocal(dm, time, label, numids, ids, Nc, comps, locU, funcs, INSERT_BC_VALUES, locX);CHKERRQ(ierr);
350   ierr = PetscFree2(funcs,ctxs);CHKERRQ(ierr);
351   PetscFunctionReturn(0);
352 }
353 
354 /*@C
355   DMPlexInsertBoundaryValuesRiemann - Insert boundary values into a local vector
356 
357   Input Parameters:
358 + dm     - The DM, with a PetscDS that matches the problem being constrained
359 . time   - The time
360 . faceGeometry - A vector with the FVM face geometry information
361 . cellGeometry - A vector with the FVM cell geometry information
362 . Grad         - A vector with the FVM cell gradient information
363 . field  - The field to constrain
364 . Nc     - The number of constrained field components, or 0 for all components
365 . comps  - An array of constrained component numbers, or NULL for all components
366 . label  - The DMLabel defining constrained points
367 . numids - The number of DMLabel ids for constrained points
368 . ids    - An array of ids for constrained points
369 . func   - A pointwise function giving boundary values
370 - ctx    - An optional user context for bcFunc
371 
372   Output Parameter:
373 . locX   - A local vector to receives the boundary values
374 
375   Note: This implementation currently ignores the numcomps/comps argument from DMAddBoundary()
376 
377   Level: developer
378 
379 .seealso: DMPlexInsertBoundaryValuesEssential(), DMPlexInsertBoundaryValuesEssentialField(), DMAddBoundary()
380 @*/
381 PetscErrorCode DMPlexInsertBoundaryValuesRiemann(DM dm, PetscReal time, Vec faceGeometry, Vec cellGeometry, Vec Grad, PetscInt field, PetscInt Nc, const PetscInt comps[], DMLabel label, PetscInt numids, const PetscInt ids[],
382                                                  PetscErrorCode (*func)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*), void *ctx, Vec locX)
383 {
384   PetscDS            prob;
385   PetscSF            sf;
386   DM                 dmFace, dmCell, dmGrad;
387   const PetscScalar *facegeom, *cellgeom = NULL, *grad;
388   const PetscInt    *leaves;
389   PetscScalar       *x, *fx;
390   PetscInt           dim, nleaves, loc, fStart, fEnd, pdim, i;
391   PetscErrorCode     ierr, ierru = 0;
392 
393   PetscFunctionBegin;
394   ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
395   ierr = PetscSFGetGraph(sf, NULL, &nleaves, &leaves, NULL);CHKERRQ(ierr);
396   nleaves = PetscMax(0, nleaves);
397   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
398   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
399   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
400   ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr);
401   ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
402   if (cellGeometry) {
403     ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr);
404     ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
405   }
406   if (Grad) {
407     PetscFV fv;
408 
409     ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fv);CHKERRQ(ierr);
410     ierr = VecGetDM(Grad, &dmGrad);CHKERRQ(ierr);
411     ierr = VecGetArrayRead(Grad, &grad);CHKERRQ(ierr);
412     ierr = PetscFVGetNumComponents(fv, &pdim);CHKERRQ(ierr);
413     ierr = DMGetWorkArray(dm, pdim, MPIU_SCALAR, &fx);CHKERRQ(ierr);
414   }
415   ierr = VecGetArray(locX, &x);CHKERRQ(ierr);
416   for (i = 0; i < numids; ++i) {
417     IS              faceIS;
418     const PetscInt *faces;
419     PetscInt        numFaces, f;
420 
421     ierr = DMLabelGetStratumIS(label, ids[i], &faceIS);CHKERRQ(ierr);
422     if (!faceIS) continue; /* No points with that id on this process */
423     ierr = ISGetLocalSize(faceIS, &numFaces);CHKERRQ(ierr);
424     ierr = ISGetIndices(faceIS, &faces);CHKERRQ(ierr);
425     for (f = 0; f < numFaces; ++f) {
426       const PetscInt         face = faces[f], *cells;
427       PetscFVFaceGeom        *fg;
428 
429       if ((face < fStart) || (face >= fEnd)) continue; /* Refinement adds non-faces to labels */
430       ierr = PetscFindInt(face, nleaves, (PetscInt *) leaves, &loc);CHKERRQ(ierr);
431       if (loc >= 0) continue;
432       ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr);
433       ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr);
434       if (Grad) {
435         PetscFVCellGeom       *cg;
436         PetscScalar           *cx, *cgrad;
437         PetscScalar           *xG;
438         PetscReal              dx[3];
439         PetscInt               d;
440 
441         ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cg);CHKERRQ(ierr);
442         ierr = DMPlexPointLocalRead(dm, cells[0], x, &cx);CHKERRQ(ierr);
443         ierr = DMPlexPointLocalRead(dmGrad, cells[0], grad, &cgrad);CHKERRQ(ierr);
444         ierr = DMPlexPointLocalFieldRef(dm, cells[1], field, x, &xG);CHKERRQ(ierr);
445         DMPlex_WaxpyD_Internal(dim, -1, cg->centroid, fg->centroid, dx);
446         for (d = 0; d < pdim; ++d) fx[d] = cx[d] + DMPlex_DotD_Internal(dim, &cgrad[d*dim], dx);
447         ierru = (*func)(time, fg->centroid, fg->normal, fx, xG, ctx);
448         if (ierru) {
449           ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr);
450           ierr = ISDestroy(&faceIS);CHKERRQ(ierr);
451           goto cleanup;
452         }
453       } else {
454         PetscScalar       *xI;
455         PetscScalar       *xG;
456 
457         ierr = DMPlexPointLocalRead(dm, cells[0], x, &xI);CHKERRQ(ierr);
458         ierr = DMPlexPointLocalFieldRef(dm, cells[1], field, x, &xG);CHKERRQ(ierr);
459         ierru = (*func)(time, fg->centroid, fg->normal, xI, xG, ctx);
460         if (ierru) {
461           ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr);
462           ierr = ISDestroy(&faceIS);CHKERRQ(ierr);
463           goto cleanup;
464         }
465       }
466     }
467     ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr);
468     ierr = ISDestroy(&faceIS);CHKERRQ(ierr);
469   }
470   cleanup:
471   ierr = VecRestoreArray(locX, &x);CHKERRQ(ierr);
472   if (Grad) {
473     ierr = DMRestoreWorkArray(dm, pdim, MPIU_SCALAR, &fx);CHKERRQ(ierr);
474     ierr = VecRestoreArrayRead(Grad, &grad);CHKERRQ(ierr);
475   }
476   if (cellGeometry) {ierr = VecRestoreArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);}
477   ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
478   CHKERRQ(ierru);
479   PetscFunctionReturn(0);
480 }
481 
482 PetscErrorCode DMPlexInsertBoundaryValues_Plex(DM dm, PetscBool insertEssential, Vec locX, PetscReal time, Vec faceGeomFVM, Vec cellGeomFVM, Vec gradFVM)
483 {
484   PetscInt       numBd, b;
485   PetscErrorCode ierr;
486 
487   PetscFunctionBegin;
488   ierr = PetscDSGetNumBoundary(dm->prob, &numBd);CHKERRQ(ierr);
489   for (b = 0; b < numBd; ++b) {
490     DMBoundaryConditionType type;
491     const char             *labelname;
492     DMLabel                 label;
493     PetscInt                field, Nc;
494     const PetscInt         *comps;
495     PetscObject             obj;
496     PetscClassId            id;
497     void                    (*func)(void);
498     PetscInt                numids;
499     const PetscInt         *ids;
500     void                   *ctx;
501 
502     ierr = DMGetBoundary(dm, b, &type, NULL, &labelname, &field, &Nc, &comps, &func, &numids, &ids, &ctx);CHKERRQ(ierr);
503     if (insertEssential != (type & DM_BC_ESSENTIAL)) continue;
504     ierr = DMGetLabel(dm, labelname, &label);CHKERRQ(ierr);
505     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
506     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
507     if (id == PETSCFE_CLASSID) {
508       switch (type) {
509         /* for FEM, there is no insertion to be done for non-essential boundary conditions */
510       case DM_BC_ESSENTIAL:
511         ierr = DMPlexLabelAddCells(dm,label);CHKERRQ(ierr);
512         ierr = DMPlexInsertBoundaryValuesEssential(dm, time, field, Nc, comps, label, numids, ids, (PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *)) func, ctx, locX);CHKERRQ(ierr);
513         ierr = DMPlexLabelClearCells(dm,label);CHKERRQ(ierr);
514         break;
515       case DM_BC_ESSENTIAL_FIELD:
516         ierr = DMPlexLabelAddCells(dm,label);CHKERRQ(ierr);
517         ierr = DMPlexInsertBoundaryValuesEssentialField(dm, time, locX, field, Nc, comps, label, numids, ids,
518                                                         (void (*)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
519                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
520                                                                   PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) func, ctx, locX);CHKERRQ(ierr);
521         ierr = DMPlexLabelClearCells(dm,label);CHKERRQ(ierr);
522         break;
523       default: break;
524       }
525     } else if (id == PETSCFV_CLASSID) {
526       if (!faceGeomFVM) continue;
527       ierr = DMPlexInsertBoundaryValuesRiemann(dm, time, faceGeomFVM, cellGeomFVM, gradFVM, field, Nc, comps, label, numids, ids,
528                                                (PetscErrorCode (*)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*)) func, ctx, locX);CHKERRQ(ierr);
529     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
530   }
531   PetscFunctionReturn(0);
532 }
533 
534 /*@
535   DMPlexInsertBoundaryValues - Puts coefficients which represent boundary values into the local solution vector
536 
537   Input Parameters:
538 + dm - The DM
539 . insertEssential - Should I insert essential (e.g. Dirichlet) or inessential (e.g. Neumann) boundary conditions
540 . time - The time
541 . faceGeomFVM - Face geometry data for FV discretizations
542 . cellGeomFVM - Cell geometry data for FV discretizations
543 - gradFVM - Gradient reconstruction data for FV discretizations
544 
545   Output Parameters:
546 . locX - Solution updated with boundary values
547 
548   Level: developer
549 
550 .seealso: DMProjectFunctionLabelLocal()
551 @*/
552 PetscErrorCode DMPlexInsertBoundaryValues(DM dm, PetscBool insertEssential, Vec locX, PetscReal time, Vec faceGeomFVM, Vec cellGeomFVM, Vec gradFVM)
553 {
554   PetscErrorCode ierr;
555 
556   PetscFunctionBegin;
557   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
558   PetscValidHeaderSpecific(locX, VEC_CLASSID, 2);
559   if (faceGeomFVM) {PetscValidHeaderSpecific(faceGeomFVM, VEC_CLASSID, 4);}
560   if (cellGeomFVM) {PetscValidHeaderSpecific(cellGeomFVM, VEC_CLASSID, 5);}
561   if (gradFVM)     {PetscValidHeaderSpecific(gradFVM, VEC_CLASSID, 6);}
562   ierr = PetscTryMethod(dm,"DMPlexInsertBoundaryValues_C",(DM,PetscBool,Vec,PetscReal,Vec,Vec,Vec),(dm,insertEssential,locX,time,faceGeomFVM,cellGeomFVM,gradFVM));CHKERRQ(ierr);
563   PetscFunctionReturn(0);
564 }
565 
566 PetscErrorCode DMComputeL2Diff_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
567 {
568   Vec              localX;
569   PetscErrorCode   ierr;
570 
571   PetscFunctionBegin;
572   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
573   ierr = DMPlexInsertBoundaryValues(dm, PETSC_TRUE, localX, time, NULL, NULL, NULL);CHKERRQ(ierr);
574   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
575   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
576   ierr = DMPlexComputeL2DiffLocal(dm, time, funcs, ctxs, localX, diff);CHKERRQ(ierr);
577   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
578   PetscFunctionReturn(0);
579 }
580 
581 /*@C
582   DMComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.
583 
584   Input Parameters:
585 + dm     - The DM
586 . time   - The time
587 . funcs  - The functions to evaluate for each field component
588 . ctxs   - Optional array of contexts to pass to each function, or NULL.
589 - localX - The coefficient vector u_h, a local vector
590 
591   Output Parameter:
592 . diff - The diff ||u - u_h||_2
593 
594   Level: developer
595 
596 .seealso: DMProjectFunction(), DMComputeL2FieldDiff(), DMComputeL2GradientDiff()
597 @*/
598 PetscErrorCode DMPlexComputeL2DiffLocal(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec localX, PetscReal *diff)
599 {
600   const PetscInt   debug = 0;
601   PetscSection     section;
602   PetscQuadrature  quad;
603   PetscScalar     *funcVal, *interpolant;
604   PetscReal       *coords, *detJ, *J;
605   PetscReal        localDiff = 0.0;
606   const PetscReal *quadWeights;
607   PetscInt         dim, coordDim, numFields, numComponents = 0, qNc, Nq, cellHeight, cStart, cEnd, cEndInterior, c, field, fieldOffset;
608   PetscErrorCode   ierr;
609 
610   PetscFunctionBegin;
611   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
612   ierr = DMGetCoordinateDim(dm, &coordDim);CHKERRQ(ierr);
613   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
614   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
615   for (field = 0; field < numFields; ++field) {
616     PetscObject  obj;
617     PetscClassId id;
618     PetscInt     Nc;
619 
620     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
621     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
622     if (id == PETSCFE_CLASSID) {
623       PetscFE fe = (PetscFE) obj;
624 
625       ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
626       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
627     } else if (id == PETSCFV_CLASSID) {
628       PetscFV fv = (PetscFV) obj;
629 
630       ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr);
631       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
632     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
633     numComponents += Nc;
634   }
635   ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, NULL, &quadWeights);CHKERRQ(ierr);
636   if ((qNc != 1) && (qNc != numComponents)) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_SIZ, "Quadrature components %D != %D field components", qNc, numComponents);
637   ierr = PetscMalloc5(numComponents,&funcVal,numComponents,&interpolant,coordDim*Nq,&coords,Nq,&detJ,coordDim*coordDim*Nq,&J);CHKERRQ(ierr);
638   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
639   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
640   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
641   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
642   for (c = cStart; c < cEnd; ++c) {
643     PetscScalar *x = NULL;
644     PetscReal    elemDiff = 0.0;
645     PetscInt     qc = 0;
646 
647     ierr = DMPlexComputeCellGeometryFEM(dm, c, quad, coords, J, NULL, detJ);CHKERRQ(ierr);
648     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
649 
650     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
651       PetscObject  obj;
652       PetscClassId id;
653       void * const ctx = ctxs ? ctxs[field] : NULL;
654       PetscInt     Nb, Nc, q, fc;
655 
656       ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
657       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
658       if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);}
659       else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;}
660       else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
661       if (debug) {
662         char title[1024];
663         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
664         ierr = DMPrintCellVector(c, title, Nb, &x[fieldOffset]);CHKERRQ(ierr);
665       }
666       for (q = 0; q < Nq; ++q) {
667         if (detJ[q] <= 0.0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D, point %D", detJ[q], c, q);
668         ierr = (*funcs[field])(coordDim, time, &coords[coordDim * q], Nc, funcVal, ctx);
669         if (ierr) {
670           PetscErrorCode ierr2;
671           ierr2 = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr2);
672           ierr2 = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr2);
673           ierr2 = PetscFree5(funcVal,interpolant,coords,detJ,J);CHKERRQ(ierr2);
674           CHKERRQ(ierr);
675         }
676         if (id == PETSCFE_CLASSID)      {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
677         else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
678         else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
679         for (fc = 0; fc < Nc; ++fc) {
680           const PetscReal wt = quadWeights[q*qNc+(qNc == 1 ? 0 : qc+fc)];
681           if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "    elem %d field %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*wt*detJ[q]);CHKERRQ(ierr);}
682           elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*wt*detJ[q];
683         }
684       }
685       fieldOffset += Nb;
686       qc += Nc;
687     }
688     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
689     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
690     localDiff += elemDiff;
691   }
692   ierr  = PetscFree5(funcVal,interpolant,coords,detJ,J);CHKERRQ(ierr);
693   ierr  = MPIU_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
694   *diff = PetscSqrtReal(*diff);
695   PetscFunctionReturn(0);
696 }
697 
698 PetscErrorCode DMComputeL2GradientDiff_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, const PetscReal n[], PetscReal *diff)
699 {
700   const PetscInt   debug = 0;
701   PetscSection     section;
702   PetscQuadrature  quad;
703   Vec              localX;
704   PetscScalar     *funcVal, *interpolant;
705   const PetscReal *quadPoints, *quadWeights;
706   PetscReal       *coords, *realSpaceDer, *J, *invJ, *detJ;
707   PetscReal        localDiff = 0.0;
708   PetscInt         dim, coordDim, qNc = 0, Nq = 0, numFields, numComponents = 0, cStart, cEnd, cEndInterior, c, field, fieldOffset;
709   PetscErrorCode   ierr;
710 
711   PetscFunctionBegin;
712   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
713   ierr = DMGetCoordinateDim(dm, &coordDim);CHKERRQ(ierr);
714   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
715   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
716   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
717   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
718   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
719   for (field = 0; field < numFields; ++field) {
720     PetscFE  fe;
721     PetscInt Nc;
722 
723     ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
724     ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
725     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
726     numComponents += Nc;
727   }
728   ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
729   if ((qNc != 1) && (qNc != numComponents)) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_SIZ, "Quadrature components %D != %D field components", qNc, numComponents);
730   /* ierr = DMProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */
731   ierr = PetscMalloc7(numComponents,&funcVal,coordDim*Nq,&coords,coordDim*Nq,&realSpaceDer,coordDim*coordDim*Nq,&J,coordDim*coordDim*Nq,&invJ,numComponents,&interpolant,Nq,&detJ);CHKERRQ(ierr);
732   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
733   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
734   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
735   for (c = cStart; c < cEnd; ++c) {
736     PetscScalar *x = NULL;
737     PetscReal    elemDiff = 0.0;
738     PetscInt     qc = 0;
739 
740     ierr = DMPlexComputeCellGeometryFEM(dm, c, quad, coords, J, invJ, detJ);CHKERRQ(ierr);
741     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
742 
743     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
744       PetscFE          fe;
745       void * const     ctx = ctxs ? ctxs[field] : NULL;
746       PetscReal       *basisDer;
747       PetscInt         Nb, Nc, q, fc;
748 
749       ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
750       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
751       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
752       ierr = PetscFEGetDefaultTabulation(fe, NULL, &basisDer, NULL);CHKERRQ(ierr);
753       if (debug) {
754         char title[1024];
755         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
756         ierr = DMPrintCellVector(c, title, Nb, &x[fieldOffset]);CHKERRQ(ierr);
757       }
758       for (q = 0; q < Nq; ++q) {
759         if (detJ[q] <= 0.0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D, quadrature points %D", detJ[q], c, q);
760         ierr = (*funcs[field])(coordDim, time, &coords[q*coordDim], n, numFields, funcVal, ctx);
761         if (ierr) {
762           PetscErrorCode ierr2;
763           ierr2 = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr2);
764           ierr2 = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr2);
765           ierr2 = PetscFree7(funcVal,coords,realSpaceDer,J,invJ,interpolant,detJ);CHKERRQ(ierr2);
766           CHKERRQ(ierr);
767         }
768         ierr = PetscFEInterpolateGradient_Static(fe, &x[fieldOffset], coordDim, invJ, n, q, interpolant);CHKERRQ(ierr);
769         for (fc = 0; fc < Nc; ++fc) {
770           const PetscReal wt = quadWeights[q*qNc+(qNc == 1 ? 0 : qc+fc)];
771           if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "    elem %d fieldDer %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*wt*detJ[q]);CHKERRQ(ierr);}
772           elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*wt*detJ[q];
773         }
774       }
775       fieldOffset += Nb;
776       qc          += Nc;
777     }
778     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
779     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
780     localDiff += elemDiff;
781   }
782   ierr  = PetscFree7(funcVal,coords,realSpaceDer,J,invJ,interpolant,detJ);CHKERRQ(ierr);
783   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
784   ierr  = MPIU_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
785   *diff = PetscSqrtReal(*diff);
786   PetscFunctionReturn(0);
787 }
788 
789 PetscErrorCode DMComputeL2FieldDiff_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
790 {
791   const PetscInt   debug = 0;
792   PetscSection     section;
793   PetscQuadrature  quad;
794   Vec              localX;
795   PetscScalar     *funcVal, *interpolant;
796   PetscReal       *coords, *detJ, *J;
797   PetscReal       *localDiff;
798   const PetscReal *quadPoints, *quadWeights;
799   PetscInt         dim, coordDim, numFields, numComponents = 0, qNc, Nq, cStart, cEnd, cEndInterior, c, field, fieldOffset;
800   PetscErrorCode   ierr;
801 
802   PetscFunctionBegin;
803   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
804   ierr = DMGetCoordinateDim(dm, &coordDim);CHKERRQ(ierr);
805   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
806   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
807   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
808   ierr = DMProjectFunctionLocal(dm, time, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
809   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
810   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
811   for (field = 0; field < numFields; ++field) {
812     PetscObject  obj;
813     PetscClassId id;
814     PetscInt     Nc;
815 
816     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
817     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
818     if (id == PETSCFE_CLASSID) {
819       PetscFE fe = (PetscFE) obj;
820 
821       ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
822       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
823     } else if (id == PETSCFV_CLASSID) {
824       PetscFV fv = (PetscFV) obj;
825 
826       ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr);
827       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
828     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
829     numComponents += Nc;
830   }
831   ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
832   if ((qNc != 1) && (qNc != numComponents)) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_SIZ, "Quadrature components %D != %D field components", qNc, numComponents);
833   ierr = PetscCalloc6(numFields,&localDiff,numComponents,&funcVal,numComponents,&interpolant,coordDim*Nq,&coords,Nq,&detJ,coordDim*coordDim*Nq,&J);CHKERRQ(ierr);
834   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
835   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
836   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
837   for (c = cStart; c < cEnd; ++c) {
838     PetscScalar *x = NULL;
839     PetscInt     qc = 0;
840 
841     ierr = DMPlexComputeCellGeometryFEM(dm, c, quad, coords, J, NULL, detJ);CHKERRQ(ierr);
842     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
843 
844     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
845       PetscObject  obj;
846       PetscClassId id;
847       void * const ctx = ctxs ? ctxs[field] : NULL;
848       PetscInt     Nb, Nc, q, fc;
849 
850       PetscReal       elemDiff = 0.0;
851 
852       ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
853       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
854       if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);}
855       else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;}
856       else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
857       if (debug) {
858         char title[1024];
859         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
860         ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr);
861       }
862       for (q = 0; q < Nq; ++q) {
863         if (detJ[q] <= 0.0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D, quadrature point %D", detJ, c, q);
864         ierr = (*funcs[field])(coordDim, time, &coords[coordDim*q], numFields, funcVal, ctx);
865         if (ierr) {
866           PetscErrorCode ierr2;
867           ierr2 = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr2);
868           ierr2 = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr2);
869           ierr2 = PetscFree6(localDiff,funcVal,interpolant,coords,detJ,J);CHKERRQ(ierr2);
870           CHKERRQ(ierr);
871         }
872         if (id == PETSCFE_CLASSID)      {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
873         else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
874         else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
875         for (fc = 0; fc < Nc; ++fc) {
876           const PetscReal wt = quadWeights[q*qNc+(qNc == 1 ? 0 : qc+fc)];
877           if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "    elem %d field %d point %g %g %g diff %g\n", c, field, coordDim > 0 ? coords[0] : 0., coordDim > 1 ? coords[1] : 0., coordDim > 2 ? coords[2] : 0., PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*wt*detJ[q]);CHKERRQ(ierr);}
878           elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*wt*detJ[q];
879         }
880       }
881       fieldOffset += Nb;
882       qc          += Nc;
883       localDiff[field] += elemDiff;
884     }
885     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
886   }
887   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
888   ierr = MPIU_Allreduce(localDiff, diff, numFields, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
889   for (field = 0; field < numFields; ++field) diff[field] = PetscSqrtReal(diff[field]);
890   ierr = PetscFree6(localDiff,funcVal,interpolant,coords,detJ,J);CHKERRQ(ierr);
891   PetscFunctionReturn(0);
892 }
893 
894 /*@C
895   DMPlexComputeL2DiffVec - This function computes the cellwise L_2 difference between a function u and an FEM interpolant solution u_h, and stores it in a Vec.
896 
897   Input Parameters:
898 + dm    - The DM
899 . time  - The time
900 . funcs - The functions to evaluate for each field component: NULL means that component does not contribute to error calculation
901 . ctxs  - Optional array of contexts to pass to each function, or NULL.
902 - X     - The coefficient vector u_h
903 
904   Output Parameter:
905 . D - A Vec which holds the difference ||u - u_h||_2 for each cell
906 
907   Level: developer
908 
909 .seealso: DMProjectFunction(), DMComputeL2Diff(), DMPlexComputeL2FieldDiff(), DMComputeL2GradientDiff()
910 @*/
911 PetscErrorCode DMPlexComputeL2DiffVec(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, Vec D)
912 {
913   PetscSection     section;
914   PetscQuadrature  quad;
915   Vec              localX;
916   PetscScalar     *funcVal, *interpolant;
917   PetscReal       *coords, *detJ, *J;
918   const PetscReal *quadPoints, *quadWeights;
919   PetscInt         dim, coordDim, numFields, numComponents = 0, qNc, Nq, cStart, cEnd, cEndInterior, c, field, fieldOffset;
920   PetscErrorCode   ierr;
921 
922   PetscFunctionBegin;
923   ierr = VecSet(D, 0.0);CHKERRQ(ierr);
924   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
925   ierr = DMGetCoordinateDim(dm, &coordDim);CHKERRQ(ierr);
926   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
927   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
928   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
929   ierr = DMProjectFunctionLocal(dm, time, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
930   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
931   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
932   for (field = 0; field < numFields; ++field) {
933     PetscObject  obj;
934     PetscClassId id;
935     PetscInt     Nc;
936 
937     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
938     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
939     if (id == PETSCFE_CLASSID) {
940       PetscFE fe = (PetscFE) obj;
941 
942       ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
943       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
944     } else if (id == PETSCFV_CLASSID) {
945       PetscFV fv = (PetscFV) obj;
946 
947       ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr);
948       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
949     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
950     numComponents += Nc;
951   }
952   ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
953   if ((qNc != 1) && (qNc != numComponents)) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_SIZ, "Quadrature components %D != %D field components", qNc, numComponents);
954   ierr = PetscMalloc5(numComponents,&funcVal,numComponents,&interpolant,coordDim*Nq,&coords,Nq,&detJ,coordDim*coordDim*Nq,&J);CHKERRQ(ierr);
955   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
956   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
957   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
958   for (c = cStart; c < cEnd; ++c) {
959     PetscScalar *x = NULL;
960     PetscScalar  elemDiff = 0.0;
961     PetscInt     qc = 0;
962 
963     ierr = DMPlexComputeCellGeometryFEM(dm, c, quad, coords, J, NULL, detJ);CHKERRQ(ierr);
964     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
965 
966     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
967       PetscObject  obj;
968       PetscClassId id;
969       void * const ctx = ctxs ? ctxs[field] : NULL;
970       PetscInt     Nb, Nc, q, fc;
971 
972       ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
973       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
974       if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);}
975       else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;}
976       else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
977       if (funcs[field]) {
978         for (q = 0; q < Nq; ++q) {
979           if (detJ[q] <= 0.0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D, quadrature points %D", (double)detJ[q], c, q);
980           ierr = (*funcs[field])(dim, time, &coords[q*coordDim], Nc, funcVal, ctx);
981           if (ierr) {
982             PetscErrorCode ierr2;
983             ierr2 = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr2);
984             ierr2 = PetscFree5(funcVal,interpolant,coords,detJ,J);CHKERRQ(ierr2);
985             ierr2 = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr2);
986             CHKERRQ(ierr);
987           }
988           if (id == PETSCFE_CLASSID)      {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
989           else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
990           else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
991           for (fc = 0; fc < Nc; ++fc) {
992             const PetscReal wt = quadWeights[q*qNc+(qNc == 1 ? 0 : qc+fc)];
993             elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*wt*detJ[q];
994           }
995         }
996       }
997       fieldOffset += Nb;
998       qc          += Nc;
999     }
1000     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
1001     ierr = VecSetValue(D, c - cStart, elemDiff, INSERT_VALUES);CHKERRQ(ierr);
1002   }
1003   ierr = PetscFree5(funcVal,interpolant,coords,detJ,J);CHKERRQ(ierr);
1004   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
1005   ierr = VecSqrtAbs(D);CHKERRQ(ierr);
1006   PetscFunctionReturn(0);
1007 }
1008 
1009 /*@C
1010   DMPlexComputeGradientClementInterpolant - This function computes the L2 projection of the cellwise gradient of a function u onto P1, and stores it in a Vec.
1011 
1012   Input Parameters:
1013 + dm - The DM
1014 - LocX  - The coefficient vector u_h
1015 
1016   Output Parameter:
1017 . locC - A Vec which holds the Clement interpolant of the gradient
1018 
1019   Notes:
1020     Add citation to (Clement, 1975) and definition of the interpolant
1021   \nabla u_h(v_i) = \sum_{T_i \in support(v_i)} |T_i| \nabla u_h(T_i) / \sum_{T_i \in support(v_i)} |T_i| where |T_i| is the cell volume
1022 
1023   Level: developer
1024 
1025 .seealso: DMProjectFunction(), DMComputeL2Diff(), DMPlexComputeL2FieldDiff(), DMComputeL2GradientDiff()
1026 @*/
1027 PetscErrorCode DMPlexComputeGradientClementInterpolant(DM dm, Vec locX, Vec locC)
1028 {
1029   DM_Plex         *mesh  = (DM_Plex *) dm->data;
1030   PetscInt         debug = mesh->printFEM;
1031   DM               dmC;
1032   PetscSection     section;
1033   PetscQuadrature  quad;
1034   PetscScalar     *interpolant, *gradsum;
1035   PetscReal       *coords, *detJ, *J, *invJ;
1036   const PetscReal *quadPoints, *quadWeights;
1037   PetscInt         dim, coordDim, numFields, numComponents = 0, qNc, Nq, cStart, cEnd, cEndInterior, vStart, vEnd, v, field, fieldOffset;
1038   PetscErrorCode   ierr;
1039 
1040   PetscFunctionBegin;
1041   ierr = VecGetDM(locC, &dmC);CHKERRQ(ierr);
1042   ierr = VecSet(locC, 0.0);CHKERRQ(ierr);
1043   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1044   ierr = DMGetCoordinateDim(dm, &coordDim);CHKERRQ(ierr);
1045   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1046   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
1047   for (field = 0; field < numFields; ++field) {
1048     PetscObject  obj;
1049     PetscClassId id;
1050     PetscInt     Nc;
1051 
1052     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
1053     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1054     if (id == PETSCFE_CLASSID) {
1055       PetscFE fe = (PetscFE) obj;
1056 
1057       ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
1058       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1059     } else if (id == PETSCFV_CLASSID) {
1060       PetscFV fv = (PetscFV) obj;
1061 
1062       ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr);
1063       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
1064     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
1065     numComponents += Nc;
1066   }
1067   ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
1068   if ((qNc != 1) && (qNc != numComponents)) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_SIZ, "Quadrature components %D != %D field components", qNc, numComponents);
1069   ierr = PetscMalloc6(coordDim*numComponents*2,&gradsum,coordDim*numComponents,&interpolant,coordDim*Nq,&coords,Nq,&detJ,coordDim*coordDim*Nq,&J,coordDim*coordDim*Nq,&invJ);CHKERRQ(ierr);
1070   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
1071   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1072   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
1073   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
1074   for (v = vStart; v < vEnd; ++v) {
1075     PetscScalar volsum = 0.0;
1076     PetscInt   *star = NULL;
1077     PetscInt    starSize, st, d, fc;
1078 
1079     ierr = PetscMemzero(gradsum, coordDim*numComponents * sizeof(PetscScalar));CHKERRQ(ierr);
1080     ierr = DMPlexGetTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
1081     for (st = 0; st < starSize*2; st += 2) {
1082       const PetscInt cell = star[st];
1083       PetscScalar   *grad = &gradsum[coordDim*numComponents];
1084       PetscScalar   *x    = NULL;
1085       PetscReal      vol  = 0.0;
1086 
1087       if ((cell < cStart) || (cell >= cEnd)) continue;
1088       ierr = DMPlexComputeCellGeometryFEM(dm, cell, quad, coords, J, invJ, detJ);CHKERRQ(ierr);
1089       ierr = DMPlexVecGetClosure(dm, NULL, locX, cell, NULL, &x);CHKERRQ(ierr);
1090       for (field = 0, fieldOffset = 0; field < numFields; ++field) {
1091         PetscObject  obj;
1092         PetscClassId id;
1093         PetscInt     Nb, Nc, q, qc = 0;
1094 
1095         ierr = PetscMemzero(grad, coordDim*numComponents * sizeof(PetscScalar));CHKERRQ(ierr);
1096         ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
1097         ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1098         if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);}
1099         else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;}
1100         else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
1101         for (q = 0; q < Nq; ++q) {
1102           if (detJ[q] <= 0.0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D, quadrature points %D", (double)detJ[q], cell, q);
1103           if (ierr) {
1104             PetscErrorCode ierr2;
1105             ierr2 = DMPlexVecRestoreClosure(dm, NULL, locX, cell, NULL, &x);CHKERRQ(ierr2);
1106             ierr2 = DMPlexRestoreTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr2);
1107             ierr2 = PetscFree4(interpolant,coords,detJ,J);CHKERRQ(ierr2);
1108             CHKERRQ(ierr);
1109           }
1110           if (id == PETSCFE_CLASSID)      {ierr = PetscFEInterpolateGradient_Static((PetscFE) obj, &x[fieldOffset], coordDim, invJ, NULL, q, interpolant);CHKERRQ(ierr);}
1111           else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
1112           for (fc = 0; fc < Nc; ++fc) {
1113             const PetscReal wt = quadWeights[q*qNc+qc+fc];
1114 
1115             for (d = 0; d < coordDim; ++d) grad[fc*coordDim+d] += interpolant[fc*dim+d]*wt*detJ[q];
1116           }
1117           vol += quadWeights[q*qNc]*detJ[q];
1118         }
1119         fieldOffset += Nb;
1120         qc          += Nc;
1121       }
1122       ierr = DMPlexVecRestoreClosure(dm, NULL, locX, cell, NULL, &x);CHKERRQ(ierr);
1123       for (fc = 0; fc < numComponents; ++fc) {
1124         for (d = 0; d < coordDim; ++d) {
1125           gradsum[fc*coordDim+d] += grad[fc*coordDim+d];
1126         }
1127       }
1128       volsum += vol;
1129       if (debug) {
1130         ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D gradient: [", cell);CHKERRQ(ierr);
1131         for (fc = 0; fc < numComponents; ++fc) {
1132           for (d = 0; d < coordDim; ++d) {
1133             if (fc || d > 0) {ierr = PetscPrintf(PETSC_COMM_SELF, ", ");CHKERRQ(ierr);}
1134             ierr = PetscPrintf(PETSC_COMM_SELF, "%g", (double)PetscRealPart(grad[fc*coordDim+d]));CHKERRQ(ierr);
1135           }
1136         }
1137         ierr = PetscPrintf(PETSC_COMM_SELF, "]\n");CHKERRQ(ierr);
1138       }
1139     }
1140     for (fc = 0; fc < numComponents; ++fc) {
1141       for (d = 0; d < coordDim; ++d) gradsum[fc*coordDim+d] /= volsum;
1142     }
1143     ierr = DMPlexRestoreTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
1144     ierr = DMPlexVecSetClosure(dmC, NULL, locC, v, gradsum, INSERT_VALUES);CHKERRQ(ierr);
1145   }
1146   ierr = PetscFree6(gradsum,interpolant,coords,detJ,J,invJ);CHKERRQ(ierr);
1147   PetscFunctionReturn(0);
1148 }
1149 
1150 static PetscErrorCode DMPlexComputeIntegral_Internal(DM dm, Vec X, PetscInt cStart, PetscInt cEnd, PetscScalar *cintegral, void *user)
1151 {
1152   DM                 dmAux = NULL;
1153   PetscDS            prob,    probAux = NULL;
1154   PetscSection       section, sectionAux;
1155   Vec                locX,    locA;
1156   PetscInt           dim, numCells = cEnd - cStart, c, f;
1157   PetscBool          useFVM = PETSC_FALSE;
1158   /* DS */
1159   PetscInt           Nf,    totDim,    *uOff, *uOff_x, numConstants;
1160   PetscInt           NfAux, totDimAux, *aOff;
1161   PetscScalar       *u, *a;
1162   const PetscScalar *constants;
1163   /* Geometry */
1164   PetscFEGeom       *cgeomFEM;
1165   DM                 dmGrad;
1166   PetscQuadrature    affineQuad = NULL;
1167   Vec                cellGeometryFVM = NULL, faceGeometryFVM = NULL, locGrad = NULL;
1168   PetscFVCellGeom   *cgeomFVM;
1169   const PetscScalar *lgrad;
1170   PetscBool          isAffine;
1171   DMField            coordField;
1172   IS                 cellIS;
1173   PetscErrorCode     ierr;
1174 
1175   PetscFunctionBegin;
1176   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
1177   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1178   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1179   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
1180   /* Determine which discretizations we have */
1181   for (f = 0; f < Nf; ++f) {
1182     PetscObject  obj;
1183     PetscClassId id;
1184 
1185     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
1186     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1187     if (id == PETSCFV_CLASSID) useFVM = PETSC_TRUE;
1188   }
1189   /* Get local solution with boundary values */
1190   ierr = DMGetLocalVector(dm, &locX);CHKERRQ(ierr);
1191   ierr = DMPlexInsertBoundaryValues(dm, PETSC_TRUE, locX, 0.0, NULL, NULL, NULL);CHKERRQ(ierr);
1192   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX);CHKERRQ(ierr);
1193   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX);CHKERRQ(ierr);
1194   /* Read DS information */
1195   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
1196   ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr);
1197   ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr);
1198   ierr = ISCreateStride(PETSC_COMM_SELF,numCells,cStart,1,&cellIS);CHKERRQ(ierr);
1199   ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr);
1200   /* Read Auxiliary DS information */
1201   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
1202   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &locA);CHKERRQ(ierr);
1203   if (dmAux) {
1204     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
1205     ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr);
1206     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
1207     ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
1208     ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr);
1209   }
1210   /* Allocate data  arrays */
1211   ierr = PetscCalloc1(numCells*totDim, &u);CHKERRQ(ierr);
1212   if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);}
1213   /* Read out geometry */
1214   ierr = DMGetCoordinateField(dm,&coordField);CHKERRQ(ierr);
1215   ierr = DMFieldGetFEInvariance(coordField,cellIS,NULL,&isAffine,NULL);CHKERRQ(ierr);
1216   if (isAffine) {
1217     ierr = DMFieldCreateDefaultQuadrature(coordField,cellIS,&affineQuad);CHKERRQ(ierr);
1218     if (affineQuad) {
1219       ierr = DMFieldCreateFEGeom(coordField,cellIS,affineQuad,PETSC_FALSE,&cgeomFEM);CHKERRQ(ierr);
1220     }
1221   }
1222   if (useFVM) {
1223     PetscFV   fv = NULL;
1224     Vec       grad;
1225     PetscInt  fStart, fEnd;
1226     PetscBool compGrad;
1227 
1228     for (f = 0; f < Nf; ++f) {
1229       PetscObject  obj;
1230       PetscClassId id;
1231 
1232       ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
1233       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1234       if (id == PETSCFV_CLASSID) {fv = (PetscFV) obj; break;}
1235     }
1236     ierr = PetscFVGetComputeGradients(fv, &compGrad);CHKERRQ(ierr);
1237     ierr = PetscFVSetComputeGradients(fv, PETSC_TRUE);CHKERRQ(ierr);
1238     ierr = DMPlexComputeGeometryFVM(dm, &cellGeometryFVM, &faceGeometryFVM);CHKERRQ(ierr);
1239     ierr = DMPlexComputeGradientFVM(dm, fv, faceGeometryFVM, cellGeometryFVM, &dmGrad);CHKERRQ(ierr);
1240     ierr = PetscFVSetComputeGradients(fv, compGrad);CHKERRQ(ierr);
1241     ierr = VecGetArrayRead(cellGeometryFVM, (const PetscScalar **) &cgeomFVM);CHKERRQ(ierr);
1242     /* Reconstruct and limit cell gradients */
1243     ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
1244     ierr = DMGetGlobalVector(dmGrad, &grad);CHKERRQ(ierr);
1245     ierr = DMPlexReconstructGradients_Internal(dm, fv, fStart, fEnd, faceGeometryFVM, cellGeometryFVM, locX, grad);CHKERRQ(ierr);
1246     /* Communicate gradient values */
1247     ierr = DMGetLocalVector(dmGrad, &locGrad);CHKERRQ(ierr);
1248     ierr = DMGlobalToLocalBegin(dmGrad, grad, INSERT_VALUES, locGrad);CHKERRQ(ierr);
1249     ierr = DMGlobalToLocalEnd(dmGrad, grad, INSERT_VALUES, locGrad);CHKERRQ(ierr);
1250     ierr = DMRestoreGlobalVector(dmGrad, &grad);CHKERRQ(ierr);
1251     /* Handle non-essential (e.g. outflow) boundary values */
1252     ierr = DMPlexInsertBoundaryValues(dm, PETSC_FALSE, locX, 0.0, faceGeometryFVM, cellGeometryFVM, locGrad);CHKERRQ(ierr);
1253     ierr = VecGetArrayRead(locGrad, &lgrad);CHKERRQ(ierr);
1254   }
1255   /* Read out data from inputs */
1256   for (c = cStart; c < cEnd; ++c) {
1257     PetscScalar *x = NULL;
1258     PetscInt     i;
1259 
1260     ierr = DMPlexVecGetClosure(dm, section, locX, c, NULL, &x);CHKERRQ(ierr);
1261     for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i];
1262     ierr = DMPlexVecRestoreClosure(dm, section, locX, c, NULL, &x);CHKERRQ(ierr);
1263     if (dmAux) {
1264       ierr = DMPlexVecGetClosure(dmAux, sectionAux, locA, c, NULL, &x);CHKERRQ(ierr);
1265       for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i];
1266       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, locA, c, NULL, &x);CHKERRQ(ierr);
1267     }
1268   }
1269   /* Do integration for each field */
1270   for (f = 0; f < Nf; ++f) {
1271     PetscObject  obj;
1272     PetscClassId id;
1273     PetscInt     numChunks, numBatches, batchSize, numBlocks, blockSize, Ne, Nr, offset;
1274 
1275     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
1276     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1277     if (id == PETSCFE_CLASSID) {
1278       PetscFE         fe = (PetscFE) obj;
1279       PetscQuadrature q;
1280       PetscFEGeom     *chunkGeom = NULL;
1281       PetscInt        Nq, Nb;
1282 
1283       ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
1284       ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
1285       ierr = PetscQuadratureGetData(q, NULL, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);
1286       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
1287       blockSize = Nb*Nq;
1288       batchSize = numBlocks * blockSize;
1289       ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1290       numChunks = numCells / (numBatches*batchSize);
1291       Ne        = numChunks*numBatches*batchSize;
1292       Nr        = numCells % (numBatches*batchSize);
1293       offset    = numCells - Nr;
1294       if (!affineQuad) {
1295         ierr = DMFieldCreateFEGeom(coordField,cellIS,q,PETSC_FALSE,&cgeomFEM);CHKERRQ(ierr);
1296       }
1297       ierr = PetscFEGeomGetChunk(cgeomFEM,0,offset,&chunkGeom);CHKERRQ(ierr);
1298       ierr = PetscFEIntegrate(fe, prob, f, Ne, chunkGeom, u, probAux, a, cintegral);CHKERRQ(ierr);
1299       ierr = PetscFEGeomGetChunk(cgeomFEM,offset,numCells,&chunkGeom);CHKERRQ(ierr);
1300       ierr = PetscFEIntegrate(fe, prob, f, Nr, chunkGeom, &u[offset*totDim], probAux, &a[offset*totDimAux], &cintegral[offset*Nf]);CHKERRQ(ierr);
1301       ierr = PetscFEGeomRestoreChunk(cgeomFEM,offset,numCells,&chunkGeom);CHKERRQ(ierr);
1302       if (!affineQuad) {
1303         ierr = PetscFEGeomDestroy(&cgeomFEM);CHKERRQ(ierr);
1304       }
1305     } else if (id == PETSCFV_CLASSID) {
1306       PetscInt       foff;
1307       PetscPointFunc obj_func;
1308       PetscScalar    lint;
1309 
1310       ierr = PetscDSGetObjective(prob, f, &obj_func);CHKERRQ(ierr);
1311       ierr = PetscDSGetFieldOffset(prob, f, &foff);CHKERRQ(ierr);
1312       if (obj_func) {
1313         for (c = 0; c < numCells; ++c) {
1314           PetscScalar *u_x;
1315 
1316           ierr = DMPlexPointLocalRead(dmGrad, c, lgrad, &u_x);CHKERRQ(ierr);
1317           obj_func(dim, Nf, NfAux, uOff, uOff_x, &u[totDim*c+foff], NULL, u_x, aOff, NULL, &a[totDimAux*c], NULL, NULL, 0.0, cgeomFVM[c].centroid, numConstants, constants, &lint);
1318           cintegral[c*Nf+f] += PetscRealPart(lint)*cgeomFVM[c].volume;
1319         }
1320       }
1321     } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
1322   }
1323   /* Cleanup data arrays */
1324   if (useFVM) {
1325     ierr = VecRestoreArrayRead(locGrad, &lgrad);CHKERRQ(ierr);
1326     ierr = VecRestoreArrayRead(cellGeometryFVM, (const PetscScalar **) &cgeomFVM);CHKERRQ(ierr);
1327     ierr = DMRestoreLocalVector(dmGrad, &locGrad);CHKERRQ(ierr);
1328     ierr = VecDestroy(&faceGeometryFVM);CHKERRQ(ierr);
1329     ierr = VecDestroy(&cellGeometryFVM);CHKERRQ(ierr);
1330     ierr = DMDestroy(&dmGrad);CHKERRQ(ierr);
1331   }
1332   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
1333   ierr = PetscFree(u);CHKERRQ(ierr);
1334   /* Cleanup */
1335   if (affineQuad) {
1336     ierr = PetscFEGeomDestroy(&cgeomFEM);CHKERRQ(ierr);
1337   }
1338   ierr = PetscQuadratureDestroy(&affineQuad);CHKERRQ(ierr);
1339   ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
1340   ierr = DMRestoreLocalVector(dm, &locX);CHKERRQ(ierr);
1341   PetscFunctionReturn(0);
1342 }
1343 
1344 /*@
1345   DMPlexComputeIntegralFEM - Form the integral over the domain from the global input X using pointwise functions specified by the user
1346 
1347   Input Parameters:
1348 + dm - The mesh
1349 . X  - Global input vector
1350 - user - The user context
1351 
1352   Output Parameter:
1353 . integral - Integral for each field
1354 
1355   Level: developer
1356 
1357 .seealso: DMPlexComputeResidualFEM()
1358 @*/
1359 PetscErrorCode DMPlexComputeIntegralFEM(DM dm, Vec X, PetscScalar *integral, void *user)
1360 {
1361   DM_Plex       *mesh = (DM_Plex *) dm->data;
1362   PetscScalar   *cintegral, *lintegral;
1363   PetscInt       Nf, f, cellHeight, cStart, cEnd, cEndInterior[4], cell;
1364   PetscErrorCode ierr;
1365 
1366   PetscFunctionBegin;
1367   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1368   PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
1369   PetscValidPointer(integral, 3);
1370   ierr = PetscLogEventBegin(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);
1371   ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr);
1372   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
1373   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
1374   ierr = DMPlexGetHybridBounds(dm, &cEndInterior[0], &cEndInterior[1], &cEndInterior[2], &cEndInterior[3]);CHKERRQ(ierr);
1375   cEnd = cEndInterior[cellHeight] < 0 ? cEnd : cEndInterior[cellHeight];
1376   /* TODO Introduce a loop over large chunks (right now this is a single chunk) */
1377   ierr = PetscCalloc2(Nf, &lintegral, (cEnd-cStart)*Nf, &cintegral);CHKERRQ(ierr);
1378   ierr = DMPlexComputeIntegral_Internal(dm, X, cStart, cEnd, cintegral, user);CHKERRQ(ierr);
1379   /* Sum up values */
1380   for (cell = cStart; cell < cEnd; ++cell) {
1381     const PetscInt c = cell - cStart;
1382 
1383     if (mesh->printFEM > 1) {ierr = DMPrintCellVector(cell, "Cell Integral", Nf, &cintegral[c*Nf]);CHKERRQ(ierr);}
1384     for (f = 0; f < Nf; ++f) lintegral[f] += cintegral[c*Nf+f];
1385   }
1386   ierr = MPIU_Allreduce(lintegral, integral, Nf, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject) dm));CHKERRQ(ierr);
1387   if (mesh->printFEM) {
1388     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Integral:");CHKERRQ(ierr);
1389     for (f = 0; f < Nf; ++f) {ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), " %g", (double) PetscRealPart(integral[f]));CHKERRQ(ierr);}
1390     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "\n");CHKERRQ(ierr);
1391   }
1392   ierr = PetscFree2(lintegral, cintegral);CHKERRQ(ierr);
1393   ierr = PetscLogEventEnd(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);
1394   PetscFunctionReturn(0);
1395 }
1396 
1397 /*@
1398   DMPlexComputeCellwiseIntegralFEM - Form the vector of cellwise integrals F from the global input X using pointwise functions specified by the user
1399 
1400   Input Parameters:
1401 + dm - The mesh
1402 . X  - Global input vector
1403 - user - The user context
1404 
1405   Output Parameter:
1406 . integral - Cellwise integrals for each field
1407 
1408   Level: developer
1409 
1410 .seealso: DMPlexComputeResidualFEM()
1411 @*/
1412 PetscErrorCode DMPlexComputeCellwiseIntegralFEM(DM dm, Vec X, Vec F, void *user)
1413 {
1414   DM_Plex       *mesh = (DM_Plex *) dm->data;
1415   DM             dmF;
1416   PetscSection   sectionF;
1417   PetscScalar   *cintegral, *af;
1418   PetscInt       Nf, f, cellHeight, cStart, cEnd, cEndInterior[4], cell;
1419   PetscErrorCode ierr;
1420 
1421   PetscFunctionBegin;
1422   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1423   PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
1424   PetscValidHeaderSpecific(F, VEC_CLASSID, 3);
1425   ierr = PetscLogEventBegin(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);
1426   ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr);
1427   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
1428   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
1429   ierr = DMPlexGetHybridBounds(dm, &cEndInterior[0], &cEndInterior[1], &cEndInterior[2], &cEndInterior[3]);CHKERRQ(ierr);
1430   cEnd = cEndInterior[cellHeight] < 0 ? cEnd : cEndInterior[cellHeight];
1431   /* TODO Introduce a loop over large chunks (right now this is a single chunk) */
1432   ierr = PetscCalloc1((cEnd-cStart)*Nf, &cintegral);CHKERRQ(ierr);
1433   ierr = DMPlexComputeIntegral_Internal(dm, X, cStart, cEnd, cintegral, user);CHKERRQ(ierr);
1434   /* Put values in F*/
1435   ierr = VecGetDM(F, &dmF);CHKERRQ(ierr);
1436   ierr = DMGetDefaultSection(dmF, &sectionF);CHKERRQ(ierr);
1437   ierr = VecGetArray(F, &af);CHKERRQ(ierr);
1438   for (cell = cStart; cell < cEnd; ++cell) {
1439     const PetscInt c = cell - cStart;
1440     PetscInt       dof, off;
1441 
1442     if (mesh->printFEM > 1) {ierr = DMPrintCellVector(cell, "Cell Integral", Nf, &cintegral[c*Nf]);CHKERRQ(ierr);}
1443     ierr = PetscSectionGetDof(sectionF, cell, &dof);CHKERRQ(ierr);
1444     ierr = PetscSectionGetOffset(sectionF, cell, &off);CHKERRQ(ierr);
1445     if (dof != Nf) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The number of cell dofs %D != %D", dof, Nf);
1446     for (f = 0; f < Nf; ++f) af[off+f] = cintegral[c*Nf+f];
1447   }
1448   ierr = VecRestoreArray(F, &af);CHKERRQ(ierr);
1449   ierr = PetscFree(cintegral);CHKERRQ(ierr);
1450   ierr = PetscLogEventEnd(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);
1451   PetscFunctionReturn(0);
1452 }
1453 
1454 /*@
1455   DMPlexComputeInterpolatorNested - Form the local portion of the interpolation matrix I from the coarse DM to the uniformly refined DM.
1456 
1457   Input Parameters:
1458 + dmf  - The fine mesh
1459 . dmc  - The coarse mesh
1460 - user - The user context
1461 
1462   Output Parameter:
1463 . In  - The interpolation matrix
1464 
1465   Level: developer
1466 
1467 .seealso: DMPlexComputeInterpolatorGeneral(), DMPlexComputeJacobianFEM()
1468 @*/
1469 PetscErrorCode DMPlexComputeInterpolatorNested(DM dmc, DM dmf, Mat In, void *user)
1470 {
1471   DM_Plex          *mesh  = (DM_Plex *) dmc->data;
1472   const char       *name  = "Interpolator";
1473   PetscDS           prob;
1474   PetscFE          *feRef;
1475   PetscFV          *fvRef;
1476   PetscSection      fsection, fglobalSection;
1477   PetscSection      csection, cglobalSection;
1478   PetscScalar      *elemMat;
1479   PetscInt          dim, Nf, f, fieldI, fieldJ, offsetI, offsetJ, cStart, cEnd, cEndInterior, c;
1480   PetscInt          cTotDim, rTotDim = 0;
1481   PetscErrorCode    ierr;
1482 
1483   PetscFunctionBegin;
1484   ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1485   ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr);
1486   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
1487   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
1488   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
1489   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
1490   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
1491   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
1492   ierr = DMPlexGetHybridBounds(dmc, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
1493   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
1494   ierr = DMGetDS(dmf, &prob);CHKERRQ(ierr);
1495   ierr = PetscCalloc2(Nf,&feRef,Nf,&fvRef);CHKERRQ(ierr);
1496   for (f = 0; f < Nf; ++f) {
1497     PetscObject  obj;
1498     PetscClassId id;
1499     PetscInt     rNb = 0, Nc = 0;
1500 
1501     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
1502     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1503     if (id == PETSCFE_CLASSID) {
1504       PetscFE fe = (PetscFE) obj;
1505 
1506       ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr);
1507       ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr);
1508       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1509     } else if (id == PETSCFV_CLASSID) {
1510       PetscFV        fv = (PetscFV) obj;
1511       PetscDualSpace Q;
1512 
1513       ierr = PetscFVRefine(fv, &fvRef[f]);CHKERRQ(ierr);
1514       ierr = PetscFVGetDualSpace(fvRef[f], &Q);CHKERRQ(ierr);
1515       ierr = PetscDualSpaceGetDimension(Q, &rNb);CHKERRQ(ierr);
1516       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
1517     }
1518     rTotDim += rNb;
1519   }
1520   ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr);
1521   ierr = PetscMalloc1(rTotDim*cTotDim,&elemMat);CHKERRQ(ierr);
1522   ierr = PetscMemzero(elemMat, rTotDim*cTotDim * sizeof(PetscScalar));CHKERRQ(ierr);
1523   for (fieldI = 0, offsetI = 0; fieldI < Nf; ++fieldI) {
1524     PetscDualSpace   Qref;
1525     PetscQuadrature  f;
1526     const PetscReal *qpoints, *qweights;
1527     PetscReal       *points;
1528     PetscInt         npoints = 0, Nc, Np, fpdim, i, k, p, d;
1529 
1530     /* Compose points from all dual basis functionals */
1531     if (feRef[fieldI]) {
1532       ierr = PetscFEGetDualSpace(feRef[fieldI], &Qref);CHKERRQ(ierr);
1533       ierr = PetscFEGetNumComponents(feRef[fieldI], &Nc);CHKERRQ(ierr);
1534     } else {
1535       ierr = PetscFVGetDualSpace(fvRef[fieldI], &Qref);CHKERRQ(ierr);
1536       ierr = PetscFVGetNumComponents(fvRef[fieldI], &Nc);CHKERRQ(ierr);
1537     }
1538     ierr = PetscDualSpaceGetDimension(Qref, &fpdim);CHKERRQ(ierr);
1539     for (i = 0; i < fpdim; ++i) {
1540       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1541       ierr = PetscQuadratureGetData(f, NULL, NULL, &Np, NULL, NULL);CHKERRQ(ierr);
1542       npoints += Np;
1543     }
1544     ierr = PetscMalloc1(npoints*dim,&points);CHKERRQ(ierr);
1545     for (i = 0, k = 0; i < fpdim; ++i) {
1546       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1547       ierr = PetscQuadratureGetData(f, NULL, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr);
1548       for (p = 0; p < Np; ++p, ++k) for (d = 0; d < dim; ++d) points[k*dim+d] = qpoints[p*dim+d];
1549     }
1550 
1551     for (fieldJ = 0, offsetJ = 0; fieldJ < Nf; ++fieldJ) {
1552       PetscObject  obj;
1553       PetscClassId id;
1554       PetscReal   *B;
1555       PetscInt     NcJ = 0, cpdim = 0, j, qNc;
1556 
1557       ierr = PetscDSGetDiscretization(prob, fieldJ, &obj);CHKERRQ(ierr);
1558       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1559       if (id == PETSCFE_CLASSID) {
1560         PetscFE fe = (PetscFE) obj;
1561 
1562         /* Evaluate basis at points */
1563         ierr = PetscFEGetNumComponents(fe, &NcJ);CHKERRQ(ierr);
1564         ierr = PetscFEGetDimension(fe, &cpdim);CHKERRQ(ierr);
1565         /* For now, fields only interpolate themselves */
1566         if (fieldI == fieldJ) {
1567           if (Nc != NcJ) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in fine space field %D does not match coarse field %D", Nc, NcJ);
1568           ierr = PetscFEGetTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);
1569           for (i = 0, k = 0; i < fpdim; ++i) {
1570             ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1571             ierr = PetscQuadratureGetData(f, NULL, &qNc, &Np, NULL, &qweights);CHKERRQ(ierr);
1572             if (qNc != NcJ) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in quadrature %D does not match coarse field %D", qNc, NcJ);
1573             for (p = 0; p < Np; ++p, ++k) {
1574               for (j = 0; j < cpdim; ++j) {
1575                 /*
1576                    cTotDim:            Total columns in element interpolation matrix, sum of number of dual basis functionals in each field
1577                    offsetI, offsetJ:   Offsets into the larger element interpolation matrix for different fields
1578                    fpdim, i, cpdim, j: Dofs for fine and coarse grids, correspond to dual space basis functionals
1579                    qNC, Nc, Ncj, c:    Number of components in this field
1580                    Np, p:              Number of quad points in the fine grid functional i
1581                    k:                  i*Np + p, overall point number for the interpolation
1582                 */
1583                 for (c = 0; c < Nc; ++c) elemMat[(offsetI + i)*cTotDim + offsetJ + j] += B[k*cpdim*NcJ+j*Nc+c]*qweights[p*qNc+c];
1584               }
1585             }
1586           }
1587           ierr = PetscFERestoreTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr);
1588         }
1589       } else if (id == PETSCFV_CLASSID) {
1590         PetscFV        fv = (PetscFV) obj;
1591 
1592         /* Evaluate constant function at points */
1593         ierr = PetscFVGetNumComponents(fv, &NcJ);CHKERRQ(ierr);
1594         cpdim = 1;
1595         /* For now, fields only interpolate themselves */
1596         if (fieldI == fieldJ) {
1597           if (Nc != NcJ) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in fine space field %d does not match coarse field %d", Nc, NcJ);
1598           for (i = 0, k = 0; i < fpdim; ++i) {
1599             ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1600             ierr = PetscQuadratureGetData(f, NULL, &qNc, &Np, NULL, &qweights);CHKERRQ(ierr);
1601             if (qNc != NcJ) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in quadrature %D does not match coarse field %D", qNc, NcJ);
1602             for (p = 0; p < Np; ++p, ++k) {
1603               for (j = 0; j < cpdim; ++j) {
1604                 for (c = 0; c < Nc; ++c) elemMat[(offsetI + i*Nc + c)*cTotDim + offsetJ + j*NcJ + c] += 1.0*qweights[p*qNc+c];
1605               }
1606             }
1607           }
1608         }
1609       }
1610       offsetJ += cpdim;
1611     }
1612     offsetI += fpdim;
1613     ierr = PetscFree(points);CHKERRQ(ierr);
1614   }
1615   if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(0, name, rTotDim, cTotDim, elemMat);CHKERRQ(ierr);}
1616   /* Preallocate matrix */
1617   {
1618     Mat          preallocator;
1619     PetscScalar *vals;
1620     PetscInt    *cellCIndices, *cellFIndices;
1621     PetscInt     locRows, locCols, cell;
1622 
1623     ierr = MatGetLocalSize(In, &locRows, &locCols);CHKERRQ(ierr);
1624     ierr = MatCreate(PetscObjectComm((PetscObject) In), &preallocator);CHKERRQ(ierr);
1625     ierr = MatSetType(preallocator, MATPREALLOCATOR);CHKERRQ(ierr);
1626     ierr = MatSetSizes(preallocator, locRows, locCols, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
1627     ierr = MatSetUp(preallocator);CHKERRQ(ierr);
1628     ierr = PetscCalloc3(rTotDim*cTotDim, &vals,cTotDim,&cellCIndices,rTotDim,&cellFIndices);CHKERRQ(ierr);
1629     for (cell = cStart; cell < cEnd; ++cell) {
1630       ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, cell, cellCIndices, cellFIndices);CHKERRQ(ierr);
1631       ierr = MatSetValues(preallocator, rTotDim, cellFIndices, cTotDim, cellCIndices, vals, INSERT_VALUES);CHKERRQ(ierr);
1632     }
1633     ierr = PetscFree3(vals,cellCIndices,cellFIndices);CHKERRQ(ierr);
1634     ierr = MatAssemblyBegin(preallocator, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1635     ierr = MatAssemblyEnd(preallocator, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1636     ierr = MatPreallocatorPreallocate(preallocator, PETSC_TRUE, In);CHKERRQ(ierr);
1637     ierr = MatDestroy(&preallocator);CHKERRQ(ierr);
1638   }
1639   /* Fill matrix */
1640   ierr = MatZeroEntries(In);CHKERRQ(ierr);
1641   for (c = cStart; c < cEnd; ++c) {
1642     ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES);CHKERRQ(ierr);
1643   }
1644   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);}
1645   ierr = PetscFree2(feRef,fvRef);CHKERRQ(ierr);
1646   ierr = PetscFree(elemMat);CHKERRQ(ierr);
1647   ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1648   ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1649   if (mesh->printFEM) {
1650     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
1651     ierr = MatChop(In, 1.0e-10);CHKERRQ(ierr);
1652     ierr = MatView(In, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1653   }
1654   ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1655   PetscFunctionReturn(0);
1656 }
1657 
1658 PetscErrorCode DMPlexComputeMassMatrixNested(DM dmc, DM dmf, Mat mass, void *user)
1659 {
1660   SETERRQ(PetscObjectComm((PetscObject) dmc), PETSC_ERR_SUP, "Laziness");
1661 }
1662 
1663 /*@
1664   DMPlexComputeInterpolatorGeneral - Form the local portion of the interpolation matrix I from the coarse DM to a non-nested fine DM.
1665 
1666   Input Parameters:
1667 + dmf  - The fine mesh
1668 . dmc  - The coarse mesh
1669 - user - The user context
1670 
1671   Output Parameter:
1672 . In  - The interpolation matrix
1673 
1674   Level: developer
1675 
1676 .seealso: DMPlexComputeInterpolatorNested(), DMPlexComputeJacobianFEM()
1677 @*/
1678 PetscErrorCode DMPlexComputeInterpolatorGeneral(DM dmc, DM dmf, Mat In, void *user)
1679 {
1680   DM_Plex       *mesh = (DM_Plex *) dmf->data;
1681   const char    *name = "Interpolator";
1682   PetscDS        prob;
1683   PetscSection   fsection, csection, globalFSection, globalCSection;
1684   PetscHashJK    ht;
1685   PetscLayout    rLayout;
1686   PetscInt      *dnz, *onz;
1687   PetscInt       locRows, rStart, rEnd;
1688   PetscReal     *x, *v0, *J, *invJ, detJ;
1689   PetscReal     *v0c, *Jc, *invJc, detJc;
1690   PetscScalar   *elemMat;
1691   PetscInt       dim, Nf, field, totDim, cStart, cEnd, cell, ccell;
1692   PetscErrorCode ierr;
1693 
1694   PetscFunctionBegin;
1695   ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1696   ierr = DMGetCoordinateDim(dmc, &dim);CHKERRQ(ierr);
1697   ierr = DMGetDS(dmc, &prob);CHKERRQ(ierr);
1698   ierr = PetscDSGetRefCoordArrays(prob, &x, NULL);CHKERRQ(ierr);
1699   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
1700   ierr = PetscMalloc3(dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
1701   ierr = PetscMalloc3(dim,&v0c,dim*dim,&Jc,dim*dim,&invJc);CHKERRQ(ierr);
1702   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
1703   ierr = DMGetDefaultGlobalSection(dmf, &globalFSection);CHKERRQ(ierr);
1704   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
1705   ierr = DMGetDefaultGlobalSection(dmc, &globalCSection);CHKERRQ(ierr);
1706   ierr = DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd);CHKERRQ(ierr);
1707   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
1708   ierr = PetscMalloc1(totDim, &elemMat);CHKERRQ(ierr);
1709 
1710   ierr = MatGetLocalSize(In, &locRows, NULL);CHKERRQ(ierr);
1711   ierr = PetscLayoutCreate(PetscObjectComm((PetscObject) In), &rLayout);CHKERRQ(ierr);
1712   ierr = PetscLayoutSetLocalSize(rLayout, locRows);CHKERRQ(ierr);
1713   ierr = PetscLayoutSetBlockSize(rLayout, 1);CHKERRQ(ierr);
1714   ierr = PetscLayoutSetUp(rLayout);CHKERRQ(ierr);
1715   ierr = PetscLayoutGetRange(rLayout, &rStart, &rEnd);CHKERRQ(ierr);
1716   ierr = PetscLayoutDestroy(&rLayout);CHKERRQ(ierr);
1717   ierr = PetscCalloc2(locRows,&dnz,locRows,&onz);CHKERRQ(ierr);
1718   ierr = PetscHashJKCreate(&ht);CHKERRQ(ierr);
1719   for (field = 0; field < Nf; ++field) {
1720     PetscObject      obj;
1721     PetscClassId     id;
1722     PetscDualSpace   Q = NULL;
1723     PetscQuadrature  f;
1724     const PetscReal *qpoints;
1725     PetscInt         Nc, Np, fpdim, i, d;
1726 
1727     ierr = PetscDSGetDiscretization(prob, field, &obj);CHKERRQ(ierr);
1728     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1729     if (id == PETSCFE_CLASSID) {
1730       PetscFE fe = (PetscFE) obj;
1731 
1732       ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr);
1733       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1734     } else if (id == PETSCFV_CLASSID) {
1735       PetscFV fv = (PetscFV) obj;
1736 
1737       ierr = PetscFVGetDualSpace(fv, &Q);CHKERRQ(ierr);
1738       Nc   = 1;
1739     }
1740     ierr = PetscDualSpaceGetDimension(Q, &fpdim);CHKERRQ(ierr);
1741     /* For each fine grid cell */
1742     for (cell = cStart; cell < cEnd; ++cell) {
1743       PetscInt *findices,   *cindices;
1744       PetscInt  numFIndices, numCIndices;
1745 
1746       ierr = DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr);
1747       ierr = DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
1748       if (numFIndices != fpdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of fine indices %d != %d dual basis vecs", numFIndices, fpdim);
1749       for (i = 0; i < fpdim; ++i) {
1750         Vec             pointVec;
1751         PetscScalar    *pV;
1752         PetscSF         coarseCellSF = NULL;
1753         const PetscSFNode *coarseCells;
1754         PetscInt        numCoarseCells, q, c;
1755 
1756         /* Get points from the dual basis functional quadrature */
1757         ierr = PetscDualSpaceGetFunctional(Q, i, &f);CHKERRQ(ierr);
1758         ierr = PetscQuadratureGetData(f, NULL, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr);
1759         ierr = VecCreateSeq(PETSC_COMM_SELF, Np*dim, &pointVec);CHKERRQ(ierr);
1760         ierr = VecSetBlockSize(pointVec, dim);CHKERRQ(ierr);
1761         ierr = VecGetArray(pointVec, &pV);CHKERRQ(ierr);
1762         for (q = 0; q < Np; ++q) {
1763           const PetscReal xi0[3] = {-1., -1., -1.};
1764 
1765           /* Transform point to real space */
1766           CoordinatesRefToReal(dim, dim, xi0, v0, J, &qpoints[q*dim], x);
1767           for (d = 0; d < dim; ++d) pV[q*dim+d] = x[d];
1768         }
1769         ierr = VecRestoreArray(pointVec, &pV);CHKERRQ(ierr);
1770         /* Get set of coarse cells that overlap points (would like to group points by coarse cell) */
1771         /* OPT: Pack all quad points from fine cell */
1772         ierr = DMLocatePoints(dmc, pointVec, DM_POINTLOCATION_NEAREST, &coarseCellSF);CHKERRQ(ierr);
1773         ierr = PetscSFViewFromOptions(coarseCellSF, NULL, "-interp_sf_view");CHKERRQ(ierr);
1774         /* Update preallocation info */
1775         ierr = PetscSFGetGraph(coarseCellSF, NULL, &numCoarseCells, NULL, &coarseCells);CHKERRQ(ierr);
1776         if (numCoarseCells != Np) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Not all closure points located");
1777         {
1778           PetscHashJKKey  key;
1779           PetscHashJKIter missing, iter;
1780 
1781           key.j = findices[i];
1782           if (key.j >= 0) {
1783             /* Get indices for coarse elements */
1784             for (ccell = 0; ccell < numCoarseCells; ++ccell) {
1785               ierr = DMPlexGetClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);CHKERRQ(ierr);
1786               for (c = 0; c < numCIndices; ++c) {
1787                 key.k = cindices[c];
1788                 if (key.k < 0) continue;
1789                 ierr = PetscHashJKPut(ht, key, &missing, &iter);CHKERRQ(ierr);
1790                 if (missing) {
1791                   ierr = PetscHashJKSet(ht, iter, 1);CHKERRQ(ierr);
1792                   if ((key.k >= rStart) && (key.k < rEnd)) ++dnz[key.j-rStart];
1793                   else                                     ++onz[key.j-rStart];
1794                 }
1795               }
1796               ierr = DMPlexRestoreClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);CHKERRQ(ierr);
1797             }
1798           }
1799         }
1800         ierr = PetscSFDestroy(&coarseCellSF);CHKERRQ(ierr);
1801         ierr = VecDestroy(&pointVec);CHKERRQ(ierr);
1802       }
1803       ierr = DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr);
1804     }
1805   }
1806   ierr = PetscHashJKDestroy(&ht);CHKERRQ(ierr);
1807   ierr = MatXAIJSetPreallocation(In, 1, dnz, onz, NULL, NULL);CHKERRQ(ierr);
1808   ierr = MatSetOption(In, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
1809   ierr = PetscFree2(dnz,onz);CHKERRQ(ierr);
1810   for (field = 0; field < Nf; ++field) {
1811     PetscObject      obj;
1812     PetscClassId     id;
1813     PetscDualSpace   Q = NULL;
1814     PetscQuadrature  f;
1815     const PetscReal *qpoints, *qweights;
1816     PetscInt         Nc, qNc, Np, fpdim, i, d;
1817 
1818     ierr = PetscDSGetDiscretization(prob, field, &obj);CHKERRQ(ierr);
1819     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1820     if (id == PETSCFE_CLASSID) {
1821       PetscFE fe = (PetscFE) obj;
1822 
1823       ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr);
1824       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1825     } else if (id == PETSCFV_CLASSID) {
1826       PetscFV fv = (PetscFV) obj;
1827 
1828       ierr = PetscFVGetDualSpace(fv, &Q);CHKERRQ(ierr);
1829       Nc   = 1;
1830     } else SETERRQ1(PetscObjectComm((PetscObject)dmc),PETSC_ERR_ARG_WRONG,"Unknown discretization type for field %d",field);
1831     ierr = PetscDualSpaceGetDimension(Q, &fpdim);CHKERRQ(ierr);
1832     /* For each fine grid cell */
1833     for (cell = cStart; cell < cEnd; ++cell) {
1834       PetscInt *findices,   *cindices;
1835       PetscInt  numFIndices, numCIndices;
1836 
1837       ierr = DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr);
1838       ierr = DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
1839       if (numFIndices != fpdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of fine indices %d != %d dual basis vecs", numFIndices, fpdim);
1840       for (i = 0; i < fpdim; ++i) {
1841         Vec             pointVec;
1842         PetscScalar    *pV;
1843         PetscSF         coarseCellSF = NULL;
1844         const PetscSFNode *coarseCells;
1845         PetscInt        numCoarseCells, cpdim, q, c, j;
1846 
1847         /* Get points from the dual basis functional quadrature */
1848         ierr = PetscDualSpaceGetFunctional(Q, i, &f);CHKERRQ(ierr);
1849         ierr = PetscQuadratureGetData(f, NULL, &qNc, &Np, &qpoints, &qweights);CHKERRQ(ierr);
1850         if (qNc != Nc) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in quadrature %D does not match coarse field %D", qNc, Nc);
1851         ierr = VecCreateSeq(PETSC_COMM_SELF, Np*dim, &pointVec);CHKERRQ(ierr);
1852         ierr = VecSetBlockSize(pointVec, dim);CHKERRQ(ierr);
1853         ierr = VecGetArray(pointVec, &pV);CHKERRQ(ierr);
1854         for (q = 0; q < Np; ++q) {
1855           const PetscReal xi0[3] = {-1., -1., -1.};
1856 
1857           /* Transform point to real space */
1858           CoordinatesRefToReal(dim, dim, xi0, v0, J, &qpoints[q*dim], x);
1859           for (d = 0; d < dim; ++d) pV[q*dim+d] = x[d];
1860         }
1861         ierr = VecRestoreArray(pointVec, &pV);CHKERRQ(ierr);
1862         /* Get set of coarse cells that overlap points (would like to group points by coarse cell) */
1863         /* OPT: Read this out from preallocation information */
1864         ierr = DMLocatePoints(dmc, pointVec, DM_POINTLOCATION_NEAREST, &coarseCellSF);CHKERRQ(ierr);
1865         /* Update preallocation info */
1866         ierr = PetscSFGetGraph(coarseCellSF, NULL, &numCoarseCells, NULL, &coarseCells);CHKERRQ(ierr);
1867         if (numCoarseCells != Np) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Not all closure points located");
1868         ierr = VecGetArray(pointVec, &pV);CHKERRQ(ierr);
1869         for (ccell = 0; ccell < numCoarseCells; ++ccell) {
1870           PetscReal pVReal[3];
1871           const PetscReal xi0[3] = {-1., -1., -1.};
1872 
1873           ierr = DMPlexGetClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);CHKERRQ(ierr);
1874           /* Transform points from real space to coarse reference space */
1875           ierr = DMPlexComputeCellGeometryFEM(dmc, coarseCells[ccell].index, NULL, v0c, Jc, invJc, &detJc);CHKERRQ(ierr);
1876           for (d = 0; d < dim; ++d) pVReal[d] = PetscRealPart(pV[ccell*dim+d]);
1877           CoordinatesRealToRef(dim, dim, xi0, v0c, invJc, pVReal, x);
1878 
1879           if (id == PETSCFE_CLASSID) {
1880             PetscFE    fe = (PetscFE) obj;
1881             PetscReal *B;
1882 
1883             /* Evaluate coarse basis on contained point */
1884             ierr = PetscFEGetDimension(fe, &cpdim);CHKERRQ(ierr);
1885             ierr = PetscFEGetTabulation(fe, 1, x, &B, NULL, NULL);CHKERRQ(ierr);
1886             ierr = PetscMemzero(elemMat, cpdim * sizeof(PetscScalar));CHKERRQ(ierr);
1887             /* Get elemMat entries by multiplying by weight */
1888             for (j = 0; j < cpdim; ++j) {
1889               for (c = 0; c < Nc; ++c) elemMat[j] += B[j*Nc + c]*qweights[ccell*qNc + c];
1890             }
1891             ierr = PetscFERestoreTabulation(fe, 1, x, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr);
1892           } else {
1893             cpdim = 1;
1894             for (j = 0; j < cpdim; ++j) {
1895               for (c = 0; c < Nc; ++c) elemMat[j] += 1.0*qweights[ccell*qNc + c];
1896             }
1897           }
1898           /* Update interpolator */
1899           if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat);CHKERRQ(ierr);}
1900           if (numCIndices != cpdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of element matrix columns %D != %D", numCIndices, cpdim);
1901           ierr = MatSetValues(In, 1, &findices[i], numCIndices, cindices, elemMat, INSERT_VALUES);CHKERRQ(ierr);
1902           ierr = DMPlexRestoreClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);CHKERRQ(ierr);
1903         }
1904         ierr = VecRestoreArray(pointVec, &pV);CHKERRQ(ierr);
1905         ierr = PetscSFDestroy(&coarseCellSF);CHKERRQ(ierr);
1906         ierr = VecDestroy(&pointVec);CHKERRQ(ierr);
1907       }
1908       ierr = DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr);
1909     }
1910   }
1911   ierr = PetscFree3(v0,J,invJ);CHKERRQ(ierr);
1912   ierr = PetscFree3(v0c,Jc,invJc);CHKERRQ(ierr);
1913   ierr = PetscFree(elemMat);CHKERRQ(ierr);
1914   ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1915   ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1916   ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1917   PetscFunctionReturn(0);
1918 }
1919 
1920 /*@
1921   DMPlexComputeMassMatrixGeneral - Form the local portion of the mass matrix M from the coarse DM to a non-nested fine DM.
1922 
1923   Input Parameters:
1924 + dmf  - The fine mesh
1925 . dmc  - The coarse mesh
1926 - user - The user context
1927 
1928   Output Parameter:
1929 . mass  - The mass matrix
1930 
1931   Level: developer
1932 
1933 .seealso: DMPlexComputeMassMatrixNested(), DMPlexComputeInterpolatorNested(), DMPlexComputeInterpolatorGeneral(), DMPlexComputeJacobianFEM()
1934 @*/
1935 PetscErrorCode DMPlexComputeMassMatrixGeneral(DM dmc, DM dmf, Mat mass, void *user)
1936 {
1937   DM_Plex       *mesh = (DM_Plex *) dmf->data;
1938   const char    *name = "Mass Matrix";
1939   PetscDS        prob;
1940   PetscSection   fsection, csection, globalFSection, globalCSection;
1941   PetscHashJK    ht;
1942   PetscLayout    rLayout;
1943   PetscInt      *dnz, *onz;
1944   PetscInt       locRows, rStart, rEnd;
1945   PetscReal     *x, *v0, *J, *invJ, detJ;
1946   PetscReal     *v0c, *Jc, *invJc, detJc;
1947   PetscScalar   *elemMat;
1948   PetscInt       dim, Nf, field, totDim, cStart, cEnd, cell, ccell;
1949   PetscErrorCode ierr;
1950 
1951   PetscFunctionBegin;
1952   ierr = DMGetCoordinateDim(dmc, &dim);CHKERRQ(ierr);
1953   ierr = DMGetDS(dmc, &prob);CHKERRQ(ierr);
1954   ierr = PetscDSGetRefCoordArrays(prob, &x, NULL);CHKERRQ(ierr);
1955   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
1956   ierr = PetscMalloc3(dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
1957   ierr = PetscMalloc3(dim,&v0c,dim*dim,&Jc,dim*dim,&invJc);CHKERRQ(ierr);
1958   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
1959   ierr = DMGetDefaultGlobalSection(dmf, &globalFSection);CHKERRQ(ierr);
1960   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
1961   ierr = DMGetDefaultGlobalSection(dmc, &globalCSection);CHKERRQ(ierr);
1962   ierr = DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd);CHKERRQ(ierr);
1963   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
1964   ierr = PetscMalloc1(totDim, &elemMat);CHKERRQ(ierr);
1965 
1966   ierr = MatGetLocalSize(mass, &locRows, NULL);CHKERRQ(ierr);
1967   ierr = PetscLayoutCreate(PetscObjectComm((PetscObject) mass), &rLayout);CHKERRQ(ierr);
1968   ierr = PetscLayoutSetLocalSize(rLayout, locRows);CHKERRQ(ierr);
1969   ierr = PetscLayoutSetBlockSize(rLayout, 1);CHKERRQ(ierr);
1970   ierr = PetscLayoutSetUp(rLayout);CHKERRQ(ierr);
1971   ierr = PetscLayoutGetRange(rLayout, &rStart, &rEnd);CHKERRQ(ierr);
1972   ierr = PetscLayoutDestroy(&rLayout);CHKERRQ(ierr);
1973   ierr = PetscCalloc2(locRows,&dnz,locRows,&onz);CHKERRQ(ierr);
1974   ierr = PetscHashJKCreate(&ht);CHKERRQ(ierr);
1975   for (field = 0; field < Nf; ++field) {
1976     PetscObject      obj;
1977     PetscClassId     id;
1978     PetscQuadrature  quad;
1979     const PetscReal *qpoints;
1980     PetscInt         Nq, Nc, i, d;
1981 
1982     ierr = PetscDSGetDiscretization(prob, field, &obj);CHKERRQ(ierr);
1983     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1984     if (id == PETSCFE_CLASSID) {ierr = PetscFEGetQuadrature((PetscFE) obj, &quad);CHKERRQ(ierr);}
1985     else                       {ierr = PetscFVGetQuadrature((PetscFV) obj, &quad);CHKERRQ(ierr);}
1986     ierr = PetscQuadratureGetData(quad, NULL, &Nc, &Nq, &qpoints, NULL);CHKERRQ(ierr);
1987     /* For each fine grid cell */
1988     for (cell = cStart; cell < cEnd; ++cell) {
1989       Vec                pointVec;
1990       PetscScalar       *pV;
1991       PetscSF            coarseCellSF = NULL;
1992       const PetscSFNode *coarseCells;
1993       PetscInt           numCoarseCells, q, c;
1994       PetscInt          *findices,   *cindices;
1995       PetscInt           numFIndices, numCIndices;
1996 
1997       ierr = DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr);
1998       ierr = DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
1999       /* Get points from the quadrature */
2000       ierr = VecCreateSeq(PETSC_COMM_SELF, Nq*dim, &pointVec);CHKERRQ(ierr);
2001       ierr = VecSetBlockSize(pointVec, dim);CHKERRQ(ierr);
2002       ierr = VecGetArray(pointVec, &pV);CHKERRQ(ierr);
2003       for (q = 0; q < Nq; ++q) {
2004         const PetscReal xi0[3] = {-1., -1., -1.};
2005 
2006         /* Transform point to real space */
2007         CoordinatesRefToReal(dim, dim, xi0, v0, J, &qpoints[q*dim], x);
2008         for (d = 0; d < dim; ++d) pV[q*dim+d] = x[d];
2009       }
2010       ierr = VecRestoreArray(pointVec, &pV);CHKERRQ(ierr);
2011       /* Get set of coarse cells that overlap points (would like to group points by coarse cell) */
2012       ierr = DMLocatePoints(dmc, pointVec, DM_POINTLOCATION_NEAREST, &coarseCellSF);CHKERRQ(ierr);
2013       ierr = PetscSFViewFromOptions(coarseCellSF, NULL, "-interp_sf_view");CHKERRQ(ierr);
2014       /* Update preallocation info */
2015       ierr = PetscSFGetGraph(coarseCellSF, NULL, &numCoarseCells, NULL, &coarseCells);CHKERRQ(ierr);
2016       if (numCoarseCells != Nq) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Not all closure points located");
2017       {
2018         PetscHashJKKey  key;
2019         PetscHashJKIter missing, iter;
2020 
2021         for (i = 0; i < numFIndices; ++i) {
2022           key.j = findices[i];
2023           if (key.j >= 0) {
2024             /* Get indices for coarse elements */
2025             for (ccell = 0; ccell < numCoarseCells; ++ccell) {
2026               ierr = DMPlexGetClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);CHKERRQ(ierr);
2027               for (c = 0; c < numCIndices; ++c) {
2028                 key.k = cindices[c];
2029                 if (key.k < 0) continue;
2030                 ierr = PetscHashJKPut(ht, key, &missing, &iter);CHKERRQ(ierr);
2031                 if (missing) {
2032                   ierr = PetscHashJKSet(ht, iter, 1);CHKERRQ(ierr);
2033                   if ((key.k >= rStart) && (key.k < rEnd)) ++dnz[key.j-rStart];
2034                   else                                     ++onz[key.j-rStart];
2035                 }
2036               }
2037               ierr = DMPlexRestoreClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);CHKERRQ(ierr);
2038             }
2039           }
2040         }
2041       }
2042       ierr = PetscSFDestroy(&coarseCellSF);CHKERRQ(ierr);
2043       ierr = VecDestroy(&pointVec);CHKERRQ(ierr);
2044       ierr = DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr);
2045     }
2046   }
2047   ierr = PetscHashJKDestroy(&ht);CHKERRQ(ierr);
2048   ierr = MatXAIJSetPreallocation(mass, 1, dnz, onz, NULL, NULL);CHKERRQ(ierr);
2049   ierr = MatSetOption(mass, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
2050   ierr = PetscFree2(dnz,onz);CHKERRQ(ierr);
2051   for (field = 0; field < Nf; ++field) {
2052     PetscObject      obj;
2053     PetscClassId     id;
2054     PetscQuadrature  quad;
2055     PetscReal       *Bfine;
2056     const PetscReal *qpoints, *qweights;
2057     PetscInt         Nq, Nc, i, d;
2058 
2059     ierr = PetscDSGetDiscretization(prob, field, &obj);CHKERRQ(ierr);
2060     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
2061     if (id == PETSCFE_CLASSID) {ierr = PetscFEGetQuadrature((PetscFE) obj, &quad);CHKERRQ(ierr);ierr = PetscFEGetDefaultTabulation((PetscFE) obj, &Bfine, NULL, NULL);CHKERRQ(ierr);}
2062     else                       {ierr = PetscFVGetQuadrature((PetscFV) obj, &quad);CHKERRQ(ierr);}
2063     ierr = PetscQuadratureGetData(quad, NULL, &Nc, &Nq, &qpoints, &qweights);CHKERRQ(ierr);
2064     /* For each fine grid cell */
2065     for (cell = cStart; cell < cEnd; ++cell) {
2066       Vec                pointVec;
2067       PetscScalar       *pV;
2068       PetscSF            coarseCellSF = NULL;
2069       const PetscSFNode *coarseCells;
2070       PetscInt           numCoarseCells, cpdim, q, c, j;
2071       PetscInt          *findices,   *cindices;
2072       PetscInt           numFIndices, numCIndices;
2073 
2074       ierr = DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr);
2075       ierr = DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
2076       /* Get points from the quadrature */
2077       ierr = VecCreateSeq(PETSC_COMM_SELF, Nq*dim, &pointVec);CHKERRQ(ierr);
2078       ierr = VecSetBlockSize(pointVec, dim);CHKERRQ(ierr);
2079       ierr = VecGetArray(pointVec, &pV);CHKERRQ(ierr);
2080       for (q = 0; q < Nq; ++q) {
2081         const PetscReal xi0[3] = {-1., -1., -1.};
2082 
2083         /* Transform point to real space */
2084         CoordinatesRefToReal(dim, dim, xi0, v0, J, &qpoints[q*dim], x);
2085         for (d = 0; d < dim; ++d) pV[q*dim+d] = x[d];
2086       }
2087       ierr = VecRestoreArray(pointVec, &pV);CHKERRQ(ierr);
2088       /* Get set of coarse cells that overlap points (would like to group points by coarse cell) */
2089       ierr = DMLocatePoints(dmc, pointVec, DM_POINTLOCATION_NEAREST, &coarseCellSF);CHKERRQ(ierr);
2090       /* Update matrix */
2091       ierr = PetscSFGetGraph(coarseCellSF, NULL, &numCoarseCells, NULL, &coarseCells);CHKERRQ(ierr);
2092       if (numCoarseCells != Nq) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Not all closure points located");
2093       ierr = VecGetArray(pointVec, &pV);CHKERRQ(ierr);
2094       for (ccell = 0; ccell < numCoarseCells; ++ccell) {
2095         PetscReal pVReal[3];
2096         const PetscReal xi0[3] = {-1., -1., -1.};
2097 
2098 
2099         ierr = DMPlexGetClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);CHKERRQ(ierr);
2100         /* Transform points from real space to coarse reference space */
2101         ierr = DMPlexComputeCellGeometryFEM(dmc, coarseCells[ccell].index, NULL, v0c, Jc, invJc, &detJc);CHKERRQ(ierr);
2102         for (d = 0; d < dim; ++d) pVReal[d] = PetscRealPart(pV[ccell*dim+d]);
2103         CoordinatesRealToRef(dim, dim, xi0, v0c, invJc, pVReal, x);
2104 
2105         if (id == PETSCFE_CLASSID) {
2106           PetscFE    fe = (PetscFE) obj;
2107           PetscReal *B;
2108 
2109           /* Evaluate coarse basis on contained point */
2110           ierr = PetscFEGetDimension(fe, &cpdim);CHKERRQ(ierr);
2111           ierr = PetscFEGetTabulation(fe, 1, x, &B, NULL, NULL);CHKERRQ(ierr);
2112           /* Get elemMat entries by multiplying by weight */
2113           for (i = 0; i < numFIndices; ++i) {
2114             ierr = PetscMemzero(elemMat, cpdim * sizeof(PetscScalar));CHKERRQ(ierr);
2115             for (j = 0; j < cpdim; ++j) {
2116               for (c = 0; c < Nc; ++c) elemMat[j] += B[j*Nc + c]*Bfine[(ccell*numFIndices + i)*Nc + c]*qweights[ccell*Nc + c]*detJ;
2117             }
2118             /* Update interpolator */
2119             if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat);CHKERRQ(ierr);}
2120             if (numCIndices != cpdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of element matrix columns %D != %D", numCIndices, cpdim);
2121             ierr = MatSetValues(mass, 1, &findices[i], numCIndices, cindices, elemMat, ADD_VALUES);CHKERRQ(ierr);
2122           }
2123           ierr = PetscFERestoreTabulation(fe, 1, x, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr);
2124         } else {
2125           cpdim = 1;
2126           for (i = 0; i < numFIndices; ++i) {
2127             ierr = PetscMemzero(elemMat, cpdim * sizeof(PetscScalar));CHKERRQ(ierr);
2128             for (j = 0; j < cpdim; ++j) {
2129               for (c = 0; c < Nc; ++c) elemMat[j] += 1.0*1.0*qweights[ccell*Nc + c]*detJ;
2130             }
2131             /* Update interpolator */
2132             if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat);CHKERRQ(ierr);}
2133             ierr = PetscPrintf(PETSC_COMM_SELF, "Nq: %d %d Nf: %d %d Nc: %d %d\n", ccell, Nq, i, numFIndices, j, numCIndices);CHKERRQ(ierr);
2134             if (numCIndices != cpdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of element matrix columns %D != %D", numCIndices, cpdim);
2135             ierr = MatSetValues(mass, 1, &findices[i], numCIndices, cindices, elemMat, ADD_VALUES);CHKERRQ(ierr);
2136           }
2137         }
2138         ierr = DMPlexRestoreClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);CHKERRQ(ierr);
2139       }
2140       ierr = VecRestoreArray(pointVec, &pV);CHKERRQ(ierr);
2141       ierr = PetscSFDestroy(&coarseCellSF);CHKERRQ(ierr);
2142       ierr = VecDestroy(&pointVec);CHKERRQ(ierr);
2143       ierr = DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr);
2144     }
2145   }
2146   ierr = PetscFree3(v0,J,invJ);CHKERRQ(ierr);
2147   ierr = PetscFree3(v0c,Jc,invJc);CHKERRQ(ierr);
2148   ierr = PetscFree(elemMat);CHKERRQ(ierr);
2149   ierr = MatAssemblyBegin(mass, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2150   ierr = MatAssemblyEnd(mass, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2151   PetscFunctionReturn(0);
2152 }
2153 
2154 /*@
2155   DMPlexComputeInjectorFEM - Compute a mapping from coarse unknowns to fine unknowns
2156 
2157   Input Parameters:
2158 + dmc  - The coarse mesh
2159 - dmf  - The fine mesh
2160 - user - The user context
2161 
2162   Output Parameter:
2163 . sc   - The mapping
2164 
2165   Level: developer
2166 
2167 .seealso: DMPlexComputeInterpolatorNested(), DMPlexComputeJacobianFEM()
2168 @*/
2169 PetscErrorCode DMPlexComputeInjectorFEM(DM dmc, DM dmf, VecScatter *sc, void *user)
2170 {
2171   PetscDS        prob;
2172   PetscFE       *feRef;
2173   PetscFV       *fvRef;
2174   Vec            fv, cv;
2175   IS             fis, cis;
2176   PetscSection   fsection, fglobalSection, csection, cglobalSection;
2177   PetscInt      *cmap, *cellCIndices, *cellFIndices, *cindices, *findices;
2178   PetscInt       cTotDim, fTotDim = 0, Nf, f, field, cStart, cEnd, cEndInterior, c, dim, d, startC, endC, offsetC, offsetF, m;
2179   PetscErrorCode ierr;
2180 
2181   PetscFunctionBegin;
2182   ierr = PetscLogEventBegin(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
2183   ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr);
2184   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
2185   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
2186   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
2187   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
2188   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
2189   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
2190   ierr = DMPlexGetHybridBounds(dmc, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
2191   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
2192   ierr = DMGetDS(dmc, &prob);CHKERRQ(ierr);
2193   ierr = PetscCalloc2(Nf,&feRef,Nf,&fvRef);CHKERRQ(ierr);
2194   for (f = 0; f < Nf; ++f) {
2195     PetscObject  obj;
2196     PetscClassId id;
2197     PetscInt     fNb = 0, Nc = 0;
2198 
2199     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
2200     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
2201     if (id == PETSCFE_CLASSID) {
2202       PetscFE fe = (PetscFE) obj;
2203 
2204       ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr);
2205       ierr = PetscFEGetDimension(feRef[f], &fNb);CHKERRQ(ierr);
2206       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
2207     } else if (id == PETSCFV_CLASSID) {
2208       PetscFV        fv = (PetscFV) obj;
2209       PetscDualSpace Q;
2210 
2211       ierr = PetscFVRefine(fv, &fvRef[f]);CHKERRQ(ierr);
2212       ierr = PetscFVGetDualSpace(fvRef[f], &Q);CHKERRQ(ierr);
2213       ierr = PetscDualSpaceGetDimension(Q, &fNb);CHKERRQ(ierr);
2214       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
2215     }
2216     fTotDim += fNb;
2217   }
2218   ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr);
2219   ierr = PetscMalloc1(cTotDim,&cmap);CHKERRQ(ierr);
2220   for (field = 0, offsetC = 0, offsetF = 0; field < Nf; ++field) {
2221     PetscFE        feC;
2222     PetscFV        fvC;
2223     PetscDualSpace QF, QC;
2224     PetscInt       order = -1, NcF, NcC, fpdim, cpdim;
2225 
2226     if (feRef[field]) {
2227       ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &feC);CHKERRQ(ierr);
2228       ierr = PetscFEGetNumComponents(feC, &NcC);CHKERRQ(ierr);
2229       ierr = PetscFEGetNumComponents(feRef[field], &NcF);CHKERRQ(ierr);
2230       ierr = PetscFEGetDualSpace(feRef[field], &QF);CHKERRQ(ierr);
2231       ierr = PetscDualSpaceGetOrder(QF, &order);CHKERRQ(ierr);
2232       ierr = PetscDualSpaceGetDimension(QF, &fpdim);CHKERRQ(ierr);
2233       ierr = PetscFEGetDualSpace(feC, &QC);CHKERRQ(ierr);
2234       ierr = PetscDualSpaceGetDimension(QC, &cpdim);CHKERRQ(ierr);
2235     } else {
2236       ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fvC);CHKERRQ(ierr);
2237       ierr = PetscFVGetNumComponents(fvC, &NcC);CHKERRQ(ierr);
2238       ierr = PetscFVGetNumComponents(fvRef[field], &NcF);CHKERRQ(ierr);
2239       ierr = PetscFVGetDualSpace(fvRef[field], &QF);CHKERRQ(ierr);
2240       ierr = PetscDualSpaceGetDimension(QF, &fpdim);CHKERRQ(ierr);
2241       ierr = PetscFVGetDualSpace(fvC, &QC);CHKERRQ(ierr);
2242       ierr = PetscDualSpaceGetDimension(QC, &cpdim);CHKERRQ(ierr);
2243     }
2244     if (NcF != NcC) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in fine space field %d does not match coarse field %d", NcF, NcC);
2245     for (c = 0; c < cpdim; ++c) {
2246       PetscQuadrature  cfunc;
2247       const PetscReal *cqpoints, *cqweights;
2248       PetscInt         NqcC, NpC;
2249       PetscBool        found = PETSC_FALSE;
2250 
2251       ierr = PetscDualSpaceGetFunctional(QC, c, &cfunc);CHKERRQ(ierr);
2252       ierr = PetscQuadratureGetData(cfunc, NULL, &NqcC, &NpC, &cqpoints, &cqweights);CHKERRQ(ierr);
2253       if (NqcC != NcC) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of quadrature components %D must match number of field components", NqcC, NcC);
2254       if (NpC != 1 && feRef[field]) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not know how to do injection for moments");
2255       for (f = 0; f < fpdim; ++f) {
2256         PetscQuadrature  ffunc;
2257         const PetscReal *fqpoints, *fqweights;
2258         PetscReal        sum = 0.0;
2259         PetscInt         NqcF, NpF;
2260 
2261         ierr = PetscDualSpaceGetFunctional(QF, f, &ffunc);CHKERRQ(ierr);
2262         ierr = PetscQuadratureGetData(ffunc, NULL, &NqcF, &NpF, &fqpoints, &fqweights);CHKERRQ(ierr);
2263         if (NqcF != NcF) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of quadrature components %D must match number of field components", NqcF, NcF);
2264         if (NpC != NpF) continue;
2265         for (d = 0; d < dim; ++d) sum += PetscAbsReal(cqpoints[d] - fqpoints[d]);
2266         if (sum > 1.0e-9) continue;
2267         for (d = 0; d < NcC; ++d) sum += PetscAbsReal(cqweights[d]*fqweights[d]);
2268         if (sum < 1.0e-9) continue;
2269         cmap[offsetC+c] = offsetF+f;
2270         found = PETSC_TRUE;
2271         break;
2272       }
2273       if (!found) {
2274         /* TODO We really want the average here, but some asshole put VecScatter in the interface */
2275         if (fvRef[field] || (feRef[field] && order == 0)) {
2276           cmap[offsetC+c] = offsetF+0;
2277         } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not locate matching functional for injection");
2278       }
2279     }
2280     offsetC += cpdim;
2281     offsetF += fpdim;
2282   }
2283   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);ierr = PetscFVDestroy(&fvRef[f]);CHKERRQ(ierr);}
2284   ierr = PetscFree2(feRef,fvRef);CHKERRQ(ierr);
2285 
2286   ierr = DMGetGlobalVector(dmf, &fv);CHKERRQ(ierr);
2287   ierr = DMGetGlobalVector(dmc, &cv);CHKERRQ(ierr);
2288   ierr = VecGetOwnershipRange(cv, &startC, &endC);CHKERRQ(ierr);
2289   ierr = PetscSectionGetConstrainedStorageSize(cglobalSection, &m);CHKERRQ(ierr);
2290   ierr = PetscMalloc2(cTotDim,&cellCIndices,fTotDim,&cellFIndices);CHKERRQ(ierr);
2291   ierr = PetscMalloc1(m,&cindices);CHKERRQ(ierr);
2292   ierr = PetscMalloc1(m,&findices);CHKERRQ(ierr);
2293   for (d = 0; d < m; ++d) cindices[d] = findices[d] = -1;
2294   for (c = cStart; c < cEnd; ++c) {
2295     ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, c, cellCIndices, cellFIndices);CHKERRQ(ierr);
2296     for (d = 0; d < cTotDim; ++d) {
2297       if ((cellCIndices[d] < startC) || (cellCIndices[d] >= endC)) continue;
2298       if ((findices[cellCIndices[d]-startC] >= 0) && (findices[cellCIndices[d]-startC] != cellFIndices[cmap[d]])) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Coarse dof %d maps to both %d and %d", cindices[cellCIndices[d]-startC], findices[cellCIndices[d]-startC], cellFIndices[cmap[d]]);
2299       cindices[cellCIndices[d]-startC] = cellCIndices[d];
2300       findices[cellCIndices[d]-startC] = cellFIndices[cmap[d]];
2301     }
2302   }
2303   ierr = PetscFree(cmap);CHKERRQ(ierr);
2304   ierr = PetscFree2(cellCIndices,cellFIndices);CHKERRQ(ierr);
2305 
2306   ierr = ISCreateGeneral(PETSC_COMM_SELF, m, cindices, PETSC_OWN_POINTER, &cis);CHKERRQ(ierr);
2307   ierr = ISCreateGeneral(PETSC_COMM_SELF, m, findices, PETSC_OWN_POINTER, &fis);CHKERRQ(ierr);
2308   ierr = VecScatterCreate(cv, cis, fv, fis, sc);CHKERRQ(ierr);
2309   ierr = ISDestroy(&cis);CHKERRQ(ierr);
2310   ierr = ISDestroy(&fis);CHKERRQ(ierr);
2311   ierr = DMRestoreGlobalVector(dmf, &fv);CHKERRQ(ierr);
2312   ierr = DMRestoreGlobalVector(dmc, &cv);CHKERRQ(ierr);
2313   ierr = PetscLogEventEnd(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
2314   PetscFunctionReturn(0);
2315 }
2316