xref: /petsc/src/dm/impls/plex/plexfem.c (revision 8d3c1932cd76dc2d1ac854cce48898cc6ac885a3)
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, 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, const PetscReal *, PetscInt, PetscScalar *, void *) = DMPlexProjectRigidBody;
158 
159     ctx[0] = dimEmbed;
160     ctx[1] = d;
161     ierr = DMPlexProjectFunction(dm, &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(), DMPlexProjectFunctionLocal(), DMPlexProjectFunctionLabelLocal()
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(), DMPlexProjectFunctionLocal(), DMPlexProjectFunctionLabelLocal()
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__ "DMPlexProjectFunctionLabelLocal"
235 PetscErrorCode DMPlexProjectFunctionLabelLocal(DM dm, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscErrorCode (**funcs)(PetscInt, 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       ierr = ISGetLocalSize(pointIS, &n);CHKERRQ(ierr);
305       ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
306       for (p = 0; p < n; ++p) {
307         const PetscInt    point = points[p];
308         PetscFECellGeom   geom;
309 
310         if ((point < pStart) || (point >= pEnd)) continue;
311         ierr          = DMPlexComputeCellGeometryFEM(dm, point, NULL, geom.v0, geom.J, NULL, &geom.detJ);CHKERRQ(ierr);
312         geom.dim      = dim - h;
313         geom.dimEmbed = dimEmbed;
314         for (f = 0, v = 0; f < numFields; ++f) {
315           void * const ctx = ctxs ? ctxs[f] : NULL;
316 
317           if (!sp[f]) continue;
318           ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
319           for (d = 0; d < spDim; ++d) {
320             if (funcs[f]) {
321               ierr = PetscDualSpaceApply(sp[f], d, &geom, numComp[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr);
322             } else {
323               for (comp = 0; comp < numComp[f]; ++comp) values[v+comp] = 0.0;
324             }
325             v += numComp[f];
326           }
327         }
328         ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, point, values, mode);CHKERRQ(ierr);
329       }
330       ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
331       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
332     }
333     if (h) {
334       for (f = 0; f < numFields; ++f) {ierr = PetscDualSpaceDestroy(&sp[f]);CHKERRQ(ierr);}
335     }
336   }
337   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
338   ierr = DMRestoreWorkArray(dm, numFields, PETSC_BOOL, &fieldActive);CHKERRQ(ierr);
339   for (f = 0; f < numFields; ++f) {ierr = PetscDualSpaceDestroy(&sp[f]);CHKERRQ(ierr);}
340   ierr = PetscFree2(sp, numComp);CHKERRQ(ierr);
341   if (maxHeight > 0) {
342     ierr = PetscFree(cellsp);CHKERRQ(ierr);
343   }
344   PetscFunctionReturn(0);
345 }
346 
347 #undef __FUNCT__
348 #define __FUNCT__ "DMPlexProjectFunctionLocal"
349 PetscErrorCode DMPlexProjectFunctionLocal(DM dm, PetscErrorCode (**funcs)(PetscInt, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
350 {
351   PetscDualSpace *sp, *cellsp;
352   PetscInt       *numComp;
353   PetscSection    section;
354   PetscScalar    *values;
355   PetscInt        Nf, dim, dimEmbed, spDim, totDim = 0, numValues, pStart, pEnd, p, cStart, cEnd, cEndInterior, f, d, v, comp, h, maxHeight;
356   PetscErrorCode  ierr;
357 
358   PetscFunctionBegin;
359   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
360   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
361   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
362   if (cEnd <= cStart) PetscFunctionReturn(0);
363   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
364   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
365   ierr = PetscMalloc2(Nf, &sp, Nf, &numComp);CHKERRQ(ierr);
366   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
367   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
368   ierr = DMPlexGetMaxProjectionHeight(dm,&maxHeight);CHKERRQ(ierr);
369   if (maxHeight < 0 || maxHeight > dim) {SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"maximum projection height %d not in [0, %d)\n", maxHeight,dim);}
370   if (maxHeight > 0) {
371     ierr = PetscMalloc1(Nf,&cellsp);CHKERRQ(ierr);
372   }
373   else {
374     cellsp = sp;
375   }
376   for (h = 0; h <= maxHeight; h++) {
377     ierr = DMPlexGetHeightStratum(dm, h, &pStart, &pEnd);CHKERRQ(ierr);
378     if (!h) {pStart = cStart; pEnd = cEnd;}
379     if (pEnd <= pStart) continue;
380     totDim = 0;
381     for (f = 0; f < Nf; ++f) {
382       PetscObject  obj;
383       PetscClassId id;
384 
385       ierr = DMGetField(dm, f, &obj);CHKERRQ(ierr);
386       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
387       if (id == PETSCFE_CLASSID) {
388         PetscFE fe = (PetscFE) obj;
389 
390         ierr = PetscFEGetNumComponents(fe, &numComp[f]);CHKERRQ(ierr);
391         if (!h) {
392           ierr  = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr);
393           sp[f] = cellsp[f];
394           ierr  = PetscObjectReference((PetscObject) sp[f]);CHKERRQ(ierr);
395         }
396         else {
397           ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], h, &sp[f]);CHKERRQ(ierr);
398           if (!sp[f]) {
399             continue;
400           }
401         }
402       } else if (id == PETSCFV_CLASSID) {
403         PetscFV fv = (PetscFV) obj;
404 
405         if (!h) {
406           ierr = PetscFVGetNumComponents(fv, &numComp[f]);CHKERRQ(ierr);
407           ierr = PetscFVGetDualSpace(fv, &cellsp[f]);CHKERRQ(ierr);
408           ierr = PetscObjectReference((PetscObject) cellsp[f]);CHKERRQ(ierr);
409           sp[f] = cellsp[f];
410         }
411         else {
412           sp[f] = NULL;
413           continue; /* FV doesn't have fields that can be evaluated at faces, corners, etc. */
414         }
415       } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
416       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
417       totDim += spDim*numComp[f];
418     }
419     ierr = DMPlexVecGetClosure(dm, section, localX, pStart, &numValues, NULL);CHKERRQ(ierr);
420     if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section point closure size %d != dual space dimension %d", numValues, totDim);
421     if (!totDim) continue;
422     ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
423     for (p = pStart; p < pEnd; ++p) {
424       PetscFECellGeom geom;
425 
426       ierr          = DMPlexComputeCellGeometryFEM(dm, p, NULL, geom.v0, geom.J, NULL, &geom.detJ);CHKERRQ(ierr);
427       geom.dim      = dim - h;
428       geom.dimEmbed = dimEmbed;
429       for (f = 0, v = 0; f < Nf; ++f) {
430         void * const ctx = ctxs ? ctxs[f] : NULL;
431 
432         if (!sp[f]) continue;
433         ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
434         for (d = 0; d < spDim; ++d) {
435           if (funcs[f]) {
436             ierr = PetscDualSpaceApply(sp[f], d, &geom, numComp[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr);
437           } else {
438             for (comp = 0; comp < numComp[f]; ++comp) values[v+comp] = 0.0;
439           }
440           v += numComp[f];
441         }
442       }
443       ierr = DMPlexVecSetClosure(dm, section, localX, p, values, mode);CHKERRQ(ierr);
444     }
445     ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
446     if (h || !maxHeight) {
447       for (f = 0; f < Nf; f++) {ierr = PetscDualSpaceDestroy(&sp[f]);CHKERRQ(ierr);}
448     }
449   }
450   ierr = PetscFree2(sp, numComp);CHKERRQ(ierr);
451   if (maxHeight > 0) {
452     for (f = 0; f < Nf; f++) {ierr = PetscDualSpaceDestroy(&cellsp[f]);CHKERRQ(ierr);}
453     ierr = PetscFree(cellsp);CHKERRQ(ierr);
454   }
455   PetscFunctionReturn(0);
456 }
457 
458 #undef __FUNCT__
459 #define __FUNCT__ "DMPlexProjectFunction"
460 /*@C
461   DMPlexProjectFunction - This projects the given function into the function space provided.
462 
463   Input Parameters:
464 + dm      - The DM
465 . funcs   - The coordinate functions to evaluate, one per field
466 . ctxs    - Optional array of contexts to pass to each coordinate function.  ctxs itself may be null.
467 - mode    - The insertion mode for values
468 
469   Output Parameter:
470 . X - vector
471 
472    Calling sequence of func:
473 $    func(PetscInt dim, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx);
474 
475 +  dim - The spatial dimension
476 .  x   - The coordinates
477 .  Nf  - The number of fields
478 .  u   - The output field values
479 -  ctx - optional user-defined function context
480 
481   Level: developer
482 
483 .seealso: DMPlexComputeL2Diff()
484 @*/
485 PetscErrorCode DMPlexProjectFunction(DM dm, PetscErrorCode (**funcs)(PetscInt, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec X)
486 {
487   Vec            localX;
488   PetscErrorCode ierr;
489 
490   PetscFunctionBegin;
491   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
492   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
493   ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, mode, localX);CHKERRQ(ierr);
494   ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr);
495   ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr);
496   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
497   PetscFunctionReturn(0);
498 }
499 
500 #undef __FUNCT__
501 #define __FUNCT__ "DMPlexProjectFieldLocal"
502 PetscErrorCode DMPlexProjectFieldLocal(DM dm, Vec localU, void (**funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal [], PetscScalar []), InsertMode mode, Vec localX)
503 {
504   DM                dmAux;
505   PetscDS           prob, probAux = NULL;
506   Vec               A;
507   PetscSection      section, sectionAux = NULL;
508   PetscDualSpace   *sp;
509   PetscInt         *Ncf;
510   PetscScalar      *values, *u, *u_x, *a, *a_x;
511   PetscReal        *x, *v0, *J, *invJ, detJ;
512   PetscInt          Nf, dim, spDim, totDim, numValues, cStart, cEnd, cEndInterior, c, f, d, v, comp, maxHeight;
513   PetscErrorCode    ierr;
514 
515   PetscFunctionBegin;
516   ierr = DMPlexGetMaxProjectionHeight(dm,&maxHeight);CHKERRQ(ierr);
517   if (maxHeight > 0) {SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Field projection for height > 0 not supported yet");}
518   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
519   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
520   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
521   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
522   ierr = PetscMalloc2(Nf, &sp, Nf, &Ncf);CHKERRQ(ierr);
523   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
524   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
525   ierr = PetscDSGetEvaluationArrays(prob, &u, NULL, &u_x);CHKERRQ(ierr);
526   ierr = PetscDSGetRefCoordArrays(prob, &x, NULL);CHKERRQ(ierr);
527   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
528   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
529   if (dmAux) {
530     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
531     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
532     ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);CHKERRQ(ierr);
533   }
534   ierr = DMPlexInsertBoundaryValues(dm, localU, 0.0, NULL, NULL, NULL);CHKERRQ(ierr);
535   ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr);
536   if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section cell closure size %d != dual space dimension %d", numValues, totDim);
537   ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
538   ierr = PetscMalloc3(dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
539   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
540   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
541   for (c = cStart; c < cEnd; ++c) {
542     PetscScalar *coefficients = NULL, *coefficientsAux = NULL;
543 
544     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
545     ierr = DMPlexVecGetClosure(dm, section, localU, c, NULL, &coefficients);CHKERRQ(ierr);
546     if (dmAux) {ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &coefficientsAux);CHKERRQ(ierr);}
547     for (f = 0, v = 0; f < Nf; ++f) {
548       PetscObject  obj;
549       PetscClassId id;
550 
551       ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
552       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
553       if (id == PETSCFE_CLASSID) {
554         PetscFE fe = (PetscFE) obj;
555 
556         ierr = PetscFEGetDualSpace(fe, &sp[f]);CHKERRQ(ierr);
557         ierr  = PetscObjectReference((PetscObject) sp[f]);CHKERRQ(ierr);
558         ierr = PetscFEGetNumComponents(fe, &Ncf[f]);CHKERRQ(ierr);
559       } else if (id == PETSCFV_CLASSID) {
560         PetscFV fv = (PetscFV) obj;
561 
562         ierr = PetscFVGetNumComponents(fv, &Ncf[f]);CHKERRQ(ierr);
563         ierr = PetscFVGetDualSpace(fv, &sp[f]);CHKERRQ(ierr);
564         ierr = PetscObjectReference((PetscObject) sp[f]);CHKERRQ(ierr);
565       }
566       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
567       for (d = 0; d < spDim; ++d) {
568         PetscQuadrature  quad;
569         const PetscReal *points, *weights;
570         PetscInt         numPoints, q;
571 
572         if (funcs[f]) {
573           ierr = PetscDualSpaceGetFunctional(sp[f], d, &quad);CHKERRQ(ierr);
574           ierr = PetscQuadratureGetData(quad, NULL, &numPoints, &points, &weights);CHKERRQ(ierr);
575           for (q = 0; q < numPoints; ++q) {
576             CoordinatesRefToReal(dim, dim, v0, J, &points[q*dim], x);
577             ierr = EvaluateFieldJets(prob,    PETSC_FALSE, q, invJ, coefficients,    NULL, u, u_x, NULL);CHKERRQ(ierr);
578             ierr = EvaluateFieldJets(probAux, PETSC_FALSE, q, invJ, coefficientsAux, NULL, a, a_x, NULL);CHKERRQ(ierr);
579             (*funcs[f])(u, NULL, u_x, a, NULL, a_x, x, &values[v]);
580           }
581         } else {
582           for (comp = 0; comp < Ncf[f]; ++comp) values[v+comp] = 0.0;
583         }
584         v += Ncf[f];
585       }
586     }
587     ierr = DMPlexVecRestoreClosure(dm, section, localU, c, NULL, &coefficients);CHKERRQ(ierr);
588     if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &coefficientsAux);CHKERRQ(ierr);}
589     ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr);
590   }
591   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
592   ierr = PetscFree3(v0,J,invJ);CHKERRQ(ierr);
593   for (f = 0; f < Nf; f++) {ierr = PetscDualSpaceDestroy(&sp[f]);CHKERRQ(ierr);}
594   ierr = PetscFree2(sp, Ncf);CHKERRQ(ierr);
595   PetscFunctionReturn(0);
596 }
597 
598 #undef __FUNCT__
599 #define __FUNCT__ "DMPlexInsertBoundaryValues_FEM_Internal"
600 static PetscErrorCode DMPlexInsertBoundaryValues_FEM_Internal(DM dm, PetscInt field, DMLabel label, PetscInt numids, const PetscInt ids[], PetscErrorCode (*func)(PetscInt, const PetscReal[], PetscInt, PetscScalar *, void *), void *ctx, Vec locX)
601 {
602   PetscErrorCode (**funcs)(PetscInt, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
603   void            **ctxs;
604   PetscInt          numFields;
605   PetscErrorCode    ierr;
606 
607   PetscFunctionBegin;
608   ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr);
609   ierr = PetscCalloc2(numFields,&funcs,numFields,&ctxs);CHKERRQ(ierr);
610   funcs[field] = func;
611   ctxs[field]  = ctx;
612   ierr = DMPlexProjectFunctionLabelLocal(dm, label, numids, ids, funcs, ctxs, INSERT_BC_VALUES, locX);CHKERRQ(ierr);
613   ierr = PetscFree2(funcs,ctxs);CHKERRQ(ierr);
614   PetscFunctionReturn(0);
615 }
616 
617 #undef __FUNCT__
618 #define __FUNCT__ "DMPlexInsertBoundaryValues_FVM_Internal"
619 /* This ignores numcomps/comps */
620 static PetscErrorCode DMPlexInsertBoundaryValues_FVM_Internal(DM dm, PetscReal time, Vec faceGeometry, Vec cellGeometry, Vec Grad,
621                                                               PetscInt field, DMLabel label, PetscInt numids, const PetscInt ids[], PetscErrorCode (*func)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*), void *ctx, Vec locX)
622 {
623   PetscDS            prob;
624   PetscSF            sf;
625   DM                 dmFace, dmCell, dmGrad;
626   const PetscScalar *facegeom, *cellgeom, *grad;
627   const PetscInt    *leaves;
628   PetscScalar       *x, *fx;
629   PetscInt           dim, nleaves, loc, fStart, fEnd, pdim, i;
630   PetscErrorCode     ierr;
631 
632   PetscFunctionBegin;
633   ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
634   ierr = PetscSFGetGraph(sf, NULL, &nleaves, &leaves, NULL);CHKERRQ(ierr);
635   nleaves = PetscMax(0, nleaves);
636   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
637   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
638   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
639   ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr);
640   ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
641   ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr);
642   ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
643   if (Grad) {
644     PetscFV fv;
645 
646     ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fv);CHKERRQ(ierr);
647     ierr = VecGetDM(Grad, &dmGrad);CHKERRQ(ierr);
648     ierr = VecGetArrayRead(Grad, &grad);CHKERRQ(ierr);
649     ierr = PetscFVGetNumComponents(fv, &pdim);CHKERRQ(ierr);
650     ierr = DMGetWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr);
651   }
652   ierr = VecGetArray(locX, &x);CHKERRQ(ierr);
653   for (i = 0; i < numids; ++i) {
654     IS              faceIS;
655     const PetscInt *faces;
656     PetscInt        numFaces, f;
657 
658     ierr = DMLabelGetStratumIS(label, ids[i], &faceIS);CHKERRQ(ierr);
659     if (!faceIS) continue; /* No points with that id on this process */
660     ierr = ISGetLocalSize(faceIS, &numFaces);CHKERRQ(ierr);
661     ierr = ISGetIndices(faceIS, &faces);CHKERRQ(ierr);
662     for (f = 0; f < numFaces; ++f) {
663       const PetscInt         face = faces[f], *cells;
664       const PetscFVFaceGeom *fg;
665 
666       if ((face < fStart) || (face >= fEnd)) continue; /* Refinement adds non-faces to labels */
667       ierr = PetscFindInt(face, nleaves, (PetscInt *) leaves, &loc);CHKERRQ(ierr);
668       if (loc >= 0) continue;
669       ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr);
670       ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr);
671       if (Grad) {
672         const PetscFVCellGeom *cg;
673         const PetscScalar     *cx, *cgrad;
674         PetscScalar           *xG;
675         PetscReal              dx[3];
676         PetscInt               d;
677 
678         ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cg);CHKERRQ(ierr);
679         ierr = DMPlexPointLocalRead(dm, cells[0], x, &cx);CHKERRQ(ierr);
680         ierr = DMPlexPointLocalRead(dmGrad, cells[0], grad, &cgrad);CHKERRQ(ierr);
681         ierr = DMPlexPointLocalFieldRef(dm, cells[1], field, x, &xG);CHKERRQ(ierr);
682         DMPlex_WaxpyD_Internal(dim, -1, cg->centroid, fg->centroid, dx);
683         for (d = 0; d < pdim; ++d) fx[d] = cx[d] + DMPlex_DotD_Internal(dim, &cgrad[d*dim], dx);
684         ierr = (*func)(time, fg->centroid, fg->normal, fx, xG, ctx);CHKERRQ(ierr);
685       } else {
686         const PetscScalar *xI;
687         PetscScalar       *xG;
688 
689         ierr = DMPlexPointLocalRead(dm, cells[0], x, &xI);CHKERRQ(ierr);
690         ierr = DMPlexPointLocalFieldRef(dm, cells[1], field, x, &xG);CHKERRQ(ierr);
691         ierr = (*func)(time, fg->centroid, fg->normal, xI, xG, ctx);CHKERRQ(ierr);
692       }
693     }
694     ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr);
695     ierr = ISDestroy(&faceIS);CHKERRQ(ierr);
696   }
697   ierr = VecRestoreArray(locX, &x);CHKERRQ(ierr);
698   if (Grad) {
699     ierr = DMRestoreWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr);
700     ierr = VecRestoreArrayRead(Grad, &grad);CHKERRQ(ierr);
701   }
702   ierr = VecRestoreArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
703   ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
704   PetscFunctionReturn(0);
705 }
706 
707 #undef __FUNCT__
708 #define __FUNCT__ "DMPlexInsertBoundaryValues"
709 PetscErrorCode DMPlexInsertBoundaryValues(DM dm, Vec locX, PetscReal time, Vec faceGeomFVM, Vec cellGeomFVM, Vec gradFVM)
710 {
711   PetscInt       numBd, b;
712   PetscErrorCode ierr;
713 
714   PetscFunctionBegin;
715   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
716   PetscValidHeaderSpecific(locX, VEC_CLASSID, 2);
717   if (faceGeomFVM) {PetscValidHeaderSpecific(faceGeomFVM, VEC_CLASSID, 4);}
718   if (cellGeomFVM) {PetscValidHeaderSpecific(cellGeomFVM, VEC_CLASSID, 5);}
719   if (gradFVM)     {PetscValidHeaderSpecific(gradFVM, VEC_CLASSID, 6);}
720   ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
721   for (b = 0; b < numBd; ++b) {
722     PetscBool       isEssential;
723     const char     *labelname;
724     DMLabel         label;
725     PetscInt        field;
726     PetscObject     obj;
727     PetscClassId    id;
728     void          (*func)();
729     PetscInt        numids;
730     const PetscInt *ids;
731     void           *ctx;
732 
733     ierr = DMPlexGetBoundary(dm, b, &isEssential, NULL, &labelname, &field, NULL, NULL, &func, &numids, &ids, &ctx);CHKERRQ(ierr);
734     if (!isEssential) continue;
735     ierr = DMPlexGetLabel(dm, labelname, &label);CHKERRQ(ierr);
736     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
737     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
738     if (id == PETSCFE_CLASSID) {
739       ierr = DMPlexInsertBoundaryValues_FEM_Internal(dm, field, label, numids, ids, (PetscErrorCode (*)(PetscInt, const PetscReal[], PetscInt, PetscScalar *, void *)) func, ctx, locX);CHKERRQ(ierr);
740     } else if (id == PETSCFV_CLASSID) {
741       if (!faceGeomFVM) continue;
742       ierr = DMPlexInsertBoundaryValues_FVM_Internal(dm, time, faceGeomFVM, cellGeomFVM, gradFVM,
743                                                      field, label, numids, ids, (PetscErrorCode (*)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*)) func, ctx, locX);CHKERRQ(ierr);
744     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
745   }
746   PetscFunctionReturn(0);
747 }
748 
749 #undef __FUNCT__
750 #define __FUNCT__ "DMPlexComputeL2Diff"
751 /*@C
752   DMPlexComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.
753 
754   Input Parameters:
755 + dm    - The DM
756 . funcs - The functions to evaluate for each field component
757 . ctxs  - Optional array of contexts to pass to each function, or NULL.
758 - X     - The coefficient vector u_h
759 
760   Output Parameter:
761 . diff - The diff ||u - u_h||_2
762 
763   Level: developer
764 
765 .seealso: DMPlexProjectFunction(), DMPlexComputeL2FieldDiff(), DMPlexComputeL2GradientDiff()
766 @*/
767 PetscErrorCode DMPlexComputeL2Diff(DM dm, PetscErrorCode (**funcs)(PetscInt, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
768 {
769   const PetscInt   debug = 0;
770   PetscSection     section;
771   PetscQuadrature  quad;
772   Vec              localX;
773   PetscScalar     *funcVal, *interpolant;
774   PetscReal       *coords, *v0, *J, *invJ, detJ;
775   PetscReal        localDiff = 0.0;
776   const PetscReal *quadPoints, *quadWeights;
777   PetscInt         dim, numFields, numComponents = 0, numQuadPoints, cStart, cEnd, cEndInterior, c, field, fieldOffset;
778   PetscErrorCode   ierr;
779 
780   PetscFunctionBegin;
781   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
782   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
783   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
784   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
785   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
786   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
787   for (field = 0; field < numFields; ++field) {
788     PetscObject  obj;
789     PetscClassId id;
790     PetscInt     Nc;
791 
792     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
793     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
794     if (id == PETSCFE_CLASSID) {
795       PetscFE fe = (PetscFE) obj;
796 
797       ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
798       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
799     } else if (id == PETSCFV_CLASSID) {
800       PetscFV fv = (PetscFV) obj;
801 
802       ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr);
803       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
804     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
805     numComponents += Nc;
806   }
807   ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
808   ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
809   ierr = PetscMalloc6(numComponents,&funcVal,numComponents,&interpolant,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
810   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
811   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
812   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
813   for (c = cStart; c < cEnd; ++c) {
814     PetscScalar *x = NULL;
815     PetscReal    elemDiff = 0.0;
816 
817     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
818     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
819     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
820 
821     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
822       PetscObject  obj;
823       PetscClassId id;
824       void * const ctx = ctxs ? ctxs[field] : NULL;
825       PetscInt     Nb, Nc, q, fc;
826 
827       ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
828       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
829       if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);}
830       else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;}
831       else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
832       if (debug) {
833         char title[1024];
834         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
835         ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr);
836       }
837       for (q = 0; q < numQuadPoints; ++q) {
838         CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], coords);
839         ierr = (*funcs[field])(dim, coords, Nc, funcVal, ctx);CHKERRQ(ierr);
840         if (id == PETSCFE_CLASSID)      {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
841         else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
842         else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
843         for (fc = 0; fc < Nc; ++fc) {
844           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);}
845           elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ;
846         }
847       }
848       fieldOffset += Nb*Nc;
849     }
850     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
851     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
852     localDiff += elemDiff;
853   }
854   ierr  = PetscFree6(funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr);
855   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
856   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
857   *diff = PetscSqrtReal(*diff);
858   PetscFunctionReturn(0);
859 }
860 
861 #undef __FUNCT__
862 #define __FUNCT__ "DMPlexComputeL2GradientDiff"
863 /*@C
864   DMPlexComputeL2GradientDiff - This function computes the L_2 difference between the gradient of a function u and an FEM interpolant solution grad u_h.
865 
866   Input Parameters:
867 + dm    - The DM
868 . funcs - The gradient functions to evaluate for each field component
869 . ctxs  - Optional array of contexts to pass to each function, or NULL.
870 . X     - The coefficient vector u_h
871 - n     - The vector to project along
872 
873   Output Parameter:
874 . diff - The diff ||(grad u - grad u_h) . n||_2
875 
876   Level: developer
877 
878 .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff()
879 @*/
880 PetscErrorCode DMPlexComputeL2GradientDiff(DM dm, PetscErrorCode (**funcs)(PetscInt, const PetscReal [], const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, const PetscReal n[], PetscReal *diff)
881 {
882   const PetscInt  debug = 0;
883   PetscSection    section;
884   PetscQuadrature quad;
885   Vec             localX;
886   PetscScalar    *funcVal, *interpolantVec;
887   PetscReal      *coords, *realSpaceDer, *v0, *J, *invJ, detJ;
888   PetscReal       localDiff = 0.0;
889   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, cEndInterior, c, field, fieldOffset, comp;
890   PetscErrorCode  ierr;
891 
892   PetscFunctionBegin;
893   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
894   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
895   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
896   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
897   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
898   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
899   for (field = 0; field < numFields; ++field) {
900     PetscFE  fe;
901     PetscInt Nc;
902 
903     ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
904     ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
905     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
906     numComponents += Nc;
907   }
908   /* ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */
909   ierr = PetscMalloc7(numComponents,&funcVal,dim,&coords,dim,&realSpaceDer,dim,&v0,dim*dim,&J,dim*dim,&invJ,dim,&interpolantVec);CHKERRQ(ierr);
910   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
911   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
912   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
913   for (c = cStart; c < cEnd; ++c) {
914     PetscScalar *x = NULL;
915     PetscReal    elemDiff = 0.0;
916 
917     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
918     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
919     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
920 
921     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
922       PetscFE          fe;
923       void * const     ctx = ctxs ? ctxs[field] : NULL;
924       const PetscReal *quadPoints, *quadWeights;
925       PetscReal       *basisDer;
926       PetscInt         numQuadPoints, Nb, Ncomp, q, d, e, fc, f, g;
927 
928       ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
929       ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
930       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
931       ierr = PetscFEGetNumComponents(fe, &Ncomp);CHKERRQ(ierr);
932       ierr = PetscFEGetDefaultTabulation(fe, NULL, &basisDer, NULL);CHKERRQ(ierr);
933       if (debug) {
934         char title[1024];
935         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
936         ierr = DMPrintCellVector(c, title, Nb*Ncomp, &x[fieldOffset]);CHKERRQ(ierr);
937       }
938       for (q = 0; q < numQuadPoints; ++q) {
939         for (d = 0; d < dim; d++) {
940           coords[d] = v0[d];
941           for (e = 0; e < dim; e++) {
942             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
943           }
944         }
945         ierr = (*funcs[field])(dim, coords, n, numFields, funcVal, ctx);CHKERRQ(ierr);
946         for (fc = 0; fc < Ncomp; ++fc) {
947           PetscScalar interpolant = 0.0;
948 
949           for (d = 0; d < dim; ++d) interpolantVec[d] = 0.0;
950           for (f = 0; f < Nb; ++f) {
951             const PetscInt fidx = f*Ncomp+fc;
952 
953             for (d = 0; d < dim; ++d) {
954               realSpaceDer[d] = 0.0;
955               for (g = 0; g < dim; ++g) {
956                 realSpaceDer[d] += invJ[g*dim+d]*basisDer[(q*Nb*Ncomp+fidx)*dim+g];
957               }
958               interpolantVec[d] += x[fieldOffset+fidx]*realSpaceDer[d];
959             }
960           }
961           for (d = 0; d < dim; ++d) interpolant += interpolantVec[d]*n[d];
962           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);}
963           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ;
964         }
965       }
966       comp        += Ncomp;
967       fieldOffset += Nb*Ncomp;
968     }
969     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
970     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
971     localDiff += elemDiff;
972   }
973   ierr  = PetscFree7(funcVal,coords,realSpaceDer,v0,J,invJ,interpolantVec);CHKERRQ(ierr);
974   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
975   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
976   *diff = PetscSqrtReal(*diff);
977   PetscFunctionReturn(0);
978 }
979 
980 #undef __FUNCT__
981 #define __FUNCT__ "DMPlexComputeL2FieldDiff"
982 /*@C
983   DMPlexComputeL2FieldDiff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h, separated into field components.
984 
985   Input Parameters:
986 + dm    - The DM
987 . funcs - The functions to evaluate for each field component
988 . ctxs  - Optional array of contexts to pass to each function, or NULL.
989 - X     - The coefficient vector u_h
990 
991   Output Parameter:
992 . diff - The array of differences, ||u^f - u^f_h||_2
993 
994   Level: developer
995 
996 .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff(), DMPlexComputeL2GradientDiff()
997 @*/
998 PetscErrorCode DMPlexComputeL2FieldDiff(DM dm, PetscErrorCode (**funcs)(PetscInt, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, PetscReal diff[])
999 {
1000   const PetscInt   debug = 0;
1001   PetscSection     section;
1002   PetscQuadrature  quad;
1003   Vec              localX;
1004   PetscScalar     *funcVal, *interpolant;
1005   PetscReal       *coords, *v0, *J, *invJ, detJ;
1006   PetscReal       *localDiff;
1007   const PetscReal *quadPoints, *quadWeights;
1008   PetscInt         dim, numFields, numComponents = 0, numQuadPoints, cStart, cEnd, cEndInterior, c, field, fieldOffset;
1009   PetscErrorCode   ierr;
1010 
1011   PetscFunctionBegin;
1012   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1013   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1014   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
1015   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
1016   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
1017   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
1018   for (field = 0; field < numFields; ++field) {
1019     PetscObject  obj;
1020     PetscClassId id;
1021     PetscInt     Nc;
1022 
1023     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
1024     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1025     if (id == PETSCFE_CLASSID) {
1026       PetscFE fe = (PetscFE) obj;
1027 
1028       ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
1029       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1030     } else if (id == PETSCFV_CLASSID) {
1031       PetscFV fv = (PetscFV) obj;
1032 
1033       ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr);
1034       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
1035     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
1036     numComponents += Nc;
1037   }
1038   ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
1039   ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
1040   ierr = PetscCalloc7(numFields,&localDiff,numComponents,&funcVal,numComponents,&interpolant,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
1041   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1042   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
1043   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
1044   for (c = cStart; c < cEnd; ++c) {
1045     PetscScalar *x = NULL;
1046 
1047     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
1048     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
1049     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
1050 
1051     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
1052       PetscObject  obj;
1053       PetscClassId id;
1054       void * const ctx = ctxs ? ctxs[field] : NULL;
1055       PetscInt     Nb, Nc, q, fc;
1056 
1057       PetscReal       elemDiff = 0.0;
1058 
1059       ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
1060       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1061       if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);}
1062       else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;}
1063       else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
1064       if (debug) {
1065         char title[1024];
1066         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
1067         ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr);
1068       }
1069       for (q = 0; q < numQuadPoints; ++q) {
1070         CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], coords);
1071         ierr = (*funcs[field])(dim, coords, numFields, funcVal, ctx);CHKERRQ(ierr);
1072         if (id == PETSCFE_CLASSID)      {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
1073         else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
1074         else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
1075         for (fc = 0; fc < Nc; ++fc) {
1076           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);}
1077           elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ;
1078         }
1079       }
1080       fieldOffset += Nb*Nc;
1081       localDiff[field] += elemDiff;
1082     }
1083     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
1084   }
1085   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
1086   ierr = MPI_Allreduce(localDiff, diff, numFields, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
1087   for (field = 0; field < numFields; ++field) diff[field] = PetscSqrtReal(diff[field]);
1088   ierr = PetscFree7(localDiff,funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr);
1089   PetscFunctionReturn(0);
1090 }
1091 
1092 #undef __FUNCT__
1093 #define __FUNCT__ "DMPlexComputeIntegralFEM"
1094 /*@
1095   DMPlexComputeIntegralFEM - Form the local integral F from the local input X using pointwise functions specified by the user
1096 
1097   Input Parameters:
1098 + dm - The mesh
1099 . X  - Local input vector
1100 - user - The user context
1101 
1102   Output Parameter:
1103 . integral - Local integral for each field
1104 
1105   Level: developer
1106 
1107 .seealso: DMPlexComputeResidualFEM()
1108 @*/
1109 PetscErrorCode DMPlexComputeIntegralFEM(DM dm, Vec X, PetscReal *integral, void *user)
1110 {
1111   DM_Plex          *mesh  = (DM_Plex *) dm->data;
1112   DM                dmAux;
1113   Vec               localX, A;
1114   PetscDS           prob, probAux = NULL;
1115   PetscSection      section, sectionAux;
1116   PetscFECellGeom  *cgeom;
1117   PetscScalar      *u, *a = NULL;
1118   PetscReal        *lintegral, *vol;
1119   PetscInt          dim, Nf, f, numCells, cStart, cEnd, cEndInterior, c;
1120   PetscInt          totDim, totDimAux;
1121   PetscErrorCode    ierr;
1122 
1123   PetscFunctionBegin;
1124   ierr = PetscLogEventBegin(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);
1125   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
1126   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
1127   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
1128   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1129   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1130   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
1131   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
1132   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
1133   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1134   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
1135   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
1136   numCells = cEnd - cStart;
1137   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
1138   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
1139   if (dmAux) {
1140     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
1141     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
1142     ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
1143   }
1144   ierr = DMPlexInsertBoundaryValues(dm, localX, 0.0, NULL, NULL, NULL);CHKERRQ(ierr);
1145   ierr = PetscMalloc4(Nf,&lintegral,numCells*totDim,&u,numCells,&cgeom,numCells,&vol);CHKERRQ(ierr);
1146   if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);}
1147   for (f = 0; f < Nf; ++f) {lintegral[f] = 0.0;}
1148   for (c = cStart; c < cEnd; ++c) {
1149     PetscScalar *x = NULL;
1150     PetscInt     i;
1151 
1152     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, cgeom[c].v0, cgeom[c].J, cgeom[c].invJ, &cgeom[c].detJ);CHKERRQ(ierr);
1153     ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol[c], NULL, NULL);CHKERRQ(ierr);
1154     if (cgeom[c].detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", cgeom[c].detJ, c);
1155     ierr = DMPlexVecGetClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr);
1156     for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i];
1157     ierr = DMPlexVecRestoreClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr);
1158     if (dmAux) {
1159       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
1160       for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i];
1161       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
1162     }
1163   }
1164   for (f = 0; f < Nf; ++f) {
1165     PetscObject  obj;
1166     PetscClassId id;
1167     PetscInt     numChunks, numBatches, batchSize, numBlocks, blockSize, Ne, Nr, offset;
1168 
1169     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
1170     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1171     if (id == PETSCFE_CLASSID) {
1172       PetscFE         fe = (PetscFE) obj;
1173       PetscQuadrature q;
1174       PetscInt        Nq, Nb;
1175 
1176       ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
1177       ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
1178       ierr = PetscQuadratureGetData(q, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);
1179       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
1180       blockSize = Nb*Nq;
1181       batchSize = numBlocks * blockSize;
1182       ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1183       numChunks = numCells / (numBatches*batchSize);
1184       Ne        = numChunks*numBatches*batchSize;
1185       Nr        = numCells % (numBatches*batchSize);
1186       offset    = numCells - Nr;
1187       ierr = PetscFEIntegrate(fe, prob, f, Ne, cgeom, u, probAux, a, lintegral);CHKERRQ(ierr);
1188       ierr = PetscFEIntegrate(fe, prob, f, Nr, &cgeom[offset], &u[offset*totDim], probAux, &a[offset*totDimAux], lintegral);CHKERRQ(ierr);
1189     } else if (id == PETSCFV_CLASSID) {
1190       /* PetscFV  fv = (PetscFV) obj; */
1191       PetscInt       foff;
1192       PetscPointFunc obj_func;
1193       PetscScalar    lint;
1194 
1195       ierr = PetscDSGetObjective(prob, f, &obj_func);CHKERRQ(ierr);
1196       ierr = PetscDSGetFieldOffset(prob, f, &foff);CHKERRQ(ierr);
1197       if (obj_func) {
1198         for (c = 0; c < numCells; ++c) {
1199           /* TODO: Need full pointwise interpolation and get centroid for x argument */
1200           obj_func(dim, Nf, 0, NULL, NULL, &u[totDim*c+foff], NULL, NULL, NULL, NULL, NULL, NULL, NULL, 0.0, NULL, &lint);
1201           lintegral[f] = PetscRealPart(lint)*vol[c];
1202         }
1203       }
1204     } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
1205   }
1206   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
1207   if (mesh->printFEM) {
1208     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Local integral:");CHKERRQ(ierr);
1209     for (f = 0; f < Nf; ++f) {ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), " %g", lintegral[f]);CHKERRQ(ierr);}
1210     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "\n");CHKERRQ(ierr);
1211   }
1212   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
1213   ierr = MPI_Allreduce(lintegral, integral, Nf, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject) dm));CHKERRQ(ierr);
1214   ierr = PetscFree4(lintegral,u,cgeom,vol);CHKERRQ(ierr);
1215   ierr = PetscLogEventEnd(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);
1216   PetscFunctionReturn(0);
1217 }
1218 
1219 #undef __FUNCT__
1220 #define __FUNCT__ "DMPlexComputeInterpolatorFEM"
1221 /*@
1222   DMPlexComputeInterpolatorFEM - Form the local portion of the interpolation matrix I from the coarse DM to the uniformly refined DM.
1223 
1224   Input Parameters:
1225 + dmf  - The fine mesh
1226 . dmc  - The coarse mesh
1227 - user - The user context
1228 
1229   Output Parameter:
1230 . In  - The interpolation matrix
1231 
1232   Note:
1233   The first member of the user context must be an FEMContext.
1234 
1235   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
1236   like a GPU, or vectorize on a multicore machine.
1237 
1238   Level: developer
1239 
1240 .seealso: DMPlexComputeJacobianFEM()
1241 @*/
1242 PetscErrorCode DMPlexComputeInterpolatorFEM(DM dmc, DM dmf, Mat In, void *user)
1243 {
1244   DM_Plex          *mesh  = (DM_Plex *) dmc->data;
1245   const char       *name  = "Interpolator";
1246   PetscDS           prob;
1247   PetscFE          *feRef;
1248   PetscFV          *fvRef;
1249   PetscSection      fsection, fglobalSection;
1250   PetscSection      csection, cglobalSection;
1251   PetscScalar      *elemMat;
1252   PetscInt          dim, Nf, f, fieldI, fieldJ, offsetI, offsetJ, cStart, cEnd, cEndInterior, c;
1253   PetscInt          cTotDim, rTotDim = 0;
1254   PetscErrorCode    ierr;
1255 
1256   PetscFunctionBegin;
1257   ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1258   ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr);
1259   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
1260   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
1261   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
1262   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
1263   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
1264   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
1265   ierr = DMPlexGetHybridBounds(dmc, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
1266   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
1267   ierr = DMGetDS(dmf, &prob);CHKERRQ(ierr);
1268   ierr = PetscCalloc2(Nf,&feRef,Nf,&fvRef);CHKERRQ(ierr);
1269   for (f = 0; f < Nf; ++f) {
1270     PetscObject  obj;
1271     PetscClassId id;
1272     PetscInt     rNb = 0, Nc = 0;
1273 
1274     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
1275     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1276     if (id == PETSCFE_CLASSID) {
1277       PetscFE fe = (PetscFE) obj;
1278 
1279       ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr);
1280       ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr);
1281       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1282     } else if (id == PETSCFV_CLASSID) {
1283       PetscFV        fv = (PetscFV) obj;
1284       PetscDualSpace Q;
1285 
1286       ierr = PetscFVRefine(fv, &fvRef[f]);CHKERRQ(ierr);
1287       ierr = PetscFVGetDualSpace(fvRef[f], &Q);CHKERRQ(ierr);
1288       ierr = PetscDualSpaceGetDimension(Q, &rNb);CHKERRQ(ierr);
1289       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
1290     }
1291     rTotDim += rNb*Nc;
1292   }
1293   ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr);
1294   ierr = PetscMalloc1(rTotDim*cTotDim,&elemMat);CHKERRQ(ierr);
1295   ierr = PetscMemzero(elemMat, rTotDim*cTotDim * sizeof(PetscScalar));CHKERRQ(ierr);
1296   for (fieldI = 0, offsetI = 0; fieldI < Nf; ++fieldI) {
1297     PetscDualSpace   Qref;
1298     PetscQuadrature  f;
1299     const PetscReal *qpoints, *qweights;
1300     PetscReal       *points;
1301     PetscInt         npoints = 0, Nc, Np, fpdim, i, k, p, d;
1302 
1303     /* Compose points from all dual basis functionals */
1304     if (feRef[fieldI]) {
1305       ierr = PetscFEGetDualSpace(feRef[fieldI], &Qref);CHKERRQ(ierr);
1306       ierr = PetscFEGetNumComponents(feRef[fieldI], &Nc);CHKERRQ(ierr);
1307     } else {
1308       ierr = PetscFVGetDualSpace(fvRef[fieldI], &Qref);CHKERRQ(ierr);
1309       ierr = PetscFVGetNumComponents(fvRef[fieldI], &Nc);CHKERRQ(ierr);
1310     }
1311     ierr = PetscDualSpaceGetDimension(Qref, &fpdim);CHKERRQ(ierr);
1312     for (i = 0; i < fpdim; ++i) {
1313       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1314       ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, NULL);CHKERRQ(ierr);
1315       npoints += Np;
1316     }
1317     ierr = PetscMalloc1(npoints*dim,&points);CHKERRQ(ierr);
1318     for (i = 0, k = 0; i < fpdim; ++i) {
1319       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1320       ierr = PetscQuadratureGetData(f, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr);
1321       for (p = 0; p < Np; ++p, ++k) for (d = 0; d < dim; ++d) points[k*dim+d] = qpoints[p*dim+d];
1322     }
1323 
1324     for (fieldJ = 0, offsetJ = 0; fieldJ < Nf; ++fieldJ) {
1325       PetscObject  obj;
1326       PetscClassId id;
1327       PetscReal   *B;
1328       PetscInt     NcJ = 0, cpdim = 0, j;
1329 
1330       ierr = PetscDSGetDiscretization(prob, fieldJ, &obj);CHKERRQ(ierr);
1331       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1332       if (id == PETSCFE_CLASSID) {
1333         PetscFE fe = (PetscFE) obj;
1334 
1335         /* Evaluate basis at points */
1336         ierr = PetscFEGetNumComponents(fe, &NcJ);CHKERRQ(ierr);
1337         ierr = PetscFEGetDimension(fe, &cpdim);CHKERRQ(ierr);
1338         /* For now, fields only interpolate themselves */
1339         if (fieldI == fieldJ) {
1340           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);
1341           ierr = PetscFEGetTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);
1342           for (i = 0, k = 0; i < fpdim; ++i) {
1343             ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1344             ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, &qweights);CHKERRQ(ierr);
1345             for (p = 0; p < Np; ++p, ++k) {
1346               for (j = 0; j < cpdim; ++j) {
1347                 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];
1348               }
1349             }
1350           }
1351           ierr = PetscFERestoreTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr);
1352         }
1353       } else if (id == PETSCFV_CLASSID) {
1354         PetscFV        fv = (PetscFV) obj;
1355 
1356         /* Evaluate constant function at points */
1357         ierr = PetscFVGetNumComponents(fv, &NcJ);CHKERRQ(ierr);
1358         cpdim = 1;
1359         /* For now, fields only interpolate themselves */
1360         if (fieldI == fieldJ) {
1361           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);
1362           for (i = 0, k = 0; i < fpdim; ++i) {
1363             ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1364             ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, &qweights);CHKERRQ(ierr);
1365             for (p = 0; p < Np; ++p, ++k) {
1366               for (j = 0; j < cpdim; ++j) {
1367                 for (c = 0; c < Nc; ++c) elemMat[(offsetI + i*Nc + c)*cTotDim + offsetJ + j*NcJ + c] += 1.0*qweights[p];
1368               }
1369             }
1370           }
1371         }
1372       }
1373       offsetJ += cpdim*NcJ;
1374     }
1375     offsetI += fpdim*Nc;
1376     ierr = PetscFree(points);CHKERRQ(ierr);
1377   }
1378   if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(0, name, rTotDim, cTotDim, elemMat);CHKERRQ(ierr);}
1379   /* Preallocate matrix */
1380   {
1381     PetscHashJK ht;
1382     PetscLayout rLayout;
1383     PetscInt   *dnz, *onz, *cellCIndices, *cellFIndices;
1384     PetscInt    locRows, rStart, rEnd, cell, r;
1385 
1386     ierr = MatGetLocalSize(In, &locRows, NULL);CHKERRQ(ierr);
1387     ierr = PetscLayoutCreate(PetscObjectComm((PetscObject) In), &rLayout);CHKERRQ(ierr);
1388     ierr = PetscLayoutSetLocalSize(rLayout, locRows);CHKERRQ(ierr);
1389     ierr = PetscLayoutSetBlockSize(rLayout, 1);CHKERRQ(ierr);
1390     ierr = PetscLayoutSetUp(rLayout);CHKERRQ(ierr);
1391     ierr = PetscLayoutGetRange(rLayout, &rStart, &rEnd);CHKERRQ(ierr);
1392     ierr = PetscLayoutDestroy(&rLayout);CHKERRQ(ierr);
1393     ierr = PetscCalloc4(locRows,&dnz,locRows,&onz,cTotDim,&cellCIndices,rTotDim,&cellFIndices);CHKERRQ(ierr);
1394     ierr = PetscHashJKCreate(&ht);CHKERRQ(ierr);
1395     for (cell = cStart; cell < cEnd; ++cell) {
1396       ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, cell, cellCIndices, cellFIndices);CHKERRQ(ierr);
1397       for (r = 0; r < rTotDim; ++r) {
1398         PetscHashJKKey  key;
1399         PetscHashJKIter missing, iter;
1400 
1401         key.j = cellFIndices[r];
1402         if (key.j < 0) continue;
1403         for (c = 0; c < cTotDim; ++c) {
1404           key.k = cellCIndices[c];
1405           if (key.k < 0) continue;
1406           ierr = PetscHashJKPut(ht, key, &missing, &iter);CHKERRQ(ierr);
1407           if (missing) {
1408             ierr = PetscHashJKSet(ht, iter, 1);CHKERRQ(ierr);
1409             if ((key.k >= rStart) && (key.k < rEnd)) ++dnz[key.j-rStart];
1410             else                                     ++onz[key.j-rStart];
1411           }
1412         }
1413       }
1414     }
1415     ierr = PetscHashJKDestroy(&ht);CHKERRQ(ierr);
1416     ierr = MatXAIJSetPreallocation(In, 1, dnz, onz, NULL, NULL);CHKERRQ(ierr);
1417     ierr = MatSetOption(In, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
1418     ierr = PetscFree4(dnz,onz,cellCIndices,cellFIndices);CHKERRQ(ierr);
1419   }
1420   /* Fill matrix */
1421   ierr = MatZeroEntries(In);CHKERRQ(ierr);
1422   for (c = cStart; c < cEnd; ++c) {
1423     ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES);CHKERRQ(ierr);
1424   }
1425   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);}
1426   ierr = PetscFree2(feRef,fvRef);CHKERRQ(ierr);
1427   ierr = PetscFree(elemMat);CHKERRQ(ierr);
1428   ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1429   ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1430   if (mesh->printFEM) {
1431     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
1432     ierr = MatChop(In, 1.0e-10);CHKERRQ(ierr);
1433     ierr = MatView(In, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1434   }
1435   ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1436   PetscFunctionReturn(0);
1437 }
1438 
1439 #undef __FUNCT__
1440 #define __FUNCT__ "DMPlexComputeInjectorFEM"
1441 PetscErrorCode DMPlexComputeInjectorFEM(DM dmc, DM dmf, VecScatter *sc, void *user)
1442 {
1443   PetscDS        prob;
1444   PetscFE       *feRef;
1445   PetscFV       *fvRef;
1446   Vec            fv, cv;
1447   IS             fis, cis;
1448   PetscSection   fsection, fglobalSection, csection, cglobalSection;
1449   PetscInt      *cmap, *cellCIndices, *cellFIndices, *cindices, *findices;
1450   PetscInt       cTotDim, fTotDim = 0, Nf, f, field, cStart, cEnd, cEndInterior, c, dim, d, startC, offsetC, offsetF, m;
1451   PetscErrorCode ierr;
1452 
1453   PetscFunctionBegin;
1454   ierr = PetscLogEventBegin(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1455   ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr);
1456   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
1457   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
1458   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
1459   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
1460   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
1461   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
1462   ierr = DMPlexGetHybridBounds(dmc, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
1463   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
1464   ierr = DMGetDS(dmc, &prob);CHKERRQ(ierr);
1465   ierr = PetscCalloc2(Nf,&feRef,Nf,&fvRef);CHKERRQ(ierr);
1466   for (f = 0; f < Nf; ++f) {
1467     PetscObject  obj;
1468     PetscClassId id;
1469     PetscInt     fNb = 0, Nc = 0;
1470 
1471     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
1472     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
1473     if (id == PETSCFE_CLASSID) {
1474       PetscFE fe = (PetscFE) obj;
1475 
1476       ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr);
1477       ierr = PetscFEGetDimension(feRef[f], &fNb);CHKERRQ(ierr);
1478       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1479     } else if (id == PETSCFV_CLASSID) {
1480       PetscFV        fv = (PetscFV) obj;
1481       PetscDualSpace Q;
1482 
1483       ierr = PetscFVRefine(fv, &fvRef[f]);CHKERRQ(ierr);
1484       ierr = PetscFVGetDualSpace(fvRef[f], &Q);CHKERRQ(ierr);
1485       ierr = PetscDualSpaceGetDimension(Q, &fNb);CHKERRQ(ierr);
1486       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
1487     }
1488     fTotDim += fNb*Nc;
1489   }
1490   ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr);
1491   ierr = PetscMalloc1(cTotDim,&cmap);CHKERRQ(ierr);
1492   for (field = 0, offsetC = 0, offsetF = 0; field < Nf; ++field) {
1493     PetscFE        feC;
1494     PetscFV        fvC;
1495     PetscDualSpace QF, QC;
1496     PetscInt       NcF, NcC, fpdim, cpdim;
1497 
1498     if (feRef[field]) {
1499       ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &feC);CHKERRQ(ierr);
1500       ierr = PetscFEGetNumComponents(feC, &NcC);CHKERRQ(ierr);
1501       ierr = PetscFEGetNumComponents(feRef[field], &NcF);CHKERRQ(ierr);
1502       ierr = PetscFEGetDualSpace(feRef[field], &QF);CHKERRQ(ierr);
1503       ierr = PetscDualSpaceGetDimension(QF, &fpdim);CHKERRQ(ierr);
1504       ierr = PetscFEGetDualSpace(feC, &QC);CHKERRQ(ierr);
1505       ierr = PetscDualSpaceGetDimension(QC, &cpdim);CHKERRQ(ierr);
1506     } else {
1507       ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fvC);CHKERRQ(ierr);
1508       ierr = PetscFVGetNumComponents(fvC, &NcC);CHKERRQ(ierr);
1509       ierr = PetscFVGetNumComponents(fvRef[field], &NcF);CHKERRQ(ierr);
1510       ierr = PetscFVGetDualSpace(fvRef[field], &QF);CHKERRQ(ierr);
1511       ierr = PetscDualSpaceGetDimension(QF, &fpdim);CHKERRQ(ierr);
1512       ierr = PetscFVGetDualSpace(fvC, &QC);CHKERRQ(ierr);
1513       ierr = PetscDualSpaceGetDimension(QC, &cpdim);CHKERRQ(ierr);
1514     }
1515     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);
1516     for (c = 0; c < cpdim; ++c) {
1517       PetscQuadrature  cfunc;
1518       const PetscReal *cqpoints;
1519       PetscInt         NpC;
1520       PetscBool        found = PETSC_FALSE;
1521 
1522       ierr = PetscDualSpaceGetFunctional(QC, c, &cfunc);CHKERRQ(ierr);
1523       ierr = PetscQuadratureGetData(cfunc, NULL, &NpC, &cqpoints, NULL);CHKERRQ(ierr);
1524       if (NpC != 1 && feRef[field]) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not know how to do injection for moments");
1525       for (f = 0; f < fpdim; ++f) {
1526         PetscQuadrature  ffunc;
1527         const PetscReal *fqpoints;
1528         PetscReal        sum = 0.0;
1529         PetscInt         NpF, comp;
1530 
1531         ierr = PetscDualSpaceGetFunctional(QF, f, &ffunc);CHKERRQ(ierr);
1532         ierr = PetscQuadratureGetData(ffunc, NULL, &NpF, &fqpoints, NULL);CHKERRQ(ierr);
1533         if (NpC != NpF) continue;
1534         for (d = 0; d < dim; ++d) sum += PetscAbsReal(cqpoints[d] - fqpoints[d]);
1535         if (sum > 1.0e-9) continue;
1536         for (comp = 0; comp < NcC; ++comp) {
1537           cmap[(offsetC+c)*NcC+comp] = (offsetF+f)*NcF+comp;
1538         }
1539         found = PETSC_TRUE;
1540         break;
1541       }
1542       if (!found) {
1543         /* TODO We really want the average here, but some asshole put VecScatter in the interface */
1544         if (fvRef[field]) {
1545           PetscInt comp;
1546           for (comp = 0; comp < NcC; ++comp) {
1547             cmap[(offsetC+c)*NcC+comp] = (offsetF+0)*NcF+comp;
1548           }
1549         } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not locate matching functional for injection");
1550       }
1551     }
1552     offsetC += cpdim*NcC;
1553     offsetF += fpdim*NcF;
1554   }
1555   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);ierr = PetscFVDestroy(&fvRef[f]);CHKERRQ(ierr);}
1556   ierr = PetscFree2(feRef,fvRef);CHKERRQ(ierr);
1557 
1558   ierr = DMGetGlobalVector(dmf, &fv);CHKERRQ(ierr);
1559   ierr = DMGetGlobalVector(dmc, &cv);CHKERRQ(ierr);
1560   ierr = VecGetOwnershipRange(cv, &startC, NULL);CHKERRQ(ierr);
1561   ierr = PetscSectionGetConstrainedStorageSize(cglobalSection, &m);CHKERRQ(ierr);
1562   ierr = PetscMalloc2(cTotDim,&cellCIndices,fTotDim,&cellFIndices);CHKERRQ(ierr);
1563   ierr = PetscMalloc1(m,&cindices);CHKERRQ(ierr);
1564   ierr = PetscMalloc1(m,&findices);CHKERRQ(ierr);
1565   for (d = 0; d < m; ++d) cindices[d] = findices[d] = -1;
1566   for (c = cStart; c < cEnd; ++c) {
1567     ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, c, cellCIndices, cellFIndices);CHKERRQ(ierr);
1568     for (d = 0; d < cTotDim; ++d) {
1569       if (cellCIndices[d] < 0) continue;
1570       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]]);
1571       cindices[cellCIndices[d]-startC] = cellCIndices[d];
1572       findices[cellCIndices[d]-startC] = cellFIndices[cmap[d]];
1573     }
1574   }
1575   ierr = PetscFree(cmap);CHKERRQ(ierr);
1576   ierr = PetscFree2(cellCIndices,cellFIndices);CHKERRQ(ierr);
1577 
1578   ierr = ISCreateGeneral(PETSC_COMM_SELF, m, cindices, PETSC_OWN_POINTER, &cis);CHKERRQ(ierr);
1579   ierr = ISCreateGeneral(PETSC_COMM_SELF, m, findices, PETSC_OWN_POINTER, &fis);CHKERRQ(ierr);
1580   ierr = VecScatterCreate(cv, cis, fv, fis, sc);CHKERRQ(ierr);
1581   ierr = ISDestroy(&cis);CHKERRQ(ierr);
1582   ierr = ISDestroy(&fis);CHKERRQ(ierr);
1583   ierr = DMRestoreGlobalVector(dmf, &fv);CHKERRQ(ierr);
1584   ierr = DMRestoreGlobalVector(dmc, &cv);CHKERRQ(ierr);
1585   ierr = PetscLogEventEnd(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1586   PetscFunctionReturn(0);
1587 }
1588