xref: /petsc/src/dm/impls/plex/plexfem.c (revision cff2dfc28d16e30196df121af9e4c6a29bd6c19e)
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 /*@
1010   DMPlexComputeIntegralFEM - Form the local integral F from the local input X using pointwise functions specified by the user
1011 
1012   Input Parameters:
1013 + dm - The mesh
1014 . X  - Local input vector
1015 - user - The user context
1016 
1017   Output Parameter:
1018 . integral - Local integral for each field
1019 
1020   Level: developer
1021 
1022 .seealso: DMPlexComputeResidualFEM()
1023 @*/
1024 PetscErrorCode DMPlexComputeIntegralFEM(DM dm, Vec X, PetscReal *integral, void *user)
1025 {
1026   DM_Plex           *mesh  = (DM_Plex *) dm->data;
1027   DM                 dmAux, dmGrad;
1028   Vec                localX, A, cellGeometryFVM = NULL, faceGeometryFVM = NULL, locGrad = NULL;
1029   PetscDS            prob, probAux = NULL;
1030   PetscSection       section, sectionAux;
1031   PetscFV            fvm = NULL;
1032   PetscFECellGeom   *cgeomFEM;
1033   PetscFVFaceGeom   *fgeomFVM;
1034   PetscFVCellGeom   *cgeomFVM;
1035   PetscScalar       *u, *a = NULL;
1036   const PetscScalar *constants, *lgrad;
1037   PetscReal         *lintegral;
1038   PetscInt          *uOff, *uOff_x, *aOff = NULL;
1039   PetscInt           dim, numConstants, Nf, NfAux = 0, f, numCells, cStart, cEnd, cEndInterior, c;
1040   PetscInt           totDim, totDimAux;
1041   PetscBool          useFVM = PETSC_FALSE;
1042   PetscErrorCode     ierr;
1043 
1044   PetscFunctionBegin;
1045   ierr = PetscLogEventBegin(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);
1046   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
1047   ierr = DMPlexInsertBoundaryValues(dm, PETSC_TRUE, localX, 0.0, NULL, NULL, NULL);CHKERRQ(ierr);
1048   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
1049   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
1050   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1051   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1052   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
1053   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
1054   ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr);
1055   ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr);
1056   ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr);
1057   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
1058   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1059   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
1060   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
1061   numCells = cEnd - cStart;
1062   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
1063   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
1064   if (dmAux) {
1065     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
1066     ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr);
1067     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
1068     ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
1069     ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr);
1070     ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL, NULL);CHKERRQ(ierr);
1071   }
1072   for (f = 0; f < Nf; ++f) {
1073     PetscObject  obj;
1074     PetscClassId id;
1075 
1076     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
1077     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1078     if (id == PETSCFV_CLASSID) {useFVM = PETSC_TRUE; fvm = (PetscFV) obj;}
1079   }
1080   if (useFVM) {
1081     Vec       grad;
1082     PetscInt  fStart, fEnd;
1083     PetscBool compGrad;
1084 
1085     ierr = PetscFVGetComputeGradients(fvm, &compGrad);CHKERRQ(ierr);
1086     ierr = PetscFVSetComputeGradients(fvm, PETSC_TRUE);CHKERRQ(ierr);
1087     ierr = DMPlexComputeGeometryFVM(dm, &cellGeometryFVM, &faceGeometryFVM);CHKERRQ(ierr);
1088     ierr = DMPlexComputeGradientFVM(dm, fvm, faceGeometryFVM, cellGeometryFVM, &dmGrad);CHKERRQ(ierr);
1089     ierr = PetscFVSetComputeGradients(fvm, compGrad);CHKERRQ(ierr);
1090     ierr = VecGetArrayRead(faceGeometryFVM, (const PetscScalar **) &fgeomFVM);CHKERRQ(ierr);
1091     ierr = VecGetArrayRead(cellGeometryFVM, (const PetscScalar **) &cgeomFVM);CHKERRQ(ierr);
1092     /* Reconstruct and limit cell gradients */
1093     ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
1094     ierr = DMGetGlobalVector(dmGrad, &grad);CHKERRQ(ierr);
1095     ierr = DMPlexReconstructGradients_Internal(dm, fvm, fStart, fEnd, faceGeometryFVM, cellGeometryFVM, localX, grad);CHKERRQ(ierr);
1096     /* Communicate gradient values */
1097     ierr = DMGetLocalVector(dmGrad, &locGrad);CHKERRQ(ierr);
1098     ierr = DMGlobalToLocalBegin(dmGrad, grad, INSERT_VALUES, locGrad);CHKERRQ(ierr);
1099     ierr = DMGlobalToLocalEnd(dmGrad, grad, INSERT_VALUES, locGrad);CHKERRQ(ierr);
1100     ierr = DMRestoreGlobalVector(dmGrad, &grad);CHKERRQ(ierr);
1101     /* Handle non-essential (e.g. outflow) boundary values */
1102     ierr = DMPlexInsertBoundaryValues(dm, PETSC_FALSE, localX, 0.0, faceGeometryFVM, cellGeometryFVM, locGrad);CHKERRQ(ierr);
1103     ierr = VecGetArrayRead(locGrad, &lgrad);CHKERRQ(ierr);
1104   }
1105   ierr = PetscMalloc3(Nf,&lintegral,numCells*totDim,&u,numCells,&cgeomFEM);CHKERRQ(ierr);
1106   if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);}
1107   for (f = 0; f < Nf; ++f) {lintegral[f] = 0.0;}
1108   for (c = cStart; c < cEnd; ++c) {
1109     PetscScalar *x = NULL;
1110     PetscInt     i;
1111 
1112     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, cgeomFEM[c].v0, cgeomFEM[c].J, cgeomFEM[c].invJ, &cgeomFEM[c].detJ);CHKERRQ(ierr);
1113     if (cgeomFEM[c].detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", cgeomFEM[c].detJ, c);
1114     ierr = DMPlexVecGetClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr);
1115     for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i];
1116     ierr = DMPlexVecRestoreClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr);
1117     if (dmAux) {
1118       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
1119       for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i];
1120       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
1121     }
1122   }
1123   for (f = 0; f < Nf; ++f) {
1124     PetscObject  obj;
1125     PetscClassId id;
1126     PetscInt     numChunks, numBatches, batchSize, numBlocks, blockSize, Ne, Nr, offset;
1127 
1128     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
1129     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1130     if (id == PETSCFE_CLASSID) {
1131       PetscFE         fe = (PetscFE) obj;
1132       PetscQuadrature q;
1133       PetscInt        Nq, Nb;
1134 
1135       ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
1136       ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
1137       ierr = PetscQuadratureGetData(q, NULL, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);
1138       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
1139       blockSize = Nb*Nq;
1140       batchSize = numBlocks * blockSize;
1141       ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1142       numChunks = numCells / (numBatches*batchSize);
1143       Ne        = numChunks*numBatches*batchSize;
1144       Nr        = numCells % (numBatches*batchSize);
1145       offset    = numCells - Nr;
1146       ierr = PetscFEIntegrate(fe, prob, f, Ne, cgeomFEM, u, probAux, a, lintegral);CHKERRQ(ierr);
1147       ierr = PetscFEIntegrate(fe, prob, f, Nr, &cgeomFEM[offset], &u[offset*totDim], probAux, &a[offset*totDimAux], lintegral);CHKERRQ(ierr);
1148     } else if (id == PETSCFV_CLASSID) {
1149       /* PetscFV  fv = (PetscFV) obj; */
1150       PetscInt       foff;
1151       PetscPointFunc obj_func;
1152       PetscScalar    lint;
1153 
1154       ierr = PetscDSGetObjective(prob, f, &obj_func);CHKERRQ(ierr);
1155       ierr = PetscDSGetFieldOffset(prob, f, &foff);CHKERRQ(ierr);
1156       if (obj_func) {
1157         for (c = 0; c < numCells; ++c) {
1158           PetscScalar *u_x;
1159 
1160           ierr = DMPlexPointLocalRead(dmGrad, c, lgrad, &u_x);CHKERRQ(ierr);
1161           obj_func(dim, Nf, NfAux, uOff, uOff_x, &u[totDim*c+foff], NULL, u_x, aOff, NULL, a, NULL, NULL, 0.0, cgeomFVM[c].centroid, numConstants, constants, &lint);
1162           lintegral[f] += PetscRealPart(lint)*cgeomFVM[c].volume;
1163         }
1164       }
1165     } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
1166   }
1167   if (useFVM) {
1168     ierr = VecRestoreArrayRead(locGrad, &lgrad);CHKERRQ(ierr);
1169     ierr = VecRestoreArrayRead(faceGeometryFVM, (const PetscScalar **) &fgeomFVM);CHKERRQ(ierr);
1170     ierr = VecRestoreArrayRead(cellGeometryFVM, (const PetscScalar **) &cgeomFVM);CHKERRQ(ierr);
1171     ierr = DMRestoreLocalVector(dmGrad, &locGrad);CHKERRQ(ierr);
1172     ierr = VecDestroy(&faceGeometryFVM);CHKERRQ(ierr);
1173     ierr = VecDestroy(&cellGeometryFVM);CHKERRQ(ierr);
1174     ierr = DMDestroy(&dmGrad);CHKERRQ(ierr);
1175   }
1176   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
1177   if (mesh->printFEM) {
1178     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Local integral:");CHKERRQ(ierr);
1179     for (f = 0; f < Nf; ++f) {ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), " %g", lintegral[f]);CHKERRQ(ierr);}
1180     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "\n");CHKERRQ(ierr);
1181   }
1182   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
1183   ierr = MPIU_Allreduce(lintegral, integral, Nf, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject) dm));CHKERRQ(ierr);
1184   ierr = PetscFree3(lintegral,u,cgeomFEM);CHKERRQ(ierr);
1185   ierr = PetscLogEventEnd(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);
1186   PetscFunctionReturn(0);
1187 }
1188 
1189 /*@
1190   DMPlexComputeInterpolatorNested - Form the local portion of the interpolation matrix I from the coarse DM to the uniformly refined DM.
1191 
1192   Input Parameters:
1193 + dmf  - The fine mesh
1194 . dmc  - The coarse mesh
1195 - user - The user context
1196 
1197   Output Parameter:
1198 . In  - The interpolation matrix
1199 
1200   Level: developer
1201 
1202 .seealso: DMPlexComputeInterpolatorGeneral(), DMPlexComputeJacobianFEM()
1203 @*/
1204 PetscErrorCode DMPlexComputeInterpolatorNested(DM dmc, DM dmf, Mat In, void *user)
1205 {
1206   DM_Plex          *mesh  = (DM_Plex *) dmc->data;
1207   const char       *name  = "Interpolator";
1208   PetscDS           prob;
1209   PetscFE          *feRef;
1210   PetscFV          *fvRef;
1211   PetscSection      fsection, fglobalSection;
1212   PetscSection      csection, cglobalSection;
1213   PetscScalar      *elemMat;
1214   PetscInt          dim, Nf, f, fieldI, fieldJ, offsetI, offsetJ, cStart, cEnd, cEndInterior, c;
1215   PetscInt          cTotDim, rTotDim = 0;
1216   PetscErrorCode    ierr;
1217 
1218   PetscFunctionBegin;
1219   ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1220   ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr);
1221   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
1222   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
1223   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
1224   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
1225   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
1226   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
1227   ierr = DMPlexGetHybridBounds(dmc, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
1228   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
1229   ierr = DMGetDS(dmf, &prob);CHKERRQ(ierr);
1230   ierr = PetscCalloc2(Nf,&feRef,Nf,&fvRef);CHKERRQ(ierr);
1231   for (f = 0; f < Nf; ++f) {
1232     PetscObject  obj;
1233     PetscClassId id;
1234     PetscInt     rNb = 0, Nc = 0;
1235 
1236     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
1237     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1238     if (id == PETSCFE_CLASSID) {
1239       PetscFE fe = (PetscFE) obj;
1240 
1241       ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr);
1242       ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr);
1243       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1244     } else if (id == PETSCFV_CLASSID) {
1245       PetscFV        fv = (PetscFV) obj;
1246       PetscDualSpace Q;
1247 
1248       ierr = PetscFVRefine(fv, &fvRef[f]);CHKERRQ(ierr);
1249       ierr = PetscFVGetDualSpace(fvRef[f], &Q);CHKERRQ(ierr);
1250       ierr = PetscDualSpaceGetDimension(Q, &rNb);CHKERRQ(ierr);
1251       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
1252     }
1253     rTotDim += rNb;
1254   }
1255   ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr);
1256   ierr = PetscMalloc1(rTotDim*cTotDim,&elemMat);CHKERRQ(ierr);
1257   ierr = PetscMemzero(elemMat, rTotDim*cTotDim * sizeof(PetscScalar));CHKERRQ(ierr);
1258   for (fieldI = 0, offsetI = 0; fieldI < Nf; ++fieldI) {
1259     PetscDualSpace   Qref;
1260     PetscQuadrature  f;
1261     const PetscReal *qpoints, *qweights;
1262     PetscReal       *points;
1263     PetscInt         npoints = 0, Nc, Np, fpdim, i, k, p, d;
1264 
1265     /* Compose points from all dual basis functionals */
1266     if (feRef[fieldI]) {
1267       ierr = PetscFEGetDualSpace(feRef[fieldI], &Qref);CHKERRQ(ierr);
1268       ierr = PetscFEGetNumComponents(feRef[fieldI], &Nc);CHKERRQ(ierr);
1269     } else {
1270       ierr = PetscFVGetDualSpace(fvRef[fieldI], &Qref);CHKERRQ(ierr);
1271       ierr = PetscFVGetNumComponents(fvRef[fieldI], &Nc);CHKERRQ(ierr);
1272     }
1273     ierr = PetscDualSpaceGetDimension(Qref, &fpdim);CHKERRQ(ierr);
1274     for (i = 0; i < fpdim; ++i) {
1275       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1276       ierr = PetscQuadratureGetData(f, NULL, NULL, &Np, NULL, NULL);CHKERRQ(ierr);
1277       npoints += Np;
1278     }
1279     ierr = PetscMalloc1(npoints*dim,&points);CHKERRQ(ierr);
1280     for (i = 0, k = 0; i < fpdim; ++i) {
1281       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1282       ierr = PetscQuadratureGetData(f, NULL, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr);
1283       for (p = 0; p < Np; ++p, ++k) for (d = 0; d < dim; ++d) points[k*dim+d] = qpoints[p*dim+d];
1284     }
1285 
1286     for (fieldJ = 0, offsetJ = 0; fieldJ < Nf; ++fieldJ) {
1287       PetscObject  obj;
1288       PetscClassId id;
1289       PetscReal   *B;
1290       PetscInt     NcJ = 0, cpdim = 0, j, qNc;
1291 
1292       ierr = PetscDSGetDiscretization(prob, fieldJ, &obj);CHKERRQ(ierr);
1293       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1294       if (id == PETSCFE_CLASSID) {
1295         PetscFE fe = (PetscFE) obj;
1296 
1297         /* Evaluate basis at points */
1298         ierr = PetscFEGetNumComponents(fe, &NcJ);CHKERRQ(ierr);
1299         ierr = PetscFEGetDimension(fe, &cpdim);CHKERRQ(ierr);
1300         /* For now, fields only interpolate themselves */
1301         if (fieldI == fieldJ) {
1302           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);
1303           ierr = PetscFEGetTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);
1304           for (i = 0, k = 0; i < fpdim; ++i) {
1305             ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1306             ierr = PetscQuadratureGetData(f, NULL, &qNc, &Np, NULL, &qweights);CHKERRQ(ierr);
1307             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);
1308             for (p = 0; p < Np; ++p, ++k) {
1309               for (j = 0; j < cpdim; ++j) {
1310                 /* NOTE: This is not quite right, unless fpdim == number of fine grid functional quad points */
1311                 for (c = 0; c < Nc; ++c) elemMat[(offsetI + i)*cTotDim + offsetJ + j] += B[k*cpdim*NcJ+j*Nc+c]*qweights[p*qNc+c];
1312               }
1313             }
1314           }
1315           ierr = PetscFERestoreTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr);
1316         }
1317       } else if (id == PETSCFV_CLASSID) {
1318         PetscFV        fv = (PetscFV) obj;
1319 
1320         /* Evaluate constant function at points */
1321         ierr = PetscFVGetNumComponents(fv, &NcJ);CHKERRQ(ierr);
1322         cpdim = 1;
1323         /* For now, fields only interpolate themselves */
1324         if (fieldI == fieldJ) {
1325           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);
1326           for (i = 0, k = 0; i < fpdim; ++i) {
1327             ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1328             ierr = PetscQuadratureGetData(f, NULL, &qNc, &Np, NULL, &qweights);CHKERRQ(ierr);
1329             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);
1330             for (p = 0; p < Np; ++p, ++k) {
1331               for (j = 0; j < cpdim; ++j) {
1332                 for (c = 0; c < Nc; ++c) elemMat[(offsetI + i*Nc + c)*cTotDim + offsetJ + j*NcJ + c] += 1.0*qweights[p*qNc+c];
1333               }
1334             }
1335           }
1336         }
1337       }
1338       offsetJ += cpdim*NcJ;
1339     }
1340     offsetI += fpdim*Nc;
1341     ierr = PetscFree(points);CHKERRQ(ierr);
1342   }
1343   if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(0, name, rTotDim, cTotDim, elemMat);CHKERRQ(ierr);}
1344   /* Preallocate matrix */
1345   {
1346     Mat          preallocator;
1347     PetscScalar *vals;
1348     PetscInt    *cellCIndices, *cellFIndices;
1349     PetscInt     locRows, locCols, cell;
1350 
1351     ierr = MatGetLocalSize(In, &locRows, &locCols);CHKERRQ(ierr);
1352     ierr = MatCreate(PetscObjectComm((PetscObject) In), &preallocator);CHKERRQ(ierr);
1353     ierr = MatSetType(preallocator, MATPREALLOCATOR);CHKERRQ(ierr);
1354     ierr = MatSetSizes(preallocator, locRows, locCols, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
1355     ierr = MatSetUp(preallocator);CHKERRQ(ierr);
1356     ierr = PetscCalloc3(rTotDim*cTotDim, &vals,cTotDim,&cellCIndices,rTotDim,&cellFIndices);CHKERRQ(ierr);
1357     for (cell = cStart; cell < cEnd; ++cell) {
1358       ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, cell, cellCIndices, cellFIndices);CHKERRQ(ierr);
1359       ierr = MatSetValues(preallocator, rTotDim, cellFIndices, cTotDim, cellCIndices, vals, INSERT_VALUES);CHKERRQ(ierr);
1360     }
1361     ierr = PetscFree3(vals,cellCIndices,cellFIndices);CHKERRQ(ierr);
1362     ierr = MatAssemblyBegin(preallocator, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1363     ierr = MatAssemblyEnd(preallocator, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1364     ierr = MatPreallocatorPreallocate(preallocator, PETSC_TRUE, In);CHKERRQ(ierr);
1365     ierr = MatDestroy(&preallocator);CHKERRQ(ierr);
1366   }
1367   /* Fill matrix */
1368   ierr = MatZeroEntries(In);CHKERRQ(ierr);
1369   for (c = cStart; c < cEnd; ++c) {
1370     ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES);CHKERRQ(ierr);
1371   }
1372   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);}
1373   ierr = PetscFree2(feRef,fvRef);CHKERRQ(ierr);
1374   ierr = PetscFree(elemMat);CHKERRQ(ierr);
1375   ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1376   ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1377   if (mesh->printFEM) {
1378     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
1379     ierr = MatChop(In, 1.0e-10);CHKERRQ(ierr);
1380     ierr = MatView(In, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1381   }
1382   ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1383   PetscFunctionReturn(0);
1384 }
1385 
1386 PetscErrorCode DMPlexComputeMassMatrixNested(DM dmc, DM dmf, Mat mass, void *user)
1387 {
1388   SETERRQ(PetscObjectComm((PetscObject) dmc), PETSC_ERR_SUP, "Laziness");
1389 }
1390 
1391 /*@
1392   DMPlexComputeInterpolatorGeneral - Form the local portion of the interpolation matrix I from the coarse DM to a non-nested fine DM.
1393 
1394   Input Parameters:
1395 + dmf  - The fine mesh
1396 . dmc  - The coarse mesh
1397 - user - The user context
1398 
1399   Output Parameter:
1400 . In  - The interpolation matrix
1401 
1402   Level: developer
1403 
1404 .seealso: DMPlexComputeInterpolatorNested(), DMPlexComputeJacobianFEM()
1405 @*/
1406 PetscErrorCode DMPlexComputeInterpolatorGeneral(DM dmc, DM dmf, Mat In, void *user)
1407 {
1408   DM_Plex       *mesh = (DM_Plex *) dmf->data;
1409   const char    *name = "Interpolator";
1410   PetscDS        prob;
1411   PetscSection   fsection, csection, globalFSection, globalCSection;
1412   PetscHashJK    ht;
1413   PetscLayout    rLayout;
1414   PetscInt      *dnz, *onz;
1415   PetscInt       locRows, rStart, rEnd;
1416   PetscReal     *x, *v0, *J, *invJ, detJ;
1417   PetscReal     *v0c, *Jc, *invJc, detJc;
1418   PetscScalar   *elemMat;
1419   PetscInt       dim, Nf, field, totDim, cStart, cEnd, cell, ccell;
1420   PetscErrorCode ierr;
1421 
1422   PetscFunctionBegin;
1423   ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1424   ierr = DMGetCoordinateDim(dmc, &dim);CHKERRQ(ierr);
1425   ierr = DMGetDS(dmc, &prob);CHKERRQ(ierr);
1426   ierr = PetscDSGetRefCoordArrays(prob, &x, NULL);CHKERRQ(ierr);
1427   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
1428   ierr = PetscMalloc3(dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
1429   ierr = PetscMalloc3(dim,&v0c,dim*dim,&Jc,dim*dim,&invJc);CHKERRQ(ierr);
1430   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
1431   ierr = DMGetDefaultGlobalSection(dmf, &globalFSection);CHKERRQ(ierr);
1432   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
1433   ierr = DMGetDefaultGlobalSection(dmc, &globalCSection);CHKERRQ(ierr);
1434   ierr = DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd);CHKERRQ(ierr);
1435   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
1436   ierr = PetscMalloc1(totDim, &elemMat);CHKERRQ(ierr);
1437 
1438   ierr = MatGetLocalSize(In, &locRows, NULL);CHKERRQ(ierr);
1439   ierr = PetscLayoutCreate(PetscObjectComm((PetscObject) In), &rLayout);CHKERRQ(ierr);
1440   ierr = PetscLayoutSetLocalSize(rLayout, locRows);CHKERRQ(ierr);
1441   ierr = PetscLayoutSetBlockSize(rLayout, 1);CHKERRQ(ierr);
1442   ierr = PetscLayoutSetUp(rLayout);CHKERRQ(ierr);
1443   ierr = PetscLayoutGetRange(rLayout, &rStart, &rEnd);CHKERRQ(ierr);
1444   ierr = PetscLayoutDestroy(&rLayout);CHKERRQ(ierr);
1445   ierr = PetscCalloc2(locRows,&dnz,locRows,&onz);CHKERRQ(ierr);
1446   ierr = PetscHashJKCreate(&ht);CHKERRQ(ierr);
1447   for (field = 0; field < Nf; ++field) {
1448     PetscObject      obj;
1449     PetscClassId     id;
1450     PetscDualSpace   Q = NULL;
1451     PetscQuadrature  f;
1452     const PetscReal *qpoints;
1453     PetscInt         Nc, Np, fpdim, i, d;
1454 
1455     ierr = PetscDSGetDiscretization(prob, field, &obj);CHKERRQ(ierr);
1456     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1457     if (id == PETSCFE_CLASSID) {
1458       PetscFE fe = (PetscFE) obj;
1459 
1460       ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr);
1461       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1462     } else if (id == PETSCFV_CLASSID) {
1463       PetscFV fv = (PetscFV) obj;
1464 
1465       ierr = PetscFVGetDualSpace(fv, &Q);CHKERRQ(ierr);
1466       Nc   = 1;
1467     }
1468     ierr = PetscDualSpaceGetDimension(Q, &fpdim);CHKERRQ(ierr);
1469     /* For each fine grid cell */
1470     for (cell = cStart; cell < cEnd; ++cell) {
1471       PetscInt *findices,   *cindices;
1472       PetscInt  numFIndices, numCIndices;
1473 
1474       ierr = DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr);
1475       ierr = DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
1476       if (numFIndices != fpdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of fine indices %d != %d dual basis vecs", numFIndices, fpdim);
1477       for (i = 0; i < fpdim; ++i) {
1478         Vec             pointVec;
1479         PetscScalar    *pV;
1480         PetscSF         coarseCellSF = NULL;
1481         const PetscSFNode *coarseCells;
1482         PetscInt        numCoarseCells, q, c;
1483 
1484         /* Get points from the dual basis functional quadrature */
1485         ierr = PetscDualSpaceGetFunctional(Q, i, &f);CHKERRQ(ierr);
1486         ierr = PetscQuadratureGetData(f, NULL, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr);
1487         ierr = VecCreateSeq(PETSC_COMM_SELF, Np*dim, &pointVec);CHKERRQ(ierr);
1488         ierr = VecSetBlockSize(pointVec, dim);CHKERRQ(ierr);
1489         ierr = VecGetArray(pointVec, &pV);CHKERRQ(ierr);
1490         for (q = 0; q < Np; ++q) {
1491           /* Transform point to real space */
1492           CoordinatesRefToReal(dim, dim, v0, J, &qpoints[q*dim], x);
1493           for (d = 0; d < dim; ++d) pV[q*dim+d] = x[d];
1494         }
1495         ierr = VecRestoreArray(pointVec, &pV);CHKERRQ(ierr);
1496         /* Get set of coarse cells that overlap points (would like to group points by coarse cell) */
1497         ierr = DMLocatePoints(dmc, pointVec, DM_POINTLOCATION_NEAREST, &coarseCellSF);CHKERRQ(ierr);
1498         ierr = PetscSFViewFromOptions(coarseCellSF, NULL, "-interp_sf_view");CHKERRQ(ierr);
1499         /* Update preallocation info */
1500         ierr = PetscSFGetGraph(coarseCellSF, NULL, &numCoarseCells, NULL, &coarseCells);CHKERRQ(ierr);
1501         if (numCoarseCells != Np) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Not all closure points located");
1502         {
1503           PetscHashJKKey  key;
1504           PetscHashJKIter missing, iter;
1505 
1506           key.j = findices[i];
1507           if (key.j >= 0) {
1508             /* Get indices for coarse elements */
1509             for (ccell = 0; ccell < numCoarseCells; ++ccell) {
1510               ierr = DMPlexGetClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);CHKERRQ(ierr);
1511               for (c = 0; c < numCIndices; ++c) {
1512                 key.k = cindices[c];
1513                 if (key.k < 0) continue;
1514                 ierr = PetscHashJKPut(ht, key, &missing, &iter);CHKERRQ(ierr);
1515                 if (missing) {
1516                   ierr = PetscHashJKSet(ht, iter, 1);CHKERRQ(ierr);
1517                   if ((key.k >= rStart) && (key.k < rEnd)) ++dnz[key.j-rStart];
1518                   else                                     ++onz[key.j-rStart];
1519                 }
1520               }
1521               ierr = DMPlexRestoreClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);CHKERRQ(ierr);
1522             }
1523           }
1524         }
1525         ierr = PetscSFDestroy(&coarseCellSF);CHKERRQ(ierr);
1526         ierr = VecDestroy(&pointVec);CHKERRQ(ierr);
1527       }
1528       ierr = DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr);
1529     }
1530   }
1531   ierr = PetscHashJKDestroy(&ht);CHKERRQ(ierr);
1532   ierr = MatXAIJSetPreallocation(In, 1, dnz, onz, NULL, NULL);CHKERRQ(ierr);
1533   ierr = MatSetOption(In, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
1534   ierr = PetscFree2(dnz,onz);CHKERRQ(ierr);
1535   for (field = 0; field < Nf; ++field) {
1536     PetscObject      obj;
1537     PetscClassId     id;
1538     PetscDualSpace   Q = NULL;
1539     PetscQuadrature  f;
1540     const PetscReal *qpoints, *qweights;
1541     PetscInt         Nc, qNc, Np, fpdim, i, d;
1542 
1543     ierr = PetscDSGetDiscretization(prob, field, &obj);CHKERRQ(ierr);
1544     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1545     if (id == PETSCFE_CLASSID) {
1546       PetscFE fe = (PetscFE) obj;
1547 
1548       ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr);
1549       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1550     } else if (id == PETSCFV_CLASSID) {
1551       PetscFV fv = (PetscFV) obj;
1552 
1553       ierr = PetscFVGetDualSpace(fv, &Q);CHKERRQ(ierr);
1554       Nc   = 1;
1555     }
1556     ierr = PetscDualSpaceGetDimension(Q, &fpdim);CHKERRQ(ierr);
1557     /* For each fine grid cell */
1558     for (cell = cStart; cell < cEnd; ++cell) {
1559       PetscInt *findices,   *cindices;
1560       PetscInt  numFIndices, numCIndices;
1561 
1562       ierr = DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr);
1563       ierr = DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
1564       if (numFIndices != fpdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of fine indices %d != %d dual basis vecs", numFIndices, fpdim);
1565       for (i = 0; i < fpdim; ++i) {
1566         Vec             pointVec;
1567         PetscScalar    *pV;
1568         PetscSF         coarseCellSF = NULL;
1569         const PetscSFNode *coarseCells;
1570         PetscInt        numCoarseCells, cpdim, q, c, j;
1571 
1572         /* Get points from the dual basis functional quadrature */
1573         ierr = PetscDualSpaceGetFunctional(Q, i, &f);CHKERRQ(ierr);
1574         ierr = PetscQuadratureGetData(f, NULL, &qNc, &Np, &qpoints, &qweights);CHKERRQ(ierr);
1575         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);
1576         ierr = VecCreateSeq(PETSC_COMM_SELF, Np*dim, &pointVec);CHKERRQ(ierr);
1577         ierr = VecSetBlockSize(pointVec, dim);CHKERRQ(ierr);
1578         ierr = VecGetArray(pointVec, &pV);CHKERRQ(ierr);
1579         for (q = 0; q < Np; ++q) {
1580           /* Transform point to real space */
1581           CoordinatesRefToReal(dim, dim, v0, J, &qpoints[q*dim], x);
1582           for (d = 0; d < dim; ++d) pV[q*dim+d] = x[d];
1583         }
1584         ierr = VecRestoreArray(pointVec, &pV);CHKERRQ(ierr);
1585         /* Get set of coarse cells that overlap points (would like to group points by coarse cell) */
1586         ierr = DMLocatePoints(dmc, pointVec, DM_POINTLOCATION_NEAREST, &coarseCellSF);CHKERRQ(ierr);
1587         /* Update preallocation info */
1588         ierr = PetscSFGetGraph(coarseCellSF, NULL, &numCoarseCells, NULL, &coarseCells);CHKERRQ(ierr);
1589         if (numCoarseCells != Np) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Not all closure points located");
1590         ierr = VecGetArray(pointVec, &pV);CHKERRQ(ierr);
1591         for (ccell = 0; ccell < numCoarseCells; ++ccell) {
1592           PetscReal pVReal[3];
1593 
1594           ierr = DMPlexGetClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);CHKERRQ(ierr);
1595           /* Transform points from real space to coarse reference space */
1596           ierr = DMPlexComputeCellGeometryFEM(dmc, coarseCells[ccell].index, NULL, v0c, Jc, invJc, &detJc);CHKERRQ(ierr);
1597           for (d = 0; d < dim; ++d) pVReal[d] = PetscRealPart(pV[ccell*dim+d]);
1598           CoordinatesRealToRef(dim, dim, v0c, invJc, pVReal, x);
1599 
1600           if (id == PETSCFE_CLASSID) {
1601             PetscFE    fe = (PetscFE) obj;
1602             PetscReal *B;
1603 
1604             /* Evaluate coarse basis on contained point */
1605             ierr = PetscFEGetDimension(fe, &cpdim);CHKERRQ(ierr);
1606             ierr = PetscFEGetTabulation(fe, 1, x, &B, NULL, NULL);CHKERRQ(ierr);
1607             ierr = PetscMemzero(elemMat, cpdim * sizeof(PetscScalar));CHKERRQ(ierr);
1608             /* Get elemMat entries by multiplying by weight */
1609             for (j = 0; j < cpdim; ++j) {
1610               for (c = 0; c < Nc; ++c) elemMat[j] += B[j*Nc + c]*qweights[ccell*qNc + c];
1611             }
1612             ierr = PetscFERestoreTabulation(fe, 1, x, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr);
1613           } else {
1614             cpdim = 1;
1615             for (j = 0; j < cpdim; ++j) {
1616               for (c = 0; c < Nc; ++c) elemMat[j] += 1.0*qweights[ccell*qNc + c];
1617             }
1618           }
1619           /* Update interpolator */
1620           if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat);CHKERRQ(ierr);}
1621           if (numCIndices != cpdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of element matrix columns %D != %D", numCIndices, cpdim);
1622           ierr = MatSetValues(In, 1, &findices[i], numCIndices, cindices, elemMat, INSERT_VALUES);CHKERRQ(ierr);
1623           ierr = DMPlexRestoreClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);CHKERRQ(ierr);
1624         }
1625         ierr = VecRestoreArray(pointVec, &pV);CHKERRQ(ierr);
1626         ierr = PetscSFDestroy(&coarseCellSF);CHKERRQ(ierr);
1627         ierr = VecDestroy(&pointVec);CHKERRQ(ierr);
1628       }
1629       ierr = DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr);
1630     }
1631   }
1632   ierr = PetscFree3(v0,J,invJ);CHKERRQ(ierr);
1633   ierr = PetscFree3(v0c,Jc,invJc);CHKERRQ(ierr);
1634   ierr = PetscFree(elemMat);CHKERRQ(ierr);
1635   ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1636   ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1637   ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1638   PetscFunctionReturn(0);
1639 }
1640 
1641 /*@
1642   DMPlexComputeMassMatrixGeneral - Form the local portion of the mass matrix M from the coarse DM to a non-nested fine DM.
1643 
1644   Input Parameters:
1645 + dmf  - The fine mesh
1646 . dmc  - The coarse mesh
1647 - user - The user context
1648 
1649   Output Parameter:
1650 . mass  - The mass matrix
1651 
1652   Level: developer
1653 
1654 .seealso: DMPlexComputeMassMatrixNested(), DMPlexComputeInterpolatorNested(), DMPlexComputeInterpolatorGeneral(), DMPlexComputeJacobianFEM()
1655 @*/
1656 PetscErrorCode DMPlexComputeMassMatrixGeneral(DM dmc, DM dmf, Mat mass, void *user)
1657 {
1658   DM_Plex       *mesh = (DM_Plex *) dmf->data;
1659   const char    *name = "Mass Matrix";
1660   PetscDS        prob;
1661   PetscSection   fsection, csection, globalFSection, globalCSection;
1662   PetscHashJK    ht;
1663   PetscLayout    rLayout;
1664   PetscInt      *dnz, *onz;
1665   PetscInt       locRows, rStart, rEnd;
1666   PetscReal     *x, *v0, *J, *invJ, detJ;
1667   PetscReal     *v0c, *Jc, *invJc, detJc;
1668   PetscScalar   *elemMat;
1669   PetscInt       dim, Nf, field, totDim, cStart, cEnd, cell, ccell;
1670   PetscErrorCode ierr;
1671 
1672   PetscFunctionBegin;
1673   ierr = DMGetCoordinateDim(dmc, &dim);CHKERRQ(ierr);
1674   ierr = DMGetDS(dmc, &prob);CHKERRQ(ierr);
1675   ierr = PetscDSGetRefCoordArrays(prob, &x, NULL);CHKERRQ(ierr);
1676   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
1677   ierr = PetscMalloc3(dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
1678   ierr = PetscMalloc3(dim,&v0c,dim*dim,&Jc,dim*dim,&invJc);CHKERRQ(ierr);
1679   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
1680   ierr = DMGetDefaultGlobalSection(dmf, &globalFSection);CHKERRQ(ierr);
1681   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
1682   ierr = DMGetDefaultGlobalSection(dmc, &globalCSection);CHKERRQ(ierr);
1683   ierr = DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd);CHKERRQ(ierr);
1684   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
1685   ierr = PetscMalloc1(totDim, &elemMat);CHKERRQ(ierr);
1686 
1687   ierr = MatGetLocalSize(mass, &locRows, NULL);CHKERRQ(ierr);
1688   ierr = PetscLayoutCreate(PetscObjectComm((PetscObject) mass), &rLayout);CHKERRQ(ierr);
1689   ierr = PetscLayoutSetLocalSize(rLayout, locRows);CHKERRQ(ierr);
1690   ierr = PetscLayoutSetBlockSize(rLayout, 1);CHKERRQ(ierr);
1691   ierr = PetscLayoutSetUp(rLayout);CHKERRQ(ierr);
1692   ierr = PetscLayoutGetRange(rLayout, &rStart, &rEnd);CHKERRQ(ierr);
1693   ierr = PetscLayoutDestroy(&rLayout);CHKERRQ(ierr);
1694   ierr = PetscCalloc2(locRows,&dnz,locRows,&onz);CHKERRQ(ierr);
1695   ierr = PetscHashJKCreate(&ht);CHKERRQ(ierr);
1696   for (field = 0; field < Nf; ++field) {
1697     PetscObject      obj;
1698     PetscClassId     id;
1699     PetscQuadrature  quad;
1700     const PetscReal *qpoints;
1701     PetscInt         Nq, Nc, i, d;
1702 
1703     ierr = PetscDSGetDiscretization(prob, field, &obj);CHKERRQ(ierr);
1704     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1705     if (id == PETSCFE_CLASSID) {ierr = PetscFEGetQuadrature((PetscFE) obj, &quad);CHKERRQ(ierr);}
1706     else                       {ierr = PetscFVGetQuadrature((PetscFV) obj, &quad);CHKERRQ(ierr);}
1707     ierr = PetscQuadratureGetData(quad, NULL, &Nc, &Nq, &qpoints, NULL);CHKERRQ(ierr);
1708     /* For each fine grid cell */
1709     for (cell = cStart; cell < cEnd; ++cell) {
1710       Vec                pointVec;
1711       PetscScalar       *pV;
1712       PetscSF            coarseCellSF = NULL;
1713       const PetscSFNode *coarseCells;
1714       PetscInt           numCoarseCells, q, c;
1715       PetscInt          *findices,   *cindices;
1716       PetscInt           numFIndices, numCIndices;
1717 
1718       ierr = DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr);
1719       ierr = DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
1720       /* Get points from the quadrature */
1721       ierr = VecCreateSeq(PETSC_COMM_SELF, Nq*dim, &pointVec);CHKERRQ(ierr);
1722       ierr = VecSetBlockSize(pointVec, dim);CHKERRQ(ierr);
1723       ierr = VecGetArray(pointVec, &pV);CHKERRQ(ierr);
1724       for (q = 0; q < Nq; ++q) {
1725         /* Transform point to real space */
1726         CoordinatesRefToReal(dim, dim, v0, J, &qpoints[q*dim], x);
1727         for (d = 0; d < dim; ++d) pV[q*dim+d] = x[d];
1728       }
1729       ierr = VecRestoreArray(pointVec, &pV);CHKERRQ(ierr);
1730       /* Get set of coarse cells that overlap points (would like to group points by coarse cell) */
1731       ierr = DMLocatePoints(dmc, pointVec, DM_POINTLOCATION_NEAREST, &coarseCellSF);CHKERRQ(ierr);
1732       ierr = PetscSFViewFromOptions(coarseCellSF, NULL, "-interp_sf_view");CHKERRQ(ierr);
1733       /* Update preallocation info */
1734       ierr = PetscSFGetGraph(coarseCellSF, NULL, &numCoarseCells, NULL, &coarseCells);CHKERRQ(ierr);
1735       if (numCoarseCells != Nq) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Not all closure points located");
1736       {
1737         PetscHashJKKey  key;
1738         PetscHashJKIter missing, iter;
1739 
1740         for (i = 0; i < numFIndices; ++i) {
1741           key.j = findices[i];
1742           if (key.j >= 0) {
1743             /* Get indices for coarse elements */
1744             for (ccell = 0; ccell < numCoarseCells; ++ccell) {
1745               ierr = DMPlexGetClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);CHKERRQ(ierr);
1746               for (c = 0; c < numCIndices; ++c) {
1747                 key.k = cindices[c];
1748                 if (key.k < 0) continue;
1749                 ierr = PetscHashJKPut(ht, key, &missing, &iter);CHKERRQ(ierr);
1750                 if (missing) {
1751                   ierr = PetscHashJKSet(ht, iter, 1);CHKERRQ(ierr);
1752                   if ((key.k >= rStart) && (key.k < rEnd)) ++dnz[key.j-rStart];
1753                   else                                     ++onz[key.j-rStart];
1754                 }
1755               }
1756               ierr = DMPlexRestoreClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);CHKERRQ(ierr);
1757             }
1758           }
1759         }
1760       }
1761       ierr = PetscSFDestroy(&coarseCellSF);CHKERRQ(ierr);
1762       ierr = VecDestroy(&pointVec);CHKERRQ(ierr);
1763       ierr = DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr);
1764     }
1765   }
1766   ierr = PetscHashJKDestroy(&ht);CHKERRQ(ierr);
1767   ierr = MatXAIJSetPreallocation(mass, 1, dnz, onz, NULL, NULL);CHKERRQ(ierr);
1768   ierr = MatSetOption(mass, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
1769   ierr = PetscFree2(dnz,onz);CHKERRQ(ierr);
1770   for (field = 0; field < Nf; ++field) {
1771     PetscObject      obj;
1772     PetscClassId     id;
1773     PetscQuadrature  quad;
1774     PetscReal       *Bfine;
1775     const PetscReal *qpoints, *qweights;
1776     PetscInt         Nq, Nc, i, d;
1777 
1778     ierr = PetscDSGetDiscretization(prob, field, &obj);CHKERRQ(ierr);
1779     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1780     if (id == PETSCFE_CLASSID) {ierr = PetscFEGetQuadrature((PetscFE) obj, &quad);CHKERRQ(ierr);ierr = PetscFEGetDefaultTabulation((PetscFE) obj, &Bfine, NULL, NULL);CHKERRQ(ierr);}
1781     else                       {ierr = PetscFVGetQuadrature((PetscFV) obj, &quad);CHKERRQ(ierr);}
1782     ierr = PetscQuadratureGetData(quad, NULL, &Nc, &Nq, &qpoints, &qweights);CHKERRQ(ierr);
1783     /* For each fine grid cell */
1784     for (cell = cStart; cell < cEnd; ++cell) {
1785       Vec                pointVec;
1786       PetscScalar       *pV;
1787       PetscSF            coarseCellSF = NULL;
1788       const PetscSFNode *coarseCells;
1789       PetscInt           numCoarseCells, cpdim, q, c, j;
1790       PetscInt          *findices,   *cindices;
1791       PetscInt           numFIndices, numCIndices;
1792 
1793       ierr = DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr);
1794       ierr = DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
1795       /* Get points from the quadrature */
1796       ierr = VecCreateSeq(PETSC_COMM_SELF, Nq*dim, &pointVec);CHKERRQ(ierr);
1797       ierr = VecSetBlockSize(pointVec, dim);CHKERRQ(ierr);
1798       ierr = VecGetArray(pointVec, &pV);CHKERRQ(ierr);
1799       for (q = 0; q < Nq; ++q) {
1800         /* Transform point to real space */
1801         CoordinatesRefToReal(dim, dim, v0, J, &qpoints[q*dim], x);
1802         for (d = 0; d < dim; ++d) pV[q*dim+d] = x[d];
1803       }
1804       ierr = VecRestoreArray(pointVec, &pV);CHKERRQ(ierr);
1805       /* Get set of coarse cells that overlap points (would like to group points by coarse cell) */
1806       ierr = DMLocatePoints(dmc, pointVec, DM_POINTLOCATION_NEAREST, &coarseCellSF);CHKERRQ(ierr);
1807       /* Update matrix */
1808       ierr = PetscSFGetGraph(coarseCellSF, NULL, &numCoarseCells, NULL, &coarseCells);CHKERRQ(ierr);
1809       if (numCoarseCells != Nq) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Not all closure points located");
1810       ierr = VecGetArray(pointVec, &pV);CHKERRQ(ierr);
1811       for (ccell = 0; ccell < numCoarseCells; ++ccell) {
1812         PetscReal pVReal[3];
1813 
1814         ierr = DMPlexGetClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);CHKERRQ(ierr);
1815         /* Transform points from real space to coarse reference space */
1816         ierr = DMPlexComputeCellGeometryFEM(dmc, coarseCells[ccell].index, NULL, v0c, Jc, invJc, &detJc);CHKERRQ(ierr);
1817         for (d = 0; d < dim; ++d) pVReal[d] = PetscRealPart(pV[ccell*dim+d]);
1818         CoordinatesRealToRef(dim, dim, v0c, invJc, pVReal, x);
1819 
1820         if (id == PETSCFE_CLASSID) {
1821           PetscFE    fe = (PetscFE) obj;
1822           PetscReal *B;
1823 
1824           /* Evaluate coarse basis on contained point */
1825           ierr = PetscFEGetDimension(fe, &cpdim);CHKERRQ(ierr);
1826           ierr = PetscFEGetTabulation(fe, 1, x, &B, NULL, NULL);CHKERRQ(ierr);
1827           /* Get elemMat entries by multiplying by weight */
1828           for (i = 0; i < numFIndices; ++i) {
1829             ierr = PetscMemzero(elemMat, cpdim * sizeof(PetscScalar));CHKERRQ(ierr);
1830             for (j = 0; j < cpdim; ++j) {
1831               for (c = 0; c < Nc; ++c) elemMat[j] += B[j*Nc + c]*Bfine[(ccell*numFIndices + i)*Nc + c]*qweights[ccell*Nc + c]*detJ;
1832             }
1833             /* Update interpolator */
1834             if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat);CHKERRQ(ierr);}
1835             if (numCIndices != cpdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of element matrix columns %D != %D", numCIndices, cpdim);
1836             ierr = MatSetValues(mass, 1, &findices[i], numCIndices, cindices, elemMat, ADD_VALUES);CHKERRQ(ierr);
1837           }
1838           ierr = PetscFERestoreTabulation(fe, 1, x, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr);
1839         } else {
1840           cpdim = 1;
1841           for (i = 0; i < numFIndices; ++i) {
1842             ierr = PetscMemzero(elemMat, cpdim * sizeof(PetscScalar));CHKERRQ(ierr);
1843             for (j = 0; j < cpdim; ++j) {
1844               for (c = 0; c < Nc; ++c) elemMat[j] += 1.0*1.0*qweights[ccell*Nc + c]*detJ;
1845             }
1846             /* Update interpolator */
1847             if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat);CHKERRQ(ierr);}
1848             ierr = PetscPrintf(PETSC_COMM_SELF, "Nq: %d %d Nf: %d %d Nc: %d %d\n", ccell, Nq, i, numFIndices, j, numCIndices);CHKERRQ(ierr);
1849             if (numCIndices != cpdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of element matrix columns %D != %D", numCIndices, cpdim);
1850             ierr = MatSetValues(mass, 1, &findices[i], numCIndices, cindices, elemMat, ADD_VALUES);CHKERRQ(ierr);
1851           }
1852         }
1853         ierr = DMPlexRestoreClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);CHKERRQ(ierr);
1854       }
1855       ierr = VecRestoreArray(pointVec, &pV);CHKERRQ(ierr);
1856       ierr = PetscSFDestroy(&coarseCellSF);CHKERRQ(ierr);
1857       ierr = VecDestroy(&pointVec);CHKERRQ(ierr);
1858       ierr = DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr);
1859     }
1860   }
1861   ierr = PetscFree3(v0,J,invJ);CHKERRQ(ierr);
1862   ierr = PetscFree3(v0c,Jc,invJc);CHKERRQ(ierr);
1863   ierr = PetscFree(elemMat);CHKERRQ(ierr);
1864   ierr = MatAssemblyBegin(mass, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1865   ierr = MatAssemblyEnd(mass, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1866   PetscFunctionReturn(0);
1867 }
1868 
1869 /*@
1870   DMPlexComputeInjectorFEM - Compute a mapping from coarse unknowns to fine unknowns
1871 
1872   Input Parameters:
1873 + dmc  - The coarse mesh
1874 - dmf  - The fine mesh
1875 - user - The user context
1876 
1877   Output Parameter:
1878 . sc   - The mapping
1879 
1880   Level: developer
1881 
1882 .seealso: DMPlexComputeInterpolatorNested(), DMPlexComputeJacobianFEM()
1883 @*/
1884 PetscErrorCode DMPlexComputeInjectorFEM(DM dmc, DM dmf, VecScatter *sc, void *user)
1885 {
1886   PetscDS        prob;
1887   PetscFE       *feRef;
1888   PetscFV       *fvRef;
1889   Vec            fv, cv;
1890   IS             fis, cis;
1891   PetscSection   fsection, fglobalSection, csection, cglobalSection;
1892   PetscInt      *cmap, *cellCIndices, *cellFIndices, *cindices, *findices;
1893   PetscInt       cTotDim, fTotDim = 0, Nf, f, field, cStart, cEnd, cEndInterior, c, dim, d, startC, endC, offsetC, offsetF, m;
1894   PetscErrorCode ierr;
1895 
1896   PetscFunctionBegin;
1897   ierr = PetscLogEventBegin(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1898   ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr);
1899   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
1900   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
1901   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
1902   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
1903   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
1904   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
1905   ierr = DMPlexGetHybridBounds(dmc, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
1906   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
1907   ierr = DMGetDS(dmc, &prob);CHKERRQ(ierr);
1908   ierr = PetscCalloc2(Nf,&feRef,Nf,&fvRef);CHKERRQ(ierr);
1909   for (f = 0; f < Nf; ++f) {
1910     PetscObject  obj;
1911     PetscClassId id;
1912     PetscInt     fNb = 0, Nc = 0;
1913 
1914     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
1915     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1916     if (id == PETSCFE_CLASSID) {
1917       PetscFE fe = (PetscFE) obj;
1918 
1919       ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr);
1920       ierr = PetscFEGetDimension(feRef[f], &fNb);CHKERRQ(ierr);
1921       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1922     } else if (id == PETSCFV_CLASSID) {
1923       PetscFV        fv = (PetscFV) obj;
1924       PetscDualSpace Q;
1925 
1926       ierr = PetscFVRefine(fv, &fvRef[f]);CHKERRQ(ierr);
1927       ierr = PetscFVGetDualSpace(fvRef[f], &Q);CHKERRQ(ierr);
1928       ierr = PetscDualSpaceGetDimension(Q, &fNb);CHKERRQ(ierr);
1929       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
1930     }
1931     fTotDim += fNb*Nc;
1932   }
1933   ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr);
1934   ierr = PetscMalloc1(cTotDim,&cmap);CHKERRQ(ierr);
1935   for (field = 0, offsetC = 0, offsetF = 0; field < Nf; ++field) {
1936     PetscFE        feC;
1937     PetscFV        fvC;
1938     PetscDualSpace QF, QC;
1939     PetscInt       NcF, NcC, fpdim, cpdim;
1940 
1941     if (feRef[field]) {
1942       ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &feC);CHKERRQ(ierr);
1943       ierr = PetscFEGetNumComponents(feC, &NcC);CHKERRQ(ierr);
1944       ierr = PetscFEGetNumComponents(feRef[field], &NcF);CHKERRQ(ierr);
1945       ierr = PetscFEGetDualSpace(feRef[field], &QF);CHKERRQ(ierr);
1946       ierr = PetscDualSpaceGetDimension(QF, &fpdim);CHKERRQ(ierr);
1947       ierr = PetscFEGetDualSpace(feC, &QC);CHKERRQ(ierr);
1948       ierr = PetscDualSpaceGetDimension(QC, &cpdim);CHKERRQ(ierr);
1949     } else {
1950       ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fvC);CHKERRQ(ierr);
1951       ierr = PetscFVGetNumComponents(fvC, &NcC);CHKERRQ(ierr);
1952       ierr = PetscFVGetNumComponents(fvRef[field], &NcF);CHKERRQ(ierr);
1953       ierr = PetscFVGetDualSpace(fvRef[field], &QF);CHKERRQ(ierr);
1954       ierr = PetscDualSpaceGetDimension(QF, &fpdim);CHKERRQ(ierr);
1955       ierr = PetscFVGetDualSpace(fvC, &QC);CHKERRQ(ierr);
1956       ierr = PetscDualSpaceGetDimension(QC, &cpdim);CHKERRQ(ierr);
1957     }
1958     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);
1959     for (c = 0; c < cpdim; ++c) {
1960       PetscQuadrature  cfunc;
1961       const PetscReal *cqpoints;
1962       PetscInt         NpC;
1963       PetscBool        found = PETSC_FALSE;
1964 
1965       ierr = PetscDualSpaceGetFunctional(QC, c, &cfunc);CHKERRQ(ierr);
1966       ierr = PetscQuadratureGetData(cfunc, NULL, NULL, &NpC, &cqpoints, NULL);CHKERRQ(ierr);
1967       if (NpC != 1 && feRef[field]) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not know how to do injection for moments");
1968       for (f = 0; f < fpdim; ++f) {
1969         PetscQuadrature  ffunc;
1970         const PetscReal *fqpoints;
1971         PetscReal        sum = 0.0;
1972         PetscInt         NpF, comp;
1973 
1974         ierr = PetscDualSpaceGetFunctional(QF, f, &ffunc);CHKERRQ(ierr);
1975         ierr = PetscQuadratureGetData(ffunc, NULL, NULL, &NpF, &fqpoints, NULL);CHKERRQ(ierr);
1976         if (NpC != NpF) continue;
1977         for (d = 0; d < dim; ++d) sum += PetscAbsReal(cqpoints[d] - fqpoints[d]);
1978         if (sum > 1.0e-9) continue;
1979         for (comp = 0; comp < NcC; ++comp) {
1980           cmap[(offsetC+c)*NcC+comp] = (offsetF+f)*NcF+comp;
1981         }
1982         found = PETSC_TRUE;
1983         break;
1984       }
1985       if (!found) {
1986         /* TODO We really want the average here, but some asshole put VecScatter in the interface */
1987         if (fvRef[field]) {
1988           PetscInt comp;
1989           for (comp = 0; comp < NcC; ++comp) {
1990             cmap[(offsetC+c)*NcC+comp] = (offsetF+0)*NcF+comp;
1991           }
1992         } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not locate matching functional for injection");
1993       }
1994     }
1995     offsetC += cpdim*NcC;
1996     offsetF += fpdim*NcF;
1997   }
1998   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);ierr = PetscFVDestroy(&fvRef[f]);CHKERRQ(ierr);}
1999   ierr = PetscFree2(feRef,fvRef);CHKERRQ(ierr);
2000 
2001   ierr = DMGetGlobalVector(dmf, &fv);CHKERRQ(ierr);
2002   ierr = DMGetGlobalVector(dmc, &cv);CHKERRQ(ierr);
2003   ierr = VecGetOwnershipRange(cv, &startC, &endC);CHKERRQ(ierr);
2004   ierr = PetscSectionGetConstrainedStorageSize(cglobalSection, &m);CHKERRQ(ierr);
2005   ierr = PetscMalloc2(cTotDim,&cellCIndices,fTotDim,&cellFIndices);CHKERRQ(ierr);
2006   ierr = PetscMalloc1(m,&cindices);CHKERRQ(ierr);
2007   ierr = PetscMalloc1(m,&findices);CHKERRQ(ierr);
2008   for (d = 0; d < m; ++d) cindices[d] = findices[d] = -1;
2009   for (c = cStart; c < cEnd; ++c) {
2010     ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, c, cellCIndices, cellFIndices);CHKERRQ(ierr);
2011     for (d = 0; d < cTotDim; ++d) {
2012       if ((cellCIndices[d] < startC) || (cellCIndices[d] >= endC)) continue;
2013       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]]);
2014       cindices[cellCIndices[d]-startC] = cellCIndices[d];
2015       findices[cellCIndices[d]-startC] = cellFIndices[cmap[d]];
2016     }
2017   }
2018   ierr = PetscFree(cmap);CHKERRQ(ierr);
2019   ierr = PetscFree2(cellCIndices,cellFIndices);CHKERRQ(ierr);
2020 
2021   ierr = ISCreateGeneral(PETSC_COMM_SELF, m, cindices, PETSC_OWN_POINTER, &cis);CHKERRQ(ierr);
2022   ierr = ISCreateGeneral(PETSC_COMM_SELF, m, findices, PETSC_OWN_POINTER, &fis);CHKERRQ(ierr);
2023   ierr = VecScatterCreate(cv, cis, fv, fis, sc);CHKERRQ(ierr);
2024   ierr = ISDestroy(&cis);CHKERRQ(ierr);
2025   ierr = ISDestroy(&fis);CHKERRQ(ierr);
2026   ierr = DMRestoreGlobalVector(dmf, &fv);CHKERRQ(ierr);
2027   ierr = DMRestoreGlobalVector(dmc, &cv);CHKERRQ(ierr);
2028   ierr = PetscLogEventEnd(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
2029   PetscFunctionReturn(0);
2030 }
2031