xref: /petsc/src/dm/impls/plex/plexfem.c (revision 9c334d8fdb557fc53fd345d68cbb3545b09ccab8)
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 #undef __FUNCT__
8 #define __FUNCT__ "DMPlexGetScale"
9 PetscErrorCode DMPlexGetScale(DM dm, PetscUnit unit, PetscReal *scale)
10 {
11   DM_Plex *mesh = (DM_Plex*) dm->data;
12 
13   PetscFunctionBegin;
14   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
15   PetscValidPointer(scale, 3);
16   *scale = mesh->scale[unit];
17   PetscFunctionReturn(0);
18 }
19 
20 #undef __FUNCT__
21 #define __FUNCT__ "DMPlexSetScale"
22 PetscErrorCode DMPlexSetScale(DM dm, PetscUnit unit, PetscReal scale)
23 {
24   DM_Plex *mesh = (DM_Plex*) dm->data;
25 
26   PetscFunctionBegin;
27   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
28   mesh->scale[unit] = scale;
29   PetscFunctionReturn(0);
30 }
31 
32 PETSC_STATIC_INLINE PetscInt epsilon(PetscInt i, PetscInt j, PetscInt k)
33 {
34   switch (i) {
35   case 0:
36     switch (j) {
37     case 0: return 0;
38     case 1:
39       switch (k) {
40       case 0: return 0;
41       case 1: return 0;
42       case 2: return 1;
43       }
44     case 2:
45       switch (k) {
46       case 0: return 0;
47       case 1: return -1;
48       case 2: return 0;
49       }
50     }
51   case 1:
52     switch (j) {
53     case 0:
54       switch (k) {
55       case 0: return 0;
56       case 1: return 0;
57       case 2: return -1;
58       }
59     case 1: return 0;
60     case 2:
61       switch (k) {
62       case 0: return 1;
63       case 1: return 0;
64       case 2: return 0;
65       }
66     }
67   case 2:
68     switch (j) {
69     case 0:
70       switch (k) {
71       case 0: return 0;
72       case 1: return 1;
73       case 2: return 0;
74       }
75     case 1:
76       switch (k) {
77       case 0: return -1;
78       case 1: return 0;
79       case 2: return 0;
80       }
81     case 2: return 0;
82     }
83   }
84   return 0;
85 }
86 
87 #undef __FUNCT__
88 #define __FUNCT__ "DMPlexProjectRigidBody"
89 static PetscErrorCode DMPlexProjectRigidBody(PetscInt dim, PetscReal t, const PetscReal X[], PetscInt Nf, PetscScalar *mode, void *ctx)
90 {
91   PetscInt *ctxInt  = (PetscInt *) ctx;
92   PetscInt  dim2    = ctxInt[0];
93   PetscInt  d       = ctxInt[1];
94   PetscInt  i, j, k = dim > 2 ? d - dim : d;
95 
96   PetscFunctionBegin;
97   if (dim != dim2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Input dimension %d does not match context dimension %d", dim, dim2);
98   for (i = 0; i < dim; i++) mode[i] = 0.;
99   if (d < dim) {
100     mode[d] = 1.;
101   } else {
102     for (i = 0; i < dim; i++) {
103       for (j = 0; j < dim; j++) {
104         mode[j] += epsilon(i, j, k)*X[i];
105       }
106     }
107   }
108   PetscFunctionReturn(0);
109 }
110 
111 #undef __FUNCT__
112 #define __FUNCT__ "DMPlexCreateRigidBody"
113 /*@C
114   DMPlexCreateRigidBody - for the default global section, create rigid body modes from coordinates
115 
116   Collective on DM
117 
118   Input Arguments:
119 . dm - the DM
120 
121   Output Argument:
122 . sp - the null space
123 
124   Note: This is necessary to take account of Dirichlet conditions on the displacements
125 
126   Level: advanced
127 
128 .seealso: MatNullSpaceCreate()
129 @*/
130 PetscErrorCode DMPlexCreateRigidBody(DM dm, MatNullSpace *sp)
131 {
132   MPI_Comm       comm;
133   Vec            mode[6];
134   PetscSection   section, globalSection;
135   PetscInt       dim, dimEmbed, n, m, d, i, j;
136   PetscErrorCode ierr;
137 
138   PetscFunctionBegin;
139   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
140   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
141   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
142   if (dim == 1) {
143     ierr = MatNullSpaceCreate(comm, PETSC_TRUE, 0, NULL, sp);CHKERRQ(ierr);
144     PetscFunctionReturn(0);
145   }
146   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
147   ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);
148   ierr = PetscSectionGetConstrainedStorageSize(globalSection, &n);CHKERRQ(ierr);
149   m    = (dim*(dim+1))/2;
150   ierr = VecCreate(comm, &mode[0]);CHKERRQ(ierr);
151   ierr = VecSetSizes(mode[0], n, PETSC_DETERMINE);CHKERRQ(ierr);
152   ierr = VecSetUp(mode[0]);CHKERRQ(ierr);
153   for (i = 1; i < m; ++i) {ierr = VecDuplicate(mode[0], &mode[i]);CHKERRQ(ierr);}
154   for (d = 0; d < m; d++) {
155     PetscInt ctx[2];
156     void    *voidctx = (void *) (&ctx[0]);
157     PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal *, PetscInt, PetscScalar *, void *) = DMPlexProjectRigidBody;
158 
159     ctx[0] = dimEmbed;
160     ctx[1] = d;
161     ierr = DMProjectFunction(dm, 0.0, &func, &voidctx, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
162   }
163   for (i = 0; i < dim; ++i) {ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);}
164   /* Orthonormalize system */
165   for (i = dim; i < m; ++i) {
166     PetscScalar dots[6];
167 
168     ierr = VecMDot(mode[i], i, mode, dots);CHKERRQ(ierr);
169     for (j = 0; j < i; ++j) dots[j] *= -1.0;
170     ierr = VecMAXPY(mode[i], i, dots, mode);CHKERRQ(ierr);
171     ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);
172   }
173   ierr = MatNullSpaceCreate(comm, PETSC_FALSE, m, mode, sp);CHKERRQ(ierr);
174   for (i = 0; i< m; ++i) {ierr = VecDestroy(&mode[i]);CHKERRQ(ierr);}
175   PetscFunctionReturn(0);
176 }
177 
178 #undef __FUNCT__
179 #define __FUNCT__ "DMPlexSetMaxProjectionHeight"
180 /*@
181   DMPlexSetMaxProjectionHeight - In DMPlexProjectXXXLocal() functions, the projected values of a basis function's dofs
182   are computed by associating the basis function with one of the mesh points in its transitively-closed support, and
183   evaluating the dual space basis of that point.  A basis function is associated with the point in its
184   transitively-closed support whose mesh height is highest (w.r.t. DAG height), but not greater than the maximum
185   projection height, which is set with this function.  By default, the maximum projection height is zero, which means
186   that only mesh cells are used to project basis functions.  A height of one, for example, evaluates a cell-interior
187   basis functions using its cells dual space basis, but all other basis functions with the dual space basis of a face.
188 
189   Input Parameters:
190 + dm - the DMPlex object
191 - height - the maximum projection height >= 0
192 
193   Level: advanced
194 
195 .seealso: DMPlexGetMaxProjectionHeight(), DMProjectFunctionLocal(), DMProjectFunctionLabelLocal()
196 @*/
197 PetscErrorCode DMPlexSetMaxProjectionHeight(DM dm, PetscInt height)
198 {
199   DM_Plex *plex = (DM_Plex *) dm->data;
200 
201   PetscFunctionBegin;
202   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
203   plex->maxProjectionHeight = height;
204   PetscFunctionReturn(0);
205 }
206 
207 #undef __FUNCT__
208 #define __FUNCT__ "DMPlexGetMaxProjectionHeight"
209 /*@
210   DMPlexGetMaxProjectionHeight - Get the maximum height (w.r.t. DAG) of mesh points used to evaluate dual bases in
211   DMPlexProjectXXXLocal() functions.
212 
213   Input Parameters:
214 . dm - the DMPlex object
215 
216   Output Parameters:
217 . height - the maximum projection height
218 
219   Level: intermediate
220 
221 .seealso: DMPlexSetMaxProjectionHeight(), DMProjectFunctionLocal(), DMProjectFunctionLabelLocal()
222 @*/
223 PetscErrorCode DMPlexGetMaxProjectionHeight(DM dm, PetscInt *height)
224 {
225   DM_Plex *plex = (DM_Plex *) dm->data;
226 
227   PetscFunctionBegin;
228   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
229   *height = plex->maxProjectionHeight;
230   PetscFunctionReturn(0);
231 }
232 
233 #undef __FUNCT__
234 #define __FUNCT__ "DMProjectFunctionLabelLocal_Plex"
235 PetscErrorCode DMProjectFunctionLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
236 {
237   PetscDualSpace *sp, *cellsp;
238   PetscInt       *numComp;
239   PetscSection    section;
240   PetscScalar    *values;
241   PetscBool      *fieldActive;
242   PetscInt        numFields, dim, dimEmbed, spDim, totDim = 0, numValues, pStart, pEnd, cStart, cEnd, cEndInterior, f, d, v, i, comp, maxHeight, h;
243   PetscErrorCode  ierr;
244 
245   PetscFunctionBegin;
246   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
247   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
248   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
249   if (cEnd <= cStart) PetscFunctionReturn(0);
250   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
251   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
252   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
253   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
254   ierr = PetscMalloc2(numFields,&sp,numFields,&numComp);CHKERRQ(ierr);
255   ierr = DMPlexGetMaxProjectionHeight(dm,&maxHeight);CHKERRQ(ierr);
256   if (maxHeight < 0 || maxHeight > dim) {SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"maximum projection height %d not in [0, %d)\n", maxHeight,dim);}
257   if (maxHeight > 0) {ierr = PetscMalloc1(numFields,&cellsp);CHKERRQ(ierr);}
258   else               {cellsp = sp;}
259   for (h = 0; h <= maxHeight; h++) {
260     ierr = DMPlexGetHeightStratum(dm, h, &pStart, &pEnd);CHKERRQ(ierr);
261     if (!h) {pStart = cStart; pEnd = cEnd;}
262     if (pEnd <= pStart) continue;
263     totDim = 0;
264     for (f = 0; f < numFields; ++f) {
265       PetscObject  obj;
266       PetscClassId id;
267 
268       ierr = DMGetField(dm, f, &obj);CHKERRQ(ierr);
269       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
270       if (id == PETSCFE_CLASSID) {
271         PetscFE fe = (PetscFE) obj;
272 
273         ierr = PetscFEGetNumComponents(fe, &numComp[f]);CHKERRQ(ierr);
274         if (!h) {
275           ierr = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr);
276           sp[f] = cellsp[f];
277           ierr = PetscObjectReference((PetscObject) sp[f]);CHKERRQ(ierr);
278         } else {
279           ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], h, &sp[f]);CHKERRQ(ierr);
280           if (!sp[f]) continue;
281         }
282       } else if (id == PETSCFV_CLASSID) {
283         PetscFV fv = (PetscFV) obj;
284 
285         ierr = PetscFVGetNumComponents(fv, &numComp[f]);CHKERRQ(ierr);
286         ierr = PetscFVGetDualSpace(fv, &sp[f]);CHKERRQ(ierr);
287         ierr = PetscObjectReference((PetscObject) sp[f]);CHKERRQ(ierr);
288       } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
289       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
290       totDim += spDim*numComp[f];
291     }
292     ierr = DMPlexVecGetClosure(dm, section, localX, pStart, &numValues, NULL);CHKERRQ(ierr);
293     if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section point closure size %d != dual space dimension %d", numValues, totDim);
294     if (!totDim) continue;
295     ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
296     ierr = DMGetWorkArray(dm, numFields, PETSC_BOOL, &fieldActive);CHKERRQ(ierr);
297     for (f = 0; f < numFields; ++f) fieldActive[f] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE;
298     for (i = 0; i < numIds; ++i) {
299       IS              pointIS;
300       const PetscInt *points;
301       PetscInt        n, p;
302 
303       ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr);
304       if (!pointIS) continue; /* No points with that id on this process */
305       ierr = ISGetLocalSize(pointIS, &n);CHKERRQ(ierr);
306       ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
307       for (p = 0; p < n; ++p) {
308         const PetscInt    point = points[p];
309         PetscFECellGeom   geom;
310 
311         if ((point < pStart) || (point >= pEnd)) continue;
312         ierr          = DMPlexComputeCellGeometryFEM(dm, point, NULL, geom.v0, geom.J, NULL, &geom.detJ);CHKERRQ(ierr);
313         geom.dim      = dim - h;
314         geom.dimEmbed = dimEmbed;
315         for (f = 0, v = 0; f < numFields; ++f) {
316           void * const ctx = ctxs ? ctxs[f] : NULL;
317 
318           if (!sp[f]) continue;
319           ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
320           for (d = 0; d < spDim; ++d) {
321             if (funcs[f]) {
322               ierr = PetscDualSpaceApply(sp[f], d, time, &geom, numComp[f], funcs[f], ctx, &values[v]);
323               if (ierr) {
324                 PetscErrorCode ierr2;
325                 ierr2 = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr2);
326                 ierr2 = DMRestoreWorkArray(dm, numFields, PETSC_BOOL, &fieldActive);CHKERRQ(ierr2);
327                 CHKERRQ(ierr);
328               }
329             } else {
330               for (comp = 0; comp < numComp[f]; ++comp) values[v+comp] = 0.0;
331             }
332             v += numComp[f];
333           }
334         }
335         ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, point, values, mode);CHKERRQ(ierr);
336       }
337       ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
338       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
339     }
340     if (h) {
341       for (f = 0; f < numFields; ++f) {ierr = PetscDualSpaceDestroy(&sp[f]);CHKERRQ(ierr);}
342     }
343     ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
344     ierr = DMRestoreWorkArray(dm, numFields, PETSC_BOOL, &fieldActive);CHKERRQ(ierr);
345   }
346   for (f = 0; f < numFields; ++f) {ierr = PetscDualSpaceDestroy(&sp[f]);CHKERRQ(ierr);}
347   ierr = PetscFree2(sp, numComp);CHKERRQ(ierr);
348   if (maxHeight > 0) {
349     ierr = PetscFree(cellsp);CHKERRQ(ierr);
350   }
351   PetscFunctionReturn(0);
352 }
353 
354 #undef __FUNCT__
355 #define __FUNCT__ "DMProjectFunctionLocal_Plex"
356 PetscErrorCode DMProjectFunctionLocal_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
357 {
358   PetscDualSpace *sp, *cellsp;
359   PetscInt       *numComp;
360   PetscSection    section;
361   PetscScalar    *values;
362   PetscInt        Nf, dim, dimEmbed, spDim, totDim = 0, numValues, pStart, pEnd, p, cStart, cEnd, cEndInterior, f, d, v, comp, h, maxHeight;
363   PetscBool      *isFE, hasFE = PETSC_FALSE, hasFV = PETSC_FALSE;
364   PetscErrorCode  ierr;
365 
366   PetscFunctionBegin;
367   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
368   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
369   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
370   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
371   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
372   ierr = PetscMalloc3(Nf, &isFE, Nf, &sp, Nf, &numComp);CHKERRQ(ierr);
373   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
374   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
375   ierr = DMPlexGetMaxProjectionHeight(dm,&maxHeight);CHKERRQ(ierr);
376   if (maxHeight < 0 || maxHeight > dim) {SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"maximum projection height %d not in [0, %d)\n", maxHeight,dim);}
377   if (maxHeight > 0) {
378     ierr = PetscMalloc1(Nf,&cellsp);CHKERRQ(ierr);
379   }
380   else {
381     cellsp = sp;
382   }
383   for (f = 0; f < Nf; ++f) {
384     PetscObject  obj;
385     PetscClassId id;
386 
387     ierr = DMGetField(dm, f, &obj);CHKERRQ(ierr);
388     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
389     if (id == PETSCFE_CLASSID) {
390       PetscFE fe = (PetscFE) obj;
391 
392       hasFE   = PETSC_TRUE;
393       isFE[f] = PETSC_TRUE;
394       ierr = PetscFEGetNumComponents(fe, &numComp[f]);CHKERRQ(ierr);
395       ierr  = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr);
396     } else if (id == PETSCFV_CLASSID) {
397       PetscFV fv = (PetscFV) obj;
398 
399       hasFV   = PETSC_TRUE;
400       isFE[f] = PETSC_FALSE;
401       ierr = PetscFVGetNumComponents(fv, &numComp[f]);CHKERRQ(ierr);
402       ierr = PetscFVGetDualSpace(fv, &cellsp[f]);CHKERRQ(ierr);
403     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
404   }
405   for (h = 0; h <= maxHeight; h++) {
406     ierr = DMPlexGetHeightStratum(dm, h, &pStart, &pEnd);CHKERRQ(ierr);
407     if (!h) {pStart = cStart; pEnd = cEnd;}
408     if (pEnd <= pStart) continue;
409     totDim = 0;
410     for (f = 0; f < Nf; ++f) {
411       if (!h) {
412         sp[f] = cellsp[f];
413       }
414       else {
415         ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], h, &sp[f]);CHKERRQ(ierr);
416         if (!sp[f]) {
417           continue;
418         }
419       }
420       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
421       totDim += spDim*numComp[f];
422     }
423     ierr = DMPlexVecGetClosure(dm, section, localX, pStart, &numValues, NULL);CHKERRQ(ierr);
424     if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section point closure size %d != dual space dimension %d", numValues, totDim);
425     if (!totDim) continue;
426     ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
427     for (p = pStart; p < pEnd; ++p) {
428       PetscFECellGeom fegeom;
429       PetscFVCellGeom fvgeom;
430 
431       if (hasFE) {
432         ierr = DMPlexComputeCellGeometryFEM(dm, p, NULL, fegeom.v0, fegeom.J, NULL, &fegeom.detJ);CHKERRQ(ierr);
433         fegeom.dim      = dim - h;
434         fegeom.dimEmbed = dimEmbed;
435       }
436       if (hasFV) {ierr = DMPlexComputeCellGeometryFVM(dm, p, &fvgeom.volume, fvgeom.centroid, NULL);CHKERRQ(ierr);}
437       for (f = 0, v = 0; f < Nf; ++f) {
438         void * const ctx = ctxs ? ctxs[f] : NULL;
439 
440         if (!sp[f]) continue;
441         ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
442         for (d = 0; d < spDim; ++d) {
443           if (funcs[f]) {
444             if (isFE[f]) {ierr = PetscDualSpaceApply(sp[f], d, time, &fegeom, numComp[f], funcs[f], ctx, &values[v]);}
445             else         {ierr = PetscDualSpaceApplyFVM(sp[f], d, time, &fvgeom, numComp[f], funcs[f], ctx, &values[v]);}
446             if (ierr) {
447               PetscErrorCode ierr2;
448               ierr2 = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr2);
449               CHKERRQ(ierr);
450             }
451           } else {
452             for (comp = 0; comp < numComp[f]; ++comp) values[v+comp] = 0.0;
453           }
454           v += numComp[f];
455         }
456       }
457       ierr = DMPlexVecSetClosure(dm, section, localX, p, values, mode);CHKERRQ(ierr);
458     }
459     ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
460     if (h) {
461       for (f = 0; f < Nf; f++) {ierr = PetscDualSpaceDestroy(&sp[f]);CHKERRQ(ierr);}
462     }
463   }
464   ierr = PetscFree3(isFE, sp, numComp);CHKERRQ(ierr);
465   if (maxHeight > 0) {
466     ierr = PetscFree(cellsp);CHKERRQ(ierr);
467   }
468   PetscFunctionReturn(0);
469 }
470 
471 #undef __FUNCT__
472 #define __FUNCT__ "DMProjectFieldLocal_Plex"
473 PetscErrorCode DMProjectFieldLocal_Plex(DM dm, Vec localU,
474                                         void (**funcs)(PetscInt, PetscInt, PetscInt,
475                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
476                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
477                                                        PetscReal, const PetscReal[], PetscScalar[]),
478                                         InsertMode mode, Vec localX)
479 {
480   DM              dmAux;
481   PetscDS         prob, probAux = NULL;
482   Vec             A;
483   PetscSection    section, sectionAux = NULL;
484   PetscDualSpace *sp;
485   PetscInt       *Ncf;
486   PetscScalar    *values, *u, *u_x, *a, *a_x;
487   PetscReal      *x, *v0, *J, *invJ, detJ;
488   PetscInt       *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
489   PetscInt        Nf, NfAux = 0, dim, spDim, totDim, numValues, cStart, cEnd, cEndInterior, c, f, d, v, comp, maxHeight;
490   PetscErrorCode  ierr;
491 
492   PetscFunctionBegin;
493   ierr = DMPlexGetMaxProjectionHeight(dm,&maxHeight);CHKERRQ(ierr);
494   if (maxHeight > 0) {SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Field projection for height > 0 not supported yet");}
495   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
496   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
497   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
498   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
499   ierr = PetscMalloc2(Nf, &sp, Nf, &Ncf);CHKERRQ(ierr);
500   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
501   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
502   ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr);
503   ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr);
504   ierr = PetscDSGetEvaluationArrays(prob, &u, NULL, &u_x);CHKERRQ(ierr);
505   ierr = PetscDSGetRefCoordArrays(prob, &x, NULL);CHKERRQ(ierr);
506   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
507   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
508   if (dmAux) {
509     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
510     ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr);
511     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
512     ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr);
513     ierr = PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);CHKERRQ(ierr);
514     ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);CHKERRQ(ierr);
515   }
516   ierr = DMPlexInsertBoundaryValues(dm, PETSC_TRUE, localU, 0.0, NULL, NULL, NULL);CHKERRQ(ierr);
517   ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr);
518   if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section cell closure size %d != dual space dimension %d", numValues, totDim);
519   ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
520   ierr = PetscMalloc3(dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
521   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
522   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
523   for (c = cStart; c < cEnd; ++c) {
524     PetscScalar *coefficients = NULL, *coefficientsAux = NULL;
525 
526     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
527     ierr = DMPlexVecGetClosure(dm, section, localU, c, NULL, &coefficients);CHKERRQ(ierr);
528     if (dmAux) {ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &coefficientsAux);CHKERRQ(ierr);}
529     for (f = 0, v = 0; f < Nf; ++f) {
530       PetscObject  obj;
531       PetscClassId id;
532 
533       ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
534       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
535       if (id == PETSCFE_CLASSID) {
536         PetscFE fe = (PetscFE) obj;
537 
538         ierr = PetscFEGetDualSpace(fe, &sp[f]);CHKERRQ(ierr);
539         ierr = PetscFEGetNumComponents(fe, &Ncf[f]);CHKERRQ(ierr);
540       } else if (id == PETSCFV_CLASSID) {
541         PetscFV fv = (PetscFV) obj;
542 
543         ierr = PetscFVGetNumComponents(fv, &Ncf[f]);CHKERRQ(ierr);
544         ierr = PetscObjectReference((PetscObject) sp[f]);CHKERRQ(ierr);
545       }
546       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
547       for (d = 0; d < spDim; ++d) {
548         PetscQuadrature  quad;
549         const PetscReal *points, *weights;
550         PetscInt         numPoints, q;
551 
552         if (funcs[f]) {
553           ierr = PetscDualSpaceGetFunctional(sp[f], d, &quad);CHKERRQ(ierr);
554           ierr = PetscQuadratureGetData(quad, NULL, &numPoints, &points, &weights);CHKERRQ(ierr);
555           for (q = 0; q < numPoints; ++q) {
556             CoordinatesRefToReal(dim, dim, v0, J, &points[q*dim], x);
557             ierr = EvaluateFieldJets(prob,    PETSC_FALSE, q, invJ, coefficients,    NULL, u, u_x, NULL);CHKERRQ(ierr);
558             ierr = EvaluateFieldJets(probAux, PETSC_FALSE, q, invJ, coefficientsAux, NULL, a, a_x, NULL);CHKERRQ(ierr);
559             (*funcs[f])(dim, Nf, NfAux, uOff, uOff_x, u, NULL, u_x, aOff, aOff_x, a, NULL, a_x, 0.0, x, &values[v]);
560           }
561         } else {
562           for (comp = 0; comp < Ncf[f]; ++comp) values[v+comp] = 0.0;
563         }
564         v += Ncf[f];
565       }
566     }
567     ierr = DMPlexVecRestoreClosure(dm, section, localU, c, NULL, &coefficients);CHKERRQ(ierr);
568     if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &coefficientsAux);CHKERRQ(ierr);}
569     ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr);
570   }
571   ierr = PetscFree3(v0,J,invJ);CHKERRQ(ierr);
572   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
573   ierr = PetscFree2(sp, Ncf);CHKERRQ(ierr);
574   PetscFunctionReturn(0);
575 }
576 
577 #undef __FUNCT__
578 #define __FUNCT__ "DMPlexInsertBoundaryValues_FEM_Internal"
579 static PetscErrorCode DMPlexInsertBoundaryValues_FEM_Internal(DM dm, PetscReal time, PetscInt field, DMLabel label, PetscInt numids, const PetscInt ids[], PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void *ctx, Vec locX)
580 {
581   PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
582   void            **ctxs;
583   PetscInt          numFields;
584   PetscErrorCode    ierr;
585 
586   PetscFunctionBegin;
587   ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr);
588   ierr = PetscCalloc2(numFields,&funcs,numFields,&ctxs);CHKERRQ(ierr);
589   funcs[field] = func;
590   ctxs[field]  = ctx;
591   ierr = DMProjectFunctionLabelLocal(dm, time, label, numids, ids, funcs, ctxs, INSERT_BC_VALUES, locX);CHKERRQ(ierr);
592   ierr = PetscFree2(funcs,ctxs);CHKERRQ(ierr);
593   PetscFunctionReturn(0);
594 }
595 
596 #undef __FUNCT__
597 #define __FUNCT__ "DMPlexInsertBoundaryValues_FVM_Internal"
598 /* This ignores numcomps/comps */
599 static PetscErrorCode DMPlexInsertBoundaryValues_FVM_Internal(DM dm, PetscReal time, Vec faceGeometry, Vec cellGeometry, Vec Grad,
600                                                               PetscInt field, DMLabel label, PetscInt numids, const PetscInt ids[], PetscErrorCode (*func)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*), void *ctx, Vec locX)
601 {
602   PetscDS            prob;
603   PetscSF            sf;
604   DM                 dmFace, dmCell, dmGrad;
605   const PetscScalar *facegeom, *cellgeom = NULL, *grad;
606   const PetscInt    *leaves;
607   PetscScalar       *x, *fx;
608   PetscInt           dim, nleaves, loc, fStart, fEnd, pdim, i;
609   PetscErrorCode     ierr, ierru = 0;
610 
611   PetscFunctionBegin;
612   ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
613   ierr = PetscSFGetGraph(sf, NULL, &nleaves, &leaves, NULL);CHKERRQ(ierr);
614   nleaves = PetscMax(0, nleaves);
615   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
616   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
617   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
618   ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr);
619   ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
620   if (cellGeometry) {
621     ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr);
622     ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
623   }
624   if (Grad) {
625     PetscFV fv;
626 
627     ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fv);CHKERRQ(ierr);
628     ierr = VecGetDM(Grad, &dmGrad);CHKERRQ(ierr);
629     ierr = VecGetArrayRead(Grad, &grad);CHKERRQ(ierr);
630     ierr = PetscFVGetNumComponents(fv, &pdim);CHKERRQ(ierr);
631     ierr = DMGetWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr);
632   }
633   ierr = VecGetArray(locX, &x);CHKERRQ(ierr);
634   for (i = 0; i < numids; ++i) {
635     IS              faceIS;
636     const PetscInt *faces;
637     PetscInt        numFaces, f;
638 
639     ierr = DMLabelGetStratumIS(label, ids[i], &faceIS);CHKERRQ(ierr);
640     if (!faceIS) continue; /* No points with that id on this process */
641     ierr = ISGetLocalSize(faceIS, &numFaces);CHKERRQ(ierr);
642     ierr = ISGetIndices(faceIS, &faces);CHKERRQ(ierr);
643     for (f = 0; f < numFaces; ++f) {
644       const PetscInt         face = faces[f], *cells;
645       PetscFVFaceGeom        *fg;
646 
647       if ((face < fStart) || (face >= fEnd)) continue; /* Refinement adds non-faces to labels */
648       ierr = PetscFindInt(face, nleaves, (PetscInt *) leaves, &loc);CHKERRQ(ierr);
649       if (loc >= 0) continue;
650       ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr);
651       ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr);
652       if (Grad) {
653         PetscFVCellGeom       *cg;
654         PetscScalar           *cx, *cgrad;
655         PetscScalar           *xG;
656         PetscReal              dx[3];
657         PetscInt               d;
658 
659         ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cg);CHKERRQ(ierr);
660         ierr = DMPlexPointLocalRead(dm, cells[0], x, &cx);CHKERRQ(ierr);
661         ierr = DMPlexPointLocalRead(dmGrad, cells[0], grad, &cgrad);CHKERRQ(ierr);
662         ierr = DMPlexPointLocalFieldRef(dm, cells[1], field, x, &xG);CHKERRQ(ierr);
663         DMPlex_WaxpyD_Internal(dim, -1, cg->centroid, fg->centroid, dx);
664         for (d = 0; d < pdim; ++d) fx[d] = cx[d] + DMPlex_DotD_Internal(dim, &cgrad[d*dim], dx);
665         ierru = (*func)(time, fg->centroid, fg->normal, fx, xG, ctx);
666         if (ierru) {
667           ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr);
668           ierr = ISDestroy(&faceIS);CHKERRQ(ierr);
669           goto cleanup;
670         }
671       } else {
672         PetscScalar       *xI;
673         PetscScalar       *xG;
674 
675         ierr = DMPlexPointLocalRead(dm, cells[0], x, &xI);CHKERRQ(ierr);
676         ierr = DMPlexPointLocalFieldRef(dm, cells[1], field, x, &xG);CHKERRQ(ierr);
677         ierru = (*func)(time, fg->centroid, fg->normal, xI, xG, ctx);
678         if (ierru) {
679           ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr);
680           ierr = ISDestroy(&faceIS);CHKERRQ(ierr);
681           goto cleanup;
682         }
683       }
684     }
685     ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr);
686     ierr = ISDestroy(&faceIS);CHKERRQ(ierr);
687   }
688   cleanup:
689   ierr = VecRestoreArray(locX, &x);CHKERRQ(ierr);
690   if (Grad) {
691     ierr = DMRestoreWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr);
692     ierr = VecRestoreArrayRead(Grad, &grad);CHKERRQ(ierr);
693   }
694   if (cellGeometry) {ierr = VecRestoreArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);}
695   ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
696   CHKERRQ(ierru);
697   PetscFunctionReturn(0);
698 }
699 
700 #undef __FUNCT__
701 #define __FUNCT__ "DMPlexInsertBoundaryValues"
702 PetscErrorCode DMPlexInsertBoundaryValues(DM dm, PetscBool insertEssential, Vec locX, PetscReal time, Vec faceGeomFVM, Vec cellGeomFVM, Vec gradFVM)
703 {
704   PetscInt       numBd, b;
705   PetscErrorCode ierr;
706 
707   PetscFunctionBegin;
708   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
709   PetscValidHeaderSpecific(locX, VEC_CLASSID, 2);
710   if (faceGeomFVM) {PetscValidHeaderSpecific(faceGeomFVM, VEC_CLASSID, 4);}
711   if (cellGeomFVM) {PetscValidHeaderSpecific(cellGeomFVM, VEC_CLASSID, 5);}
712   if (gradFVM)     {PetscValidHeaderSpecific(gradFVM, VEC_CLASSID, 6);}
713   ierr = DMGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
714   for (b = 0; b < numBd; ++b) {
715     PetscBool       isEssential;
716     const char     *labelname;
717     DMLabel         label;
718     PetscInt        field;
719     PetscObject     obj;
720     PetscClassId    id;
721     void          (*func)();
722     PetscInt        numids;
723     const PetscInt *ids;
724     void           *ctx;
725 
726     ierr = DMGetBoundary(dm, b, &isEssential, NULL, &labelname, &field, NULL, NULL, &func, &numids, &ids, &ctx);CHKERRQ(ierr);
727     if (insertEssential != isEssential) continue;
728     ierr = DMGetLabel(dm, labelname, &label);CHKERRQ(ierr);
729     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
730     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
731     if (id == PETSCFE_CLASSID) {
732       if (!isEssential) continue; /* for FEM, there is no insertion to be done for non-essential boundary conditions */
733       ierr = DMPlexLabelAddCells(dm,label);CHKERRQ(ierr);
734       ierr = DMPlexInsertBoundaryValues_FEM_Internal(dm, time, field, label, numids, ids, (PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *)) func, ctx, locX);CHKERRQ(ierr);
735       ierr = DMPlexLabelClearCells(dm,label);CHKERRQ(ierr);
736     } else if (id == PETSCFV_CLASSID) {
737       if (!faceGeomFVM) continue;
738       ierr = DMPlexInsertBoundaryValues_FVM_Internal(dm, time, faceGeomFVM, cellGeomFVM, gradFVM,
739                                                      field, label, numids, ids, (PetscErrorCode (*)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*)) func, ctx, locX);CHKERRQ(ierr);
740     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
741   }
742   PetscFunctionReturn(0);
743 }
744 
745 #undef __FUNCT__
746 #define __FUNCT__ "DMComputeL2Diff_Plex"
747 PetscErrorCode DMComputeL2Diff_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
748 {
749   const PetscInt   debug = 0;
750   PetscSection     section;
751   PetscQuadrature  quad;
752   Vec              localX;
753   PetscScalar     *funcVal, *interpolant;
754   PetscReal       *coords, *v0, *J, *invJ, detJ;
755   PetscReal        localDiff = 0.0;
756   const PetscReal *quadPoints, *quadWeights;
757   PetscInt         dim, numFields, numComponents = 0, numQuadPoints, cStart, cEnd, cEndInterior, c, field, fieldOffset;
758   PetscErrorCode   ierr;
759 
760   PetscFunctionBegin;
761   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
762   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
763   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
764   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
765   ierr = DMProjectFunctionLocal(dm, time, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
766   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
767   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
768   for (field = 0; field < numFields; ++field) {
769     PetscObject  obj;
770     PetscClassId id;
771     PetscInt     Nc;
772 
773     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
774     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
775     if (id == PETSCFE_CLASSID) {
776       PetscFE fe = (PetscFE) obj;
777 
778       ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
779       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
780     } else if (id == PETSCFV_CLASSID) {
781       PetscFV fv = (PetscFV) obj;
782 
783       ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr);
784       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
785     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
786     numComponents += Nc;
787   }
788   ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
789   ierr = PetscMalloc6(numComponents,&funcVal,numComponents,&interpolant,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
790   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
791   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
792   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
793   for (c = cStart; c < cEnd; ++c) {
794     PetscScalar *x = NULL;
795     PetscReal    elemDiff = 0.0;
796 
797     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
798     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
799     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
800 
801     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
802       PetscObject  obj;
803       PetscClassId id;
804       void * const ctx = ctxs ? ctxs[field] : NULL;
805       PetscInt     Nb, Nc, q, fc;
806 
807       ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
808       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
809       if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);}
810       else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;}
811       else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
812       if (debug) {
813         char title[1024];
814         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
815         ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr);
816       }
817       for (q = 0; q < numQuadPoints; ++q) {
818         CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], coords);
819         ierr = (*funcs[field])(dim, time, coords, Nc, funcVal, ctx);
820         if (ierr) {
821           PetscErrorCode ierr2;
822           ierr2 = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr2);
823           ierr2 = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr2);
824           ierr2 = PetscFree6(funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr2);
825           CHKERRQ(ierr);
826         }
827         if (id == PETSCFE_CLASSID)      {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
828         else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
829         else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
830         for (fc = 0; fc < Nc; ++fc) {
831           if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "    elem %d field %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ);CHKERRQ(ierr);}
832           elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ;
833         }
834       }
835       fieldOffset += Nb*Nc;
836     }
837     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
838     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
839     localDiff += elemDiff;
840   }
841   ierr  = PetscFree6(funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr);
842   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
843   ierr  = MPIU_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
844   *diff = PetscSqrtReal(*diff);
845   PetscFunctionReturn(0);
846 }
847 
848 #undef __FUNCT__
849 #define __FUNCT__ "DMComputeL2GradientDiff_Plex"
850 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)
851 {
852   const PetscInt  debug = 0;
853   PetscSection    section;
854   PetscQuadrature quad;
855   Vec             localX;
856   PetscScalar    *funcVal, *interpolantVec;
857   PetscReal      *coords, *realSpaceDer, *v0, *J, *invJ, detJ;
858   PetscReal       localDiff = 0.0;
859   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, cEndInterior, c, field, fieldOffset, comp;
860   PetscErrorCode  ierr;
861 
862   PetscFunctionBegin;
863   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
864   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
865   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
866   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
867   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
868   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
869   for (field = 0; field < numFields; ++field) {
870     PetscFE  fe;
871     PetscInt Nc;
872 
873     ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
874     ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
875     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
876     numComponents += Nc;
877   }
878   /* ierr = DMProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */
879   ierr = PetscMalloc7(numComponents,&funcVal,dim,&coords,dim,&realSpaceDer,dim,&v0,dim*dim,&J,dim*dim,&invJ,dim,&interpolantVec);CHKERRQ(ierr);
880   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
881   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
882   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
883   for (c = cStart; c < cEnd; ++c) {
884     PetscScalar *x = NULL;
885     PetscReal    elemDiff = 0.0;
886 
887     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
888     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
889     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
890 
891     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
892       PetscFE          fe;
893       void * const     ctx = ctxs ? ctxs[field] : NULL;
894       const PetscReal *quadPoints, *quadWeights;
895       PetscReal       *basisDer;
896       PetscInt         numQuadPoints, Nb, Ncomp, q, d, e, fc, f, g;
897 
898       ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
899       ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
900       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
901       ierr = PetscFEGetNumComponents(fe, &Ncomp);CHKERRQ(ierr);
902       ierr = PetscFEGetDefaultTabulation(fe, NULL, &basisDer, NULL);CHKERRQ(ierr);
903       if (debug) {
904         char title[1024];
905         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
906         ierr = DMPrintCellVector(c, title, Nb*Ncomp, &x[fieldOffset]);CHKERRQ(ierr);
907       }
908       for (q = 0; q < numQuadPoints; ++q) {
909         for (d = 0; d < dim; d++) {
910           coords[d] = v0[d];
911           for (e = 0; e < dim; e++) {
912             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
913           }
914         }
915         ierr = (*funcs[field])(dim, time, coords, n, numFields, funcVal, ctx);
916         if (ierr) {
917           PetscErrorCode ierr2;
918           ierr2 = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr2);
919           ierr2 = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr2);
920           ierr2 = PetscFree7(funcVal,coords,realSpaceDer,v0,J,invJ,interpolantVec);CHKERRQ(ierr2);
921           CHKERRQ(ierr);
922         }
923         for (fc = 0; fc < Ncomp; ++fc) {
924           PetscScalar interpolant = 0.0;
925 
926           for (d = 0; d < dim; ++d) interpolantVec[d] = 0.0;
927           for (f = 0; f < Nb; ++f) {
928             const PetscInt fidx = f*Ncomp+fc;
929 
930             for (d = 0; d < dim; ++d) {
931               realSpaceDer[d] = 0.0;
932               for (g = 0; g < dim; ++g) {
933                 realSpaceDer[d] += invJ[g*dim+d]*basisDer[(q*Nb*Ncomp+fidx)*dim+g];
934               }
935               interpolantVec[d] += x[fieldOffset+fidx]*realSpaceDer[d];
936             }
937           }
938           for (d = 0; d < dim; ++d) interpolant += interpolantVec[d]*n[d];
939           if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "    elem %d fieldDer %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ);CHKERRQ(ierr);}
940           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ;
941         }
942       }
943       comp        += Ncomp;
944       fieldOffset += Nb*Ncomp;
945     }
946     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
947     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
948     localDiff += elemDiff;
949   }
950   ierr  = PetscFree7(funcVal,coords,realSpaceDer,v0,J,invJ,interpolantVec);CHKERRQ(ierr);
951   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
952   ierr  = MPIU_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
953   *diff = PetscSqrtReal(*diff);
954   PetscFunctionReturn(0);
955 }
956 
957 #undef __FUNCT__
958 #define __FUNCT__ "DMComputeL2FieldDiff_Plex"
959 PetscErrorCode DMComputeL2FieldDiff_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
960 {
961   const PetscInt   debug = 0;
962   PetscSection     section;
963   PetscQuadrature  quad;
964   Vec              localX;
965   PetscScalar     *funcVal, *interpolant;
966   PetscReal       *coords, *v0, *J, *invJ, detJ;
967   PetscReal       *localDiff;
968   const PetscReal *quadPoints, *quadWeights;
969   PetscInt         dim, numFields, numComponents = 0, numQuadPoints, cStart, cEnd, cEndInterior, c, field, fieldOffset;
970   PetscErrorCode   ierr;
971 
972   PetscFunctionBegin;
973   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
974   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
975   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
976   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
977   ierr = DMProjectFunctionLocal(dm, time, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
978   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
979   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
980   for (field = 0; field < numFields; ++field) {
981     PetscObject  obj;
982     PetscClassId id;
983     PetscInt     Nc;
984 
985     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
986     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
987     if (id == PETSCFE_CLASSID) {
988       PetscFE fe = (PetscFE) obj;
989 
990       ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
991       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
992     } else if (id == PETSCFV_CLASSID) {
993       PetscFV fv = (PetscFV) obj;
994 
995       ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr);
996       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
997     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
998     numComponents += Nc;
999   }
1000   ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
1001   ierr = PetscCalloc7(numFields,&localDiff,numComponents,&funcVal,numComponents,&interpolant,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
1002   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1003   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
1004   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
1005   for (c = cStart; c < cEnd; ++c) {
1006     PetscScalar *x = NULL;
1007 
1008     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
1009     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
1010     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
1011 
1012     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
1013       PetscObject  obj;
1014       PetscClassId id;
1015       void * const ctx = ctxs ? ctxs[field] : NULL;
1016       PetscInt     Nb, Nc, q, fc;
1017 
1018       PetscReal       elemDiff = 0.0;
1019 
1020       ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
1021       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1022       if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);}
1023       else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;}
1024       else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
1025       if (debug) {
1026         char title[1024];
1027         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
1028         ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr);
1029       }
1030       for (q = 0; q < numQuadPoints; ++q) {
1031         CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], coords);
1032         ierr = (*funcs[field])(dim, time, coords, numFields, funcVal, ctx);
1033         if (ierr) {
1034           PetscErrorCode ierr2;
1035           ierr2 = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr2);
1036           ierr2 = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr2);
1037           ierr2 = PetscFree7(localDiff,funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr2);
1038           CHKERRQ(ierr);
1039         }
1040         if (id == PETSCFE_CLASSID)      {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
1041         else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
1042         else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
1043         for (fc = 0; fc < Nc; ++fc) {
1044           if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "    elem %d field %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ);CHKERRQ(ierr);}
1045           elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ;
1046         }
1047       }
1048       fieldOffset += Nb*Nc;
1049       localDiff[field] += elemDiff;
1050     }
1051     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
1052   }
1053   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
1054   ierr = MPIU_Allreduce(localDiff, diff, numFields, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
1055   for (field = 0; field < numFields; ++field) diff[field] = PetscSqrtReal(diff[field]);
1056   ierr = PetscFree7(localDiff,funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr);
1057   PetscFunctionReturn(0);
1058 }
1059 
1060 #undef __FUNCT__
1061 #define __FUNCT__ "DMPlexComputeL2DiffVec"
1062 /*@C
1063   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.
1064 
1065   Input Parameters:
1066 + dm    - The DM
1067 . time  - The time
1068 . funcs - The functions to evaluate for each field component: NULL means that component does not contribute to error calculation
1069 . ctxs  - Optional array of contexts to pass to each function, or NULL.
1070 - X     - The coefficient vector u_h
1071 
1072   Output Parameter:
1073 . D - A Vec which holds the difference ||u - u_h||_2 for each cell
1074 
1075   Level: developer
1076 
1077 .seealso: DMProjectFunction(), DMComputeL2Diff(), DMPlexComputeL2FieldDiff(), DMComputeL2GradientDiff()
1078 @*/
1079 PetscErrorCode DMPlexComputeL2DiffVec(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, Vec D)
1080 {
1081   PetscSection     section;
1082   PetscQuadrature  quad;
1083   Vec              localX;
1084   PetscScalar     *funcVal, *interpolant;
1085   PetscReal       *coords, *v0, *J, *invJ, detJ;
1086   const PetscReal *quadPoints, *quadWeights;
1087   PetscInt         dim, numFields, numComponents = 0, numQuadPoints, cStart, cEnd, cEndInterior, c, field, fieldOffset;
1088   PetscErrorCode   ierr;
1089 
1090   PetscFunctionBegin;
1091   ierr = VecSet(D, 0.0);CHKERRQ(ierr);
1092   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1093   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1094   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
1095   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
1096   ierr = DMProjectFunctionLocal(dm, time, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
1097   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
1098   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
1099   for (field = 0; field < numFields; ++field) {
1100     PetscObject  obj;
1101     PetscClassId id;
1102     PetscInt     Nc;
1103 
1104     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
1105     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1106     if (id == PETSCFE_CLASSID) {
1107       PetscFE fe = (PetscFE) obj;
1108 
1109       ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
1110       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1111     } else if (id == PETSCFV_CLASSID) {
1112       PetscFV fv = (PetscFV) obj;
1113 
1114       ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr);
1115       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
1116     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
1117     numComponents += Nc;
1118   }
1119   ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
1120   ierr = PetscMalloc6(numComponents,&funcVal,numComponents,&interpolant,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
1121   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1122   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
1123   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
1124   for (c = cStart; c < cEnd; ++c) {
1125     PetscScalar *x = NULL;
1126     PetscScalar  elemDiff = 0.0;
1127 
1128     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
1129     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
1130     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
1131 
1132     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
1133       PetscObject  obj;
1134       PetscClassId id;
1135       void * const ctx = ctxs ? ctxs[field] : NULL;
1136       PetscInt     Nb, Nc, q, fc;
1137 
1138       ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
1139       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1140       if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);}
1141       else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;}
1142       else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
1143       if (funcs[field]) {
1144         for (q = 0; q < numQuadPoints; ++q) {
1145           CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], coords);
1146           ierr = (*funcs[field])(dim, time, coords, Nc, funcVal, ctx);
1147           if (ierr) {
1148             PetscErrorCode ierr2;
1149             ierr2 = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr2);
1150             ierr2 = PetscFree6(funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr2);
1151             ierr2 = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr2);
1152             CHKERRQ(ierr);
1153           }
1154           if (id == PETSCFE_CLASSID)      {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
1155           else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
1156           else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
1157           for (fc = 0; fc < Nc; ++fc) {
1158             elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ;
1159           }
1160         }
1161       }
1162       fieldOffset += Nb*Nc;
1163     }
1164     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
1165     ierr = VecSetValue(D, c - cStart, elemDiff, INSERT_VALUES);CHKERRQ(ierr);
1166   }
1167   ierr = PetscFree6(funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr);
1168   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
1169   ierr = VecSqrtAbs(D);CHKERRQ(ierr);
1170   PetscFunctionReturn(0);
1171 }
1172 
1173 #undef __FUNCT__
1174 #define __FUNCT__ "DMPlexComputeIntegralFEM"
1175 /*@
1176   DMPlexComputeIntegralFEM - Form the local integral F from the local input X using pointwise functions specified by the user
1177 
1178   Input Parameters:
1179 + dm - The mesh
1180 . X  - Local input vector
1181 - user - The user context
1182 
1183   Output Parameter:
1184 . integral - Local integral for each field
1185 
1186   Level: developer
1187 
1188 .seealso: DMPlexComputeResidualFEM()
1189 @*/
1190 PetscErrorCode DMPlexComputeIntegralFEM(DM dm, Vec X, PetscReal *integral, void *user)
1191 {
1192   DM_Plex          *mesh  = (DM_Plex *) dm->data;
1193   DM                dmAux;
1194   Vec               localX, A;
1195   PetscDS           prob, probAux = NULL;
1196   PetscSection      section, sectionAux;
1197   PetscFECellGeom  *cgeom;
1198   PetscScalar      *u, *a = NULL;
1199   PetscReal        *lintegral, *vol;
1200   PetscInt          dim, Nf, f, numCells, cStart, cEnd, cEndInterior, c;
1201   PetscInt          totDim, totDimAux;
1202   PetscErrorCode    ierr;
1203 
1204   PetscFunctionBegin;
1205   ierr = PetscLogEventBegin(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);
1206   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
1207   ierr = DMPlexInsertBoundaryValues(dm, PETSC_TRUE, localX, 0.0, NULL, NULL, NULL);CHKERRQ(ierr);
1208   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
1209   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
1210   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1211   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1212   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
1213   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
1214   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
1215   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1216   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
1217   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
1218   numCells = cEnd - cStart;
1219   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
1220   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
1221   if (dmAux) {
1222     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
1223     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
1224     ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
1225   }
1226   ierr = PetscMalloc4(Nf,&lintegral,numCells*totDim,&u,numCells,&cgeom,numCells,&vol);CHKERRQ(ierr);
1227   if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);}
1228   for (f = 0; f < Nf; ++f) {lintegral[f] = 0.0;}
1229   for (c = cStart; c < cEnd; ++c) {
1230     PetscScalar *x = NULL;
1231     PetscInt     i;
1232 
1233     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, cgeom[c].v0, cgeom[c].J, cgeom[c].invJ, &cgeom[c].detJ);CHKERRQ(ierr);
1234     ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol[c], NULL, NULL);CHKERRQ(ierr);
1235     if (cgeom[c].detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", cgeom[c].detJ, c);
1236     ierr = DMPlexVecGetClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr);
1237     for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i];
1238     ierr = DMPlexVecRestoreClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr);
1239     if (dmAux) {
1240       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
1241       for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i];
1242       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
1243     }
1244   }
1245   for (f = 0; f < Nf; ++f) {
1246     PetscObject  obj;
1247     PetscClassId id;
1248     PetscInt     numChunks, numBatches, batchSize, numBlocks, blockSize, Ne, Nr, offset;
1249 
1250     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
1251     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1252     if (id == PETSCFE_CLASSID) {
1253       PetscFE         fe = (PetscFE) obj;
1254       PetscQuadrature q;
1255       PetscInt        Nq, Nb;
1256 
1257       ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
1258       ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
1259       ierr = PetscQuadratureGetData(q, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);
1260       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
1261       blockSize = Nb*Nq;
1262       batchSize = numBlocks * blockSize;
1263       ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1264       numChunks = numCells / (numBatches*batchSize);
1265       Ne        = numChunks*numBatches*batchSize;
1266       Nr        = numCells % (numBatches*batchSize);
1267       offset    = numCells - Nr;
1268       ierr = PetscFEIntegrate(fe, prob, f, Ne, cgeom, u, probAux, a, lintegral);CHKERRQ(ierr);
1269       ierr = PetscFEIntegrate(fe, prob, f, Nr, &cgeom[offset], &u[offset*totDim], probAux, &a[offset*totDimAux], lintegral);CHKERRQ(ierr);
1270     } else if (id == PETSCFV_CLASSID) {
1271       /* PetscFV  fv = (PetscFV) obj; */
1272       PetscInt       foff;
1273       PetscPointFunc obj_func;
1274       PetscScalar    lint;
1275 
1276       ierr = PetscDSGetObjective(prob, f, &obj_func);CHKERRQ(ierr);
1277       ierr = PetscDSGetFieldOffset(prob, f, &foff);CHKERRQ(ierr);
1278       if (obj_func) {
1279         for (c = 0; c < numCells; ++c) {
1280           /* TODO: Need full pointwise interpolation and get centroid for x argument */
1281           obj_func(dim, Nf, 0, NULL, NULL, &u[totDim*c+foff], NULL, NULL, NULL, NULL, NULL, NULL, NULL, 0.0, NULL, &lint);
1282           lintegral[f] = PetscRealPart(lint)*vol[c];
1283         }
1284       }
1285     } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
1286   }
1287   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
1288   if (mesh->printFEM) {
1289     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Local integral:");CHKERRQ(ierr);
1290     for (f = 0; f < Nf; ++f) {ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), " %g", lintegral[f]);CHKERRQ(ierr);}
1291     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "\n");CHKERRQ(ierr);
1292   }
1293   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
1294   ierr = MPIU_Allreduce(lintegral, integral, Nf, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject) dm));CHKERRQ(ierr);
1295   ierr = PetscFree4(lintegral,u,cgeom,vol);CHKERRQ(ierr);
1296   ierr = PetscLogEventEnd(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);
1297   PetscFunctionReturn(0);
1298 }
1299 
1300 #undef __FUNCT__
1301 #define __FUNCT__ "DMPlexComputeInterpolatorNested"
1302 /*@
1303   DMPlexComputeInterpolatorNested - Form the local portion of the interpolation matrix I from the coarse DM to the uniformly refined DM.
1304 
1305   Input Parameters:
1306 + dmf  - The fine mesh
1307 . dmc  - The coarse mesh
1308 - user - The user context
1309 
1310   Output Parameter:
1311 . In  - The interpolation matrix
1312 
1313   Level: developer
1314 
1315 .seealso: DMPlexComputeInterpolatorGeneral(), DMPlexComputeJacobianFEM()
1316 @*/
1317 PetscErrorCode DMPlexComputeInterpolatorNested(DM dmc, DM dmf, Mat In, void *user)
1318 {
1319   DM_Plex          *mesh  = (DM_Plex *) dmc->data;
1320   const char       *name  = "Interpolator";
1321   PetscDS           prob;
1322   PetscFE          *feRef;
1323   PetscFV          *fvRef;
1324   PetscSection      fsection, fglobalSection;
1325   PetscSection      csection, cglobalSection;
1326   PetscScalar      *elemMat;
1327   PetscInt          dim, Nf, f, fieldI, fieldJ, offsetI, offsetJ, cStart, cEnd, cEndInterior, c;
1328   PetscInt          cTotDim, rTotDim = 0;
1329   PetscErrorCode    ierr;
1330 
1331   PetscFunctionBegin;
1332   ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1333   ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr);
1334   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
1335   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
1336   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
1337   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
1338   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
1339   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
1340   ierr = DMPlexGetHybridBounds(dmc, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
1341   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
1342   ierr = DMGetDS(dmf, &prob);CHKERRQ(ierr);
1343   ierr = PetscCalloc2(Nf,&feRef,Nf,&fvRef);CHKERRQ(ierr);
1344   for (f = 0; f < Nf; ++f) {
1345     PetscObject  obj;
1346     PetscClassId id;
1347     PetscInt     rNb = 0, Nc = 0;
1348 
1349     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
1350     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1351     if (id == PETSCFE_CLASSID) {
1352       PetscFE fe = (PetscFE) obj;
1353 
1354       ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr);
1355       ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr);
1356       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1357     } else if (id == PETSCFV_CLASSID) {
1358       PetscFV        fv = (PetscFV) obj;
1359       PetscDualSpace Q;
1360 
1361       ierr = PetscFVRefine(fv, &fvRef[f]);CHKERRQ(ierr);
1362       ierr = PetscFVGetDualSpace(fvRef[f], &Q);CHKERRQ(ierr);
1363       ierr = PetscDualSpaceGetDimension(Q, &rNb);CHKERRQ(ierr);
1364       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
1365     }
1366     rTotDim += rNb*Nc;
1367   }
1368   ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr);
1369   ierr = PetscMalloc1(rTotDim*cTotDim,&elemMat);CHKERRQ(ierr);
1370   ierr = PetscMemzero(elemMat, rTotDim*cTotDim * sizeof(PetscScalar));CHKERRQ(ierr);
1371   for (fieldI = 0, offsetI = 0; fieldI < Nf; ++fieldI) {
1372     PetscDualSpace   Qref;
1373     PetscQuadrature  f;
1374     const PetscReal *qpoints, *qweights;
1375     PetscReal       *points;
1376     PetscInt         npoints = 0, Nc, Np, fpdim, i, k, p, d;
1377 
1378     /* Compose points from all dual basis functionals */
1379     if (feRef[fieldI]) {
1380       ierr = PetscFEGetDualSpace(feRef[fieldI], &Qref);CHKERRQ(ierr);
1381       ierr = PetscFEGetNumComponents(feRef[fieldI], &Nc);CHKERRQ(ierr);
1382     } else {
1383       ierr = PetscFVGetDualSpace(fvRef[fieldI], &Qref);CHKERRQ(ierr);
1384       ierr = PetscFVGetNumComponents(fvRef[fieldI], &Nc);CHKERRQ(ierr);
1385     }
1386     ierr = PetscDualSpaceGetDimension(Qref, &fpdim);CHKERRQ(ierr);
1387     for (i = 0; i < fpdim; ++i) {
1388       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1389       ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, NULL);CHKERRQ(ierr);
1390       npoints += Np;
1391     }
1392     ierr = PetscMalloc1(npoints*dim,&points);CHKERRQ(ierr);
1393     for (i = 0, k = 0; i < fpdim; ++i) {
1394       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1395       ierr = PetscQuadratureGetData(f, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr);
1396       for (p = 0; p < Np; ++p, ++k) for (d = 0; d < dim; ++d) points[k*dim+d] = qpoints[p*dim+d];
1397     }
1398 
1399     for (fieldJ = 0, offsetJ = 0; fieldJ < Nf; ++fieldJ) {
1400       PetscObject  obj;
1401       PetscClassId id;
1402       PetscReal   *B;
1403       PetscInt     NcJ = 0, cpdim = 0, j;
1404 
1405       ierr = PetscDSGetDiscretization(prob, fieldJ, &obj);CHKERRQ(ierr);
1406       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1407       if (id == PETSCFE_CLASSID) {
1408         PetscFE fe = (PetscFE) obj;
1409 
1410         /* Evaluate basis at points */
1411         ierr = PetscFEGetNumComponents(fe, &NcJ);CHKERRQ(ierr);
1412         ierr = PetscFEGetDimension(fe, &cpdim);CHKERRQ(ierr);
1413         /* For now, fields only interpolate themselves */
1414         if (fieldI == fieldJ) {
1415           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);
1416           ierr = PetscFEGetTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);
1417           for (i = 0, k = 0; i < fpdim; ++i) {
1418             ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1419             ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, &qweights);CHKERRQ(ierr);
1420             for (p = 0; p < Np; ++p, ++k) {
1421               for (j = 0; j < cpdim; ++j) {
1422                 for (c = 0; c < Nc; ++c) elemMat[(offsetI + i*Nc + c)*cTotDim + offsetJ + j*NcJ + c] += B[k*cpdim*NcJ+j*Nc+c]*qweights[p];
1423               }
1424             }
1425           }
1426           ierr = PetscFERestoreTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr);
1427         }
1428       } else if (id == PETSCFV_CLASSID) {
1429         PetscFV        fv = (PetscFV) obj;
1430 
1431         /* Evaluate constant function at points */
1432         ierr = PetscFVGetNumComponents(fv, &NcJ);CHKERRQ(ierr);
1433         cpdim = 1;
1434         /* For now, fields only interpolate themselves */
1435         if (fieldI == fieldJ) {
1436           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);
1437           for (i = 0, k = 0; i < fpdim; ++i) {
1438             ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1439             ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, &qweights);CHKERRQ(ierr);
1440             for (p = 0; p < Np; ++p, ++k) {
1441               for (j = 0; j < cpdim; ++j) {
1442                 for (c = 0; c < Nc; ++c) elemMat[(offsetI + i*Nc + c)*cTotDim + offsetJ + j*NcJ + c] += 1.0*qweights[p];
1443               }
1444             }
1445           }
1446         }
1447       }
1448       offsetJ += cpdim*NcJ;
1449     }
1450     offsetI += fpdim*Nc;
1451     ierr = PetscFree(points);CHKERRQ(ierr);
1452   }
1453   if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(0, name, rTotDim, cTotDim, elemMat);CHKERRQ(ierr);}
1454   /* Preallocate matrix */
1455   {
1456     Mat          preallocator;
1457     PetscScalar *vals;
1458     PetscInt    *cellCIndices, *cellFIndices;
1459     PetscInt     locRows, locCols, cell;
1460 
1461     ierr = MatGetLocalSize(In, &locRows, &locCols);CHKERRQ(ierr);
1462     ierr = MatCreate(PetscObjectComm((PetscObject) In), &preallocator);CHKERRQ(ierr);
1463     ierr = MatSetType(preallocator, MATPREALLOCATOR);CHKERRQ(ierr);
1464     ierr = MatSetSizes(preallocator, locRows, locCols, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
1465     ierr = MatSetUp(preallocator);CHKERRQ(ierr);
1466     ierr = PetscCalloc3(rTotDim*cTotDim, &vals,cTotDim,&cellCIndices,rTotDim,&cellFIndices);CHKERRQ(ierr);
1467     for (cell = cStart; cell < cEnd; ++cell) {
1468       ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, cell, cellCIndices, cellFIndices);CHKERRQ(ierr);
1469       ierr = MatSetValues(preallocator, rTotDim, cellFIndices, cTotDim, cellCIndices, vals, INSERT_VALUES);CHKERRQ(ierr);
1470     }
1471     ierr = PetscFree3(vals,cellCIndices,cellFIndices);CHKERRQ(ierr);
1472     ierr = MatAssemblyBegin(preallocator, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1473     ierr = MatAssemblyEnd(preallocator, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1474     ierr = MatPreallocatorPreallocate(preallocator, PETSC_TRUE, In);CHKERRQ(ierr);
1475     ierr = MatDestroy(&preallocator);CHKERRQ(ierr);
1476   }
1477   /* Fill matrix */
1478   ierr = MatZeroEntries(In);CHKERRQ(ierr);
1479   for (c = cStart; c < cEnd; ++c) {
1480     ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES);CHKERRQ(ierr);
1481   }
1482   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);}
1483   ierr = PetscFree2(feRef,fvRef);CHKERRQ(ierr);
1484   ierr = PetscFree(elemMat);CHKERRQ(ierr);
1485   ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1486   ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1487   if (mesh->printFEM) {
1488     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
1489     ierr = MatChop(In, 1.0e-10);CHKERRQ(ierr);
1490     ierr = MatView(In, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1491   }
1492   ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1493   PetscFunctionReturn(0);
1494 }
1495 
1496 #undef __FUNCT__
1497 #define __FUNCT__ "DMPlexComputeInterpolatorGeneral"
1498 /*@
1499   DMPlexComputeInterpolatorGeneral - Form the local portion of the interpolation matrix I from the coarse DM to a non-nested fine DM.
1500 
1501   Input Parameters:
1502 + dmf  - The fine mesh
1503 . dmc  - The coarse mesh
1504 - user - The user context
1505 
1506   Output Parameter:
1507 . In  - The interpolation matrix
1508 
1509   Level: developer
1510 
1511 .seealso: DMPlexComputeInterpolatorNested(), DMPlexComputeJacobianFEM()
1512 @*/
1513 PetscErrorCode DMPlexComputeInterpolatorGeneral(DM dmc, DM dmf, Mat In, void *user)
1514 {
1515   DM_Plex       *mesh = (DM_Plex *) dmf->data;
1516   const char    *name = "Interpolator";
1517   PetscDS        prob;
1518   PetscSection   fsection, csection, globalFSection, globalCSection;
1519   PetscHashJK    ht;
1520   PetscLayout    rLayout;
1521   PetscInt      *dnz, *onz;
1522   PetscInt       locRows, rStart, rEnd;
1523   PetscReal     *x, *v0, *J, *invJ, detJ;
1524   PetscReal     *v0c, *Jc, *invJc, detJc;
1525   PetscScalar   *elemMat;
1526   PetscInt       dim, Nf, field, totDim, cStart, cEnd, cell, ccell;
1527   PetscErrorCode ierr;
1528 
1529   PetscFunctionBegin;
1530   ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1531   ierr = DMGetCoordinateDim(dmc, &dim);CHKERRQ(ierr);
1532   ierr = DMGetDS(dmc, &prob);CHKERRQ(ierr);
1533   ierr = PetscDSGetRefCoordArrays(prob, &x, NULL);CHKERRQ(ierr);
1534   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
1535   ierr = PetscMalloc3(dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
1536   ierr = PetscMalloc3(dim,&v0c,dim*dim,&Jc,dim*dim,&invJc);CHKERRQ(ierr);
1537   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
1538   ierr = DMGetDefaultGlobalSection(dmf, &globalFSection);CHKERRQ(ierr);
1539   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
1540   ierr = DMGetDefaultGlobalSection(dmc, &globalCSection);CHKERRQ(ierr);
1541   ierr = DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd);CHKERRQ(ierr);
1542   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
1543   ierr = PetscMalloc1(totDim*totDim, &elemMat);CHKERRQ(ierr);
1544 
1545   ierr = MatGetLocalSize(In, &locRows, NULL);CHKERRQ(ierr);
1546   ierr = PetscLayoutCreate(PetscObjectComm((PetscObject) In), &rLayout);CHKERRQ(ierr);
1547   ierr = PetscLayoutSetLocalSize(rLayout, locRows);CHKERRQ(ierr);
1548   ierr = PetscLayoutSetBlockSize(rLayout, 1);CHKERRQ(ierr);
1549   ierr = PetscLayoutSetUp(rLayout);CHKERRQ(ierr);
1550   ierr = PetscLayoutGetRange(rLayout, &rStart, &rEnd);CHKERRQ(ierr);
1551   ierr = PetscLayoutDestroy(&rLayout);CHKERRQ(ierr);
1552   ierr = PetscCalloc2(locRows,&dnz,locRows,&onz);CHKERRQ(ierr);
1553   ierr = PetscHashJKCreate(&ht);CHKERRQ(ierr);
1554   for (field = 0; field < Nf; ++field) {
1555     PetscObject      obj;
1556     PetscClassId     id;
1557     PetscDualSpace   Q = NULL;
1558     PetscQuadrature  f;
1559     const PetscReal *qpoints;
1560     PetscInt         Nc, Np, fpdim, i, d;
1561 
1562     ierr = PetscDSGetDiscretization(prob, field, &obj);CHKERRQ(ierr);
1563     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1564     if (id == PETSCFE_CLASSID) {
1565       PetscFE fe = (PetscFE) obj;
1566 
1567       ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr);
1568       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1569     } else if (id == PETSCFV_CLASSID) {
1570       PetscFV fv = (PetscFV) obj;
1571 
1572       ierr = PetscFVGetDualSpace(fv, &Q);CHKERRQ(ierr);
1573       Nc   = 1;
1574     }
1575     ierr = PetscDualSpaceGetDimension(Q, &fpdim);CHKERRQ(ierr);
1576     /* For each fine grid cell */
1577     for (cell = cStart; cell < cEnd; ++cell) {
1578       PetscInt *findices,   *cindices;
1579       PetscInt  numFIndices, numCIndices;
1580 
1581       ierr = DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr);
1582       ierr = DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
1583       if (numFIndices != fpdim*Nc) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of fine indices %d != %d dual basis vecs", numFIndices, fpdim*Nc);
1584       for (i = 0; i < fpdim; ++i) {
1585         Vec             pointVec;
1586         PetscScalar    *pV;
1587         IS              coarseCellIS;
1588         const PetscInt *coarseCells;
1589         PetscInt        numCoarseCells, q, r, c;
1590 
1591         /* Get points from the dual basis functional quadrature */
1592         ierr = PetscDualSpaceGetFunctional(Q, i, &f);CHKERRQ(ierr);
1593         ierr = PetscQuadratureGetData(f, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr);
1594         ierr = VecCreateSeq(PETSC_COMM_SELF, Np*dim, &pointVec);CHKERRQ(ierr);
1595         ierr = VecSetBlockSize(pointVec, dim);CHKERRQ(ierr);
1596         ierr = VecGetArray(pointVec, &pV);CHKERRQ(ierr);
1597         for (q = 0; q < Np; ++q) {
1598           /* Transform point to real space */
1599           CoordinatesRefToReal(dim, dim, v0, J, &qpoints[q*dim], x);
1600           for (d = 0; d < dim; ++d) pV[q*dim+d] = x[d];
1601         }
1602         ierr = VecRestoreArray(pointVec, &pV);CHKERRQ(ierr);
1603         /* Get set of coarse cells that overlap points (would like to group points by coarse cell) */
1604         ierr = DMLocatePoints(dmc, pointVec, &coarseCellIS);CHKERRQ(ierr);
1605         ierr = ISViewFromOptions(coarseCellIS, NULL, "-interp_is_view");CHKERRQ(ierr);
1606         /* Update preallocation info */
1607         ierr = ISGetLocalSize(coarseCellIS, &numCoarseCells);CHKERRQ(ierr);
1608         ierr = ISGetIndices(coarseCellIS, &coarseCells);CHKERRQ(ierr);
1609         for (r = 0; r < Nc; ++r) {
1610           PetscHashJKKey  key;
1611           PetscHashJKIter missing, iter;
1612 
1613           key.j = findices[i*Nc+r];
1614           if (key.j < 0) continue;
1615           /* Get indices for coarse elements */
1616           for (ccell = 0; ccell < numCoarseCells; ++ccell) {
1617             ierr = DMPlexGetClosureIndices(dmc, csection, globalCSection, coarseCells[ccell], &numCIndices, &cindices, NULL);CHKERRQ(ierr);
1618             for (c = 0; c < numCIndices; ++c) {
1619               key.k = cindices[c];
1620               if (key.k < 0) continue;
1621               ierr = PetscHashJKPut(ht, key, &missing, &iter);CHKERRQ(ierr);
1622               if (missing) {
1623                 ierr = PetscHashJKSet(ht, iter, 1);CHKERRQ(ierr);
1624                 if ((key.k >= rStart) && (key.k < rEnd)) ++dnz[key.j-rStart];
1625                 else                                     ++onz[key.j-rStart];
1626               }
1627             }
1628             ierr = DMPlexRestoreClosureIndices(dmc, csection, globalCSection, coarseCells[ccell], &numCIndices, &cindices, NULL);CHKERRQ(ierr);
1629           }
1630         }
1631         ierr = ISRestoreIndices(coarseCellIS, &coarseCells);CHKERRQ(ierr);
1632         ierr = ISDestroy(&coarseCellIS);CHKERRQ(ierr);
1633         ierr = VecDestroy(&pointVec);CHKERRQ(ierr);
1634       }
1635       ierr = DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr);
1636     }
1637   }
1638   ierr = PetscHashJKDestroy(&ht);CHKERRQ(ierr);
1639   ierr = MatXAIJSetPreallocation(In, 1, dnz, onz, NULL, NULL);CHKERRQ(ierr);
1640   ierr = MatSetOption(In, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
1641   ierr = PetscFree2(dnz,onz);CHKERRQ(ierr);
1642   for (field = 0; field < Nf; ++field) {
1643     PetscObject      obj;
1644     PetscClassId     id;
1645     PetscDualSpace   Q = NULL;
1646     PetscQuadrature  f;
1647     const PetscReal *qpoints, *qweights;
1648     PetscInt         Nc, Np, fpdim, i, d;
1649 
1650     ierr = PetscDSGetDiscretization(prob, field, &obj);CHKERRQ(ierr);
1651     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1652     if (id == PETSCFE_CLASSID) {
1653       PetscFE fe = (PetscFE) obj;
1654 
1655       ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr);
1656       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1657     } else if (id == PETSCFV_CLASSID) {
1658       PetscFV fv = (PetscFV) obj;
1659 
1660       ierr = PetscFVGetDualSpace(fv, &Q);CHKERRQ(ierr);
1661       Nc   = 1;
1662     }
1663     ierr = PetscDualSpaceGetDimension(Q, &fpdim);CHKERRQ(ierr);
1664     /* For each fine grid cell */
1665     for (cell = cStart; cell < cEnd; ++cell) {
1666       PetscInt *findices,   *cindices;
1667       PetscInt  numFIndices, numCIndices;
1668 
1669       ierr = DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr);
1670       ierr = DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
1671       if (numFIndices != fpdim*Nc) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of fine indices %d != %d dual basis vecs", numFIndices, fpdim*Nc);
1672       for (i = 0; i < fpdim; ++i) {
1673         Vec             pointVec;
1674         PetscScalar    *pV;
1675         IS              coarseCellIS;
1676         const PetscInt *coarseCells;
1677         PetscInt        numCoarseCells, cpdim, q, c, j;
1678 
1679         /* Get points from the dual basis functional quadrature */
1680         ierr = PetscDualSpaceGetFunctional(Q, i, &f);CHKERRQ(ierr);
1681         ierr = PetscQuadratureGetData(f, NULL, &Np, &qpoints, &qweights);CHKERRQ(ierr);
1682         ierr = VecCreateSeq(PETSC_COMM_SELF, Np*dim, &pointVec);CHKERRQ(ierr);
1683         ierr = VecSetBlockSize(pointVec, dim);CHKERRQ(ierr);
1684         ierr = VecGetArray(pointVec, &pV);CHKERRQ(ierr);
1685         for (q = 0; q < Np; ++q) {
1686           /* Transform point to real space */
1687           CoordinatesRefToReal(dim, dim, v0, J, &qpoints[q*dim], x);
1688           for (d = 0; d < dim; ++d) pV[q*dim+d] = x[d];
1689         }
1690         ierr = VecRestoreArray(pointVec, &pV);CHKERRQ(ierr);
1691         /* Get set of coarse cells that overlap points (would like to group points by coarse cell) */
1692         ierr = DMLocatePoints(dmc, pointVec, &coarseCellIS);CHKERRQ(ierr);
1693         /* Update preallocation info */
1694         ierr = ISGetLocalSize(coarseCellIS, &numCoarseCells);CHKERRQ(ierr);
1695         ierr = ISGetIndices(coarseCellIS, &coarseCells);CHKERRQ(ierr);
1696         ierr = VecGetArray(pointVec, &pV);CHKERRQ(ierr);
1697         for (ccell = 0; ccell < numCoarseCells; ++ccell) {
1698           PetscReal pVReal[3];
1699 
1700           ierr = DMPlexGetClosureIndices(dmc, csection, globalCSection, coarseCells[ccell], &numCIndices, &cindices, NULL);CHKERRQ(ierr);
1701           /* Transform points from real space to coarse reference space */
1702           ierr = DMPlexComputeCellGeometryFEM(dmc, coarseCells[ccell], NULL, v0c, Jc, invJc, &detJc);CHKERRQ(ierr);
1703           for (d = 0; d < dim; ++d) pVReal[d] = PetscRealPart(pV[ccell*dim+d]);
1704           CoordinatesRealToRef(dim, dim, v0c, invJc, pVReal, x);
1705 
1706           if (id == PETSCFE_CLASSID) {
1707             PetscFE    fe = (PetscFE) obj;
1708             PetscReal *B;
1709 
1710             /* Evaluate coarse basis on contained point */
1711             ierr = PetscFEGetDimension(fe, &cpdim);CHKERRQ(ierr);
1712             ierr = PetscFEGetTabulation(fe, 1, x, &B, NULL, NULL);CHKERRQ(ierr);
1713             /* Get elemMat entries by multiplying by weight */
1714             for (j = 0; j < cpdim; ++j) {
1715               for (c = 0; c < Nc; ++c) elemMat[(c*cpdim + j)*Nc + c] = B[j*Nc + c]*qweights[ccell];
1716             }
1717             ierr = PetscFERestoreTabulation(fe, 1, x, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr);
1718           } else {
1719             cpdim = 1;
1720             for (j = 0; j < cpdim; ++j) {
1721               for (c = 0; c < Nc; ++c) elemMat[(c*cpdim + j)*Nc + c] = 1.0*qweights[ccell];
1722             }
1723           }
1724           /* Update interpolator */
1725           if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(cell, name, Nc, numCIndices, elemMat);CHKERRQ(ierr);}
1726           ierr = MatSetValues(In, Nc, &findices[i*Nc], numCIndices, cindices, elemMat, INSERT_VALUES);CHKERRQ(ierr);
1727           ierr = DMPlexRestoreClosureIndices(dmc, csection, globalCSection, coarseCells[ccell], &numCIndices, &cindices, NULL);CHKERRQ(ierr);
1728         }
1729         ierr = VecRestoreArray(pointVec, &pV);CHKERRQ(ierr);
1730         ierr = ISRestoreIndices(coarseCellIS, &coarseCells);CHKERRQ(ierr);
1731         ierr = ISDestroy(&coarseCellIS);CHKERRQ(ierr);
1732         ierr = VecDestroy(&pointVec);CHKERRQ(ierr);
1733       }
1734       ierr = DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr);
1735     }
1736   }
1737   ierr = PetscFree3(v0,J,invJ);CHKERRQ(ierr);
1738   ierr = PetscFree3(v0c,Jc,invJc);CHKERRQ(ierr);
1739   ierr = PetscFree(elemMat);CHKERRQ(ierr);
1740   ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1741   ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1742   ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1743   PetscFunctionReturn(0);
1744 }
1745 
1746 #undef __FUNCT__
1747 #define __FUNCT__ "DMPlexComputeInjectorFEM"
1748 PetscErrorCode DMPlexComputeInjectorFEM(DM dmc, DM dmf, VecScatter *sc, void *user)
1749 {
1750   PetscDS        prob;
1751   PetscFE       *feRef;
1752   PetscFV       *fvRef;
1753   Vec            fv, cv;
1754   IS             fis, cis;
1755   PetscSection   fsection, fglobalSection, csection, cglobalSection;
1756   PetscInt      *cmap, *cellCIndices, *cellFIndices, *cindices, *findices;
1757   PetscInt       cTotDim, fTotDim = 0, Nf, f, field, cStart, cEnd, cEndInterior, c, dim, d, startC, endC, offsetC, offsetF, m;
1758   PetscErrorCode ierr;
1759 
1760   PetscFunctionBegin;
1761   ierr = PetscLogEventBegin(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1762   ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr);
1763   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
1764   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
1765   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
1766   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
1767   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
1768   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
1769   ierr = DMPlexGetHybridBounds(dmc, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
1770   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
1771   ierr = DMGetDS(dmc, &prob);CHKERRQ(ierr);
1772   ierr = PetscCalloc2(Nf,&feRef,Nf,&fvRef);CHKERRQ(ierr);
1773   for (f = 0; f < Nf; ++f) {
1774     PetscObject  obj;
1775     PetscClassId id;
1776     PetscInt     fNb = 0, Nc = 0;
1777 
1778     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
1779     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1780     if (id == PETSCFE_CLASSID) {
1781       PetscFE fe = (PetscFE) obj;
1782 
1783       ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr);
1784       ierr = PetscFEGetDimension(feRef[f], &fNb);CHKERRQ(ierr);
1785       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1786     } else if (id == PETSCFV_CLASSID) {
1787       PetscFV        fv = (PetscFV) obj;
1788       PetscDualSpace Q;
1789 
1790       ierr = PetscFVRefine(fv, &fvRef[f]);CHKERRQ(ierr);
1791       ierr = PetscFVGetDualSpace(fvRef[f], &Q);CHKERRQ(ierr);
1792       ierr = PetscDualSpaceGetDimension(Q, &fNb);CHKERRQ(ierr);
1793       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
1794     }
1795     fTotDim += fNb*Nc;
1796   }
1797   ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr);
1798   ierr = PetscMalloc1(cTotDim,&cmap);CHKERRQ(ierr);
1799   for (field = 0, offsetC = 0, offsetF = 0; field < Nf; ++field) {
1800     PetscFE        feC;
1801     PetscFV        fvC;
1802     PetscDualSpace QF, QC;
1803     PetscInt       NcF, NcC, fpdim, cpdim;
1804 
1805     if (feRef[field]) {
1806       ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &feC);CHKERRQ(ierr);
1807       ierr = PetscFEGetNumComponents(feC, &NcC);CHKERRQ(ierr);
1808       ierr = PetscFEGetNumComponents(feRef[field], &NcF);CHKERRQ(ierr);
1809       ierr = PetscFEGetDualSpace(feRef[field], &QF);CHKERRQ(ierr);
1810       ierr = PetscDualSpaceGetDimension(QF, &fpdim);CHKERRQ(ierr);
1811       ierr = PetscFEGetDualSpace(feC, &QC);CHKERRQ(ierr);
1812       ierr = PetscDualSpaceGetDimension(QC, &cpdim);CHKERRQ(ierr);
1813     } else {
1814       ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fvC);CHKERRQ(ierr);
1815       ierr = PetscFVGetNumComponents(fvC, &NcC);CHKERRQ(ierr);
1816       ierr = PetscFVGetNumComponents(fvRef[field], &NcF);CHKERRQ(ierr);
1817       ierr = PetscFVGetDualSpace(fvRef[field], &QF);CHKERRQ(ierr);
1818       ierr = PetscDualSpaceGetDimension(QF, &fpdim);CHKERRQ(ierr);
1819       ierr = PetscFVGetDualSpace(fvC, &QC);CHKERRQ(ierr);
1820       ierr = PetscDualSpaceGetDimension(QC, &cpdim);CHKERRQ(ierr);
1821     }
1822     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);
1823     for (c = 0; c < cpdim; ++c) {
1824       PetscQuadrature  cfunc;
1825       const PetscReal *cqpoints;
1826       PetscInt         NpC;
1827       PetscBool        found = PETSC_FALSE;
1828 
1829       ierr = PetscDualSpaceGetFunctional(QC, c, &cfunc);CHKERRQ(ierr);
1830       ierr = PetscQuadratureGetData(cfunc, NULL, &NpC, &cqpoints, NULL);CHKERRQ(ierr);
1831       if (NpC != 1 && feRef[field]) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not know how to do injection for moments");
1832       for (f = 0; f < fpdim; ++f) {
1833         PetscQuadrature  ffunc;
1834         const PetscReal *fqpoints;
1835         PetscReal        sum = 0.0;
1836         PetscInt         NpF, comp;
1837 
1838         ierr = PetscDualSpaceGetFunctional(QF, f, &ffunc);CHKERRQ(ierr);
1839         ierr = PetscQuadratureGetData(ffunc, NULL, &NpF, &fqpoints, NULL);CHKERRQ(ierr);
1840         if (NpC != NpF) continue;
1841         for (d = 0; d < dim; ++d) sum += PetscAbsReal(cqpoints[d] - fqpoints[d]);
1842         if (sum > 1.0e-9) continue;
1843         for (comp = 0; comp < NcC; ++comp) {
1844           cmap[(offsetC+c)*NcC+comp] = (offsetF+f)*NcF+comp;
1845         }
1846         found = PETSC_TRUE;
1847         break;
1848       }
1849       if (!found) {
1850         /* TODO We really want the average here, but some asshole put VecScatter in the interface */
1851         if (fvRef[field]) {
1852           PetscInt comp;
1853           for (comp = 0; comp < NcC; ++comp) {
1854             cmap[(offsetC+c)*NcC+comp] = (offsetF+0)*NcF+comp;
1855           }
1856         } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not locate matching functional for injection");
1857       }
1858     }
1859     offsetC += cpdim*NcC;
1860     offsetF += fpdim*NcF;
1861   }
1862   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);ierr = PetscFVDestroy(&fvRef[f]);CHKERRQ(ierr);}
1863   ierr = PetscFree2(feRef,fvRef);CHKERRQ(ierr);
1864 
1865   ierr = DMGetGlobalVector(dmf, &fv);CHKERRQ(ierr);
1866   ierr = DMGetGlobalVector(dmc, &cv);CHKERRQ(ierr);
1867   ierr = VecGetOwnershipRange(cv, &startC, &endC);CHKERRQ(ierr);
1868   ierr = PetscSectionGetConstrainedStorageSize(cglobalSection, &m);CHKERRQ(ierr);
1869   ierr = PetscMalloc2(cTotDim,&cellCIndices,fTotDim,&cellFIndices);CHKERRQ(ierr);
1870   ierr = PetscMalloc1(m,&cindices);CHKERRQ(ierr);
1871   ierr = PetscMalloc1(m,&findices);CHKERRQ(ierr);
1872   for (d = 0; d < m; ++d) cindices[d] = findices[d] = -1;
1873   for (c = cStart; c < cEnd; ++c) {
1874     ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, c, cellCIndices, cellFIndices);CHKERRQ(ierr);
1875     for (d = 0; d < cTotDim; ++d) {
1876       if ((cellCIndices[d] < startC) || (cellCIndices[d] >= endC)) continue;
1877       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]]);
1878       cindices[cellCIndices[d]-startC] = cellCIndices[d];
1879       findices[cellCIndices[d]-startC] = cellFIndices[cmap[d]];
1880     }
1881   }
1882   ierr = PetscFree(cmap);CHKERRQ(ierr);
1883   ierr = PetscFree2(cellCIndices,cellFIndices);CHKERRQ(ierr);
1884 
1885   ierr = ISCreateGeneral(PETSC_COMM_SELF, m, cindices, PETSC_OWN_POINTER, &cis);CHKERRQ(ierr);
1886   ierr = ISCreateGeneral(PETSC_COMM_SELF, m, findices, PETSC_OWN_POINTER, &fis);CHKERRQ(ierr);
1887   ierr = VecScatterCreate(cv, cis, fv, fis, sc);CHKERRQ(ierr);
1888   ierr = ISDestroy(&cis);CHKERRQ(ierr);
1889   ierr = ISDestroy(&fis);CHKERRQ(ierr);
1890   ierr = DMRestoreGlobalVector(dmf, &fv);CHKERRQ(ierr);
1891   ierr = DMRestoreGlobalVector(dmc, &cv);CHKERRQ(ierr);
1892   ierr = PetscLogEventEnd(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1893   PetscFunctionReturn(0);
1894 }
1895