xref: /petsc/src/dm/impls/plex/plexfem.c (revision ee2838f6c52a678b46227921cbfd65cc7c49a7bb)
1 #include <petsc-private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2 
3 #include <petsc-private/petscfeimpl.h>
4 #include <petsc-private/petscfvimpl.h>
5 
6 #undef __FUNCT__
7 #define __FUNCT__ "DMPlexGetScale"
8 PetscErrorCode DMPlexGetScale(DM dm, PetscUnit unit, PetscReal *scale)
9 {
10   DM_Plex *mesh = (DM_Plex*) dm->data;
11 
12   PetscFunctionBegin;
13   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
14   PetscValidPointer(scale, 3);
15   *scale = mesh->scale[unit];
16   PetscFunctionReturn(0);
17 }
18 
19 #undef __FUNCT__
20 #define __FUNCT__ "DMPlexSetScale"
21 PetscErrorCode DMPlexSetScale(DM dm, PetscUnit unit, PetscReal scale)
22 {
23   DM_Plex *mesh = (DM_Plex*) dm->data;
24 
25   PetscFunctionBegin;
26   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
27   mesh->scale[unit] = scale;
28   PetscFunctionReturn(0);
29 }
30 
31 PETSC_STATIC_INLINE PetscInt epsilon(PetscInt i, PetscInt j, PetscInt k)
32 {
33   switch (i) {
34   case 0:
35     switch (j) {
36     case 0: return 0;
37     case 1:
38       switch (k) {
39       case 0: return 0;
40       case 1: return 0;
41       case 2: return 1;
42       }
43     case 2:
44       switch (k) {
45       case 0: return 0;
46       case 1: return -1;
47       case 2: return 0;
48       }
49     }
50   case 1:
51     switch (j) {
52     case 0:
53       switch (k) {
54       case 0: return 0;
55       case 1: return 0;
56       case 2: return -1;
57       }
58     case 1: return 0;
59     case 2:
60       switch (k) {
61       case 0: return 1;
62       case 1: return 0;
63       case 2: return 0;
64       }
65     }
66   case 2:
67     switch (j) {
68     case 0:
69       switch (k) {
70       case 0: return 0;
71       case 1: return 1;
72       case 2: return 0;
73       }
74     case 1:
75       switch (k) {
76       case 0: return -1;
77       case 1: return 0;
78       case 2: return 0;
79       }
80     case 2: return 0;
81     }
82   }
83   return 0;
84 }
85 
86 #undef __FUNCT__
87 #define __FUNCT__ "DMPlexCreateRigidBody"
88 /*@C
89   DMPlexCreateRigidBody - create rigid body modes from coordinates
90 
91   Collective on DM
92 
93   Input Arguments:
94 + dm - the DM
95 . section - the local section associated with the rigid field, or NULL for the default section
96 - globalSection - the global section associated with the rigid field, or NULL for the default section
97 
98   Output Argument:
99 . sp - the null space
100 
101   Note: This is necessary to take account of Dirichlet conditions on the displacements
102 
103   Level: advanced
104 
105 .seealso: MatNullSpaceCreate()
106 @*/
107 PetscErrorCode DMPlexCreateRigidBody(DM dm, PetscSection section, PetscSection globalSection, MatNullSpace *sp)
108 {
109   MPI_Comm       comm;
110   Vec            coordinates, localMode, mode[6];
111   PetscSection   coordSection;
112   PetscScalar   *coords;
113   PetscInt       dim, vStart, vEnd, v, n, m, d, i, j;
114   PetscErrorCode ierr;
115 
116   PetscFunctionBegin;
117   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
118   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
119   if (dim == 1) {
120     ierr = MatNullSpaceCreate(comm, PETSC_TRUE, 0, NULL, sp);CHKERRQ(ierr);
121     PetscFunctionReturn(0);
122   }
123   if (!section)       {ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);}
124   if (!globalSection) {ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);}
125   ierr = PetscSectionGetConstrainedStorageSize(globalSection, &n);CHKERRQ(ierr);
126   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
127   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
128   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
129   m    = (dim*(dim+1))/2;
130   ierr = VecCreate(comm, &mode[0]);CHKERRQ(ierr);
131   ierr = VecSetSizes(mode[0], n, PETSC_DETERMINE);CHKERRQ(ierr);
132   ierr = VecSetUp(mode[0]);CHKERRQ(ierr);
133   for (i = 1; i < m; ++i) {ierr = VecDuplicate(mode[0], &mode[i]);CHKERRQ(ierr);}
134   /* Assume P1 */
135   ierr = DMGetLocalVector(dm, &localMode);CHKERRQ(ierr);
136   for (d = 0; d < dim; ++d) {
137     PetscScalar values[3] = {0.0, 0.0, 0.0};
138 
139     values[d] = 1.0;
140     ierr      = VecSet(localMode, 0.0);CHKERRQ(ierr);
141     for (v = vStart; v < vEnd; ++v) {
142       ierr = DMPlexVecSetClosure(dm, section, localMode, v, values, INSERT_VALUES);CHKERRQ(ierr);
143     }
144     ierr = DMLocalToGlobalBegin(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
145     ierr = DMLocalToGlobalEnd(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
146   }
147   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
148   for (d = dim; d < dim*(dim+1)/2; ++d) {
149     PetscInt i, j, k = dim > 2 ? d - dim : d;
150 
151     ierr = VecSet(localMode, 0.0);CHKERRQ(ierr);
152     for (v = vStart; v < vEnd; ++v) {
153       PetscScalar values[3] = {0.0, 0.0, 0.0};
154       PetscInt    off;
155 
156       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
157       for (i = 0; i < dim; ++i) {
158         for (j = 0; j < dim; ++j) {
159           values[j] += epsilon(i, j, k)*PetscRealPart(coords[off+i]);
160         }
161       }
162       ierr = DMPlexVecSetClosure(dm, section, localMode, v, values, INSERT_VALUES);CHKERRQ(ierr);
163     }
164     ierr = DMLocalToGlobalBegin(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
165     ierr = DMLocalToGlobalEnd(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
166   }
167   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
168   ierr = DMRestoreLocalVector(dm, &localMode);CHKERRQ(ierr);
169   for (i = 0; i < dim; ++i) {ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);}
170   /* Orthonormalize system */
171   for (i = dim; i < m; ++i) {
172     PetscScalar dots[6];
173 
174     ierr = VecMDot(mode[i], i, mode, dots);CHKERRQ(ierr);
175     for (j = 0; j < i; ++j) dots[j] *= -1.0;
176     ierr = VecMAXPY(mode[i], i, dots, mode);CHKERRQ(ierr);
177     ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);
178   }
179   ierr = MatNullSpaceCreate(comm, PETSC_FALSE, m, mode, sp);CHKERRQ(ierr);
180   for (i = 0; i< m; ++i) {ierr = VecDestroy(&mode[i]);CHKERRQ(ierr);}
181   PetscFunctionReturn(0);
182 }
183 
184 #undef __FUNCT__
185 #define __FUNCT__ "DMPlexSetMaxProjectionHeight"
186 /*@
187   DMPlexSetMaxProjectionHeight - In DMPlexProjectXXXLocal() functions, the projected values of a basis function's dofs
188   are computed by associating the basis function with one of the mesh points in its transitively-closed support, and
189   evaluating the dual space basis of that point.  A basis function is associated with the point in its
190   transitively-closed support whose mesh height is highest (w.r.t. DAG height), but not greater than the maximum
191   projection height, which is set with this function.  By default, the maximum projection height is zero, which means
192   that only mesh cells are used to project basis functions.  A height of one, for example, evaluates a cell-interior
193   basis functions using its cells dual space basis, but all other basis functions with the dual space basis of a face.
194 
195   Input Parameters:
196 + dm - the DMPlex object
197 - height - the maximum projection height >= 0
198 
199   Level: advanced
200 
201 .seealso: DMPlexGetMaxProjectionHeight(), DMPlexProjectFunctionLocal(), DMPlexProjectFunctionLabelLocal()
202 @*/
203 PetscErrorCode DMPlexSetMaxProjectionHeight(DM dm, PetscInt height)
204 {
205   DM_Plex *plex = (DM_Plex *) dm->data;
206 
207   PetscFunctionBegin;
208   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
209   plex->maxProjectionHeight = height;
210   PetscFunctionReturn(0);
211 }
212 
213 #undef __FUNCT__
214 #define __FUNCT__ "DMPlexGetMaxProjectionHeight"
215 /*@
216   DMPlexGetMaxProjectionHeight - Get the maximum height (w.r.t. DAG) of mesh points used to evaluate dual bases in
217   DMPlexProjectXXXLocal() functions.
218 
219   Input Parameters:
220 . dm - the DMPlex object
221 
222   Output Parameters:
223 . height - the maximum projection height
224 
225   Level: intermediate
226 
227 .seealso: DMPlexSetMaxProjectionHeight(), DMPlexProjectFunctionLocal(), DMPlexProjectFunctionLabelLocal()
228 @*/
229 PetscErrorCode DMPlexGetMaxProjectionHeight(DM dm, PetscInt *height)
230 {
231   DM_Plex *plex = (DM_Plex *) dm->data;
232 
233   PetscFunctionBegin;
234   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
235   *height = plex->maxProjectionHeight;
236   PetscFunctionReturn(0);
237 }
238 
239 #undef __FUNCT__
240 #define __FUNCT__ "DMPlexProjectFunctionLabelLocal"
241 PetscErrorCode DMPlexProjectFunctionLabelLocal(DM dm, DMLabel label, PetscInt numIds, const PetscInt ids[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
242 {
243   PetscFE         fe;
244   PetscDualSpace *sp, *cellsp;
245   PetscSection    section;
246   PetscScalar    *values;
247   PetscBool      *fieldActive;
248   PetscInt        numFields, numComp, dim, dimEmbed, spDim, totDim, numValues, pStart, pEnd, f, d, v, i, comp, maxHeight, h;
249   PetscErrorCode  ierr;
250 
251   PetscFunctionBegin;
252   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
253   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
254   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
255   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
256   ierr = PetscMalloc1(numFields,&sp);CHKERRQ(ierr);
257   ierr = DMPlexGetMaxProjectionHeight(dm,&maxHeight);CHKERRQ(ierr);
258   if (maxHeight < 0 || maxHeight > dim) {SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"maximum projection height %d not in [0, %d)\n", maxHeight,dim);}
259   if (maxHeight > 0) {
260     ierr = PetscMalloc1(numFields,&cellsp);CHKERRQ(ierr);
261   }
262   else {
263     cellsp = sp;
264   }
265   for (h = 0; h <= maxHeight; h++) {
266     ierr = DMPlexGetHeightStratum(dm, h, &pStart, &pEnd);CHKERRQ(ierr);
267     if (pEnd <= pStart) continue;
268     totDim = 0;
269     for (f = 0; f < numFields; ++f) {
270       ierr = DMGetField(dm, f, (PetscObject *) &fe);CHKERRQ(ierr);
271       if (!h) {
272         ierr = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr);
273         sp[f] = cellsp[f];
274       }
275       else {
276         ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], h, &sp[f]);CHKERRQ(ierr);
277         if (!sp[f]) continue;
278       }
279       ierr = PetscFEGetNumComponents(fe, &numComp);CHKERRQ(ierr);
280       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
281       totDim += spDim*numComp;
282     }
283     ierr = DMPlexVecGetClosure(dm, section, localX, pStart, &numValues, NULL);CHKERRQ(ierr);
284     if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section point closure size %d != dual space dimension %d", numValues, totDim);
285     if (!totDim) continue;
286     ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
287     ierr = DMGetWorkArray(dm, numFields, PETSC_BOOL, &fieldActive);CHKERRQ(ierr);
288     for (f = 0; f < numFields; ++f) fieldActive[f] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE;
289     for (i = 0; i < numIds; ++i) {
290       IS              pointIS;
291       const PetscInt *points;
292       PetscInt        n, p;
293 
294       ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr);
295       ierr = ISGetLocalSize(pointIS, &n);CHKERRQ(ierr);
296       ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
297       for (p = 0; p < n; ++p) {
298         const PetscInt    point = points[p];
299         PetscFECellGeom   geom;
300 
301         if ((point < pStart) || (point >= pEnd)) continue;
302         ierr          = DMPlexComputeCellGeometryFEM(dm, point, NULL, geom.v0, geom.J, NULL, &geom.detJ);CHKERRQ(ierr);
303         geom.dim      = dim - h;
304         geom.dimEmbed = dimEmbed;
305         for (f = 0, v = 0; f < numFields; ++f) {
306           void * const ctx = ctxs ? ctxs[f] : NULL;
307 
308           if (!sp[f]) continue;
309           ierr = DMGetField(dm, f, (PetscObject *) &fe);CHKERRQ(ierr);
310           ierr = PetscFEGetNumComponents(fe, &numComp);CHKERRQ(ierr);
311           ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
312           for (d = 0; d < spDim; ++d) {
313             if (funcs[f]) {
314               ierr = PetscDualSpaceApply(sp[f], d, &geom, numComp, funcs[f], ctx, &values[v]);CHKERRQ(ierr);
315             } else {
316               for (comp = 0; comp < numComp; ++comp) values[v+comp] = 0.0;
317             }
318             v += numComp;
319           }
320         }
321         ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, point, values, mode);CHKERRQ(ierr);
322       }
323       ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
324       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
325     }
326     if (h) {
327       for (f = 0; f < numFields; ++f) {ierr = PetscDualSpaceDestroy(&sp[f]);CHKERRQ(ierr);}
328     }
329   }
330   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
331   ierr = DMRestoreWorkArray(dm, numFields, PETSC_BOOL, &fieldActive);CHKERRQ(ierr);
332   ierr = PetscFree(sp);CHKERRQ(ierr);
333   if (maxHeight > 0) {
334     ierr = PetscFree(cellsp);CHKERRQ(ierr);
335   }
336   PetscFunctionReturn(0);
337 }
338 
339 #undef __FUNCT__
340 #define __FUNCT__ "DMPlexProjectFunctionLocal"
341 PetscErrorCode DMPlexProjectFunctionLocal(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
342 {
343   PetscDualSpace *sp, *cellsp;
344   PetscInt       *numComp;
345   PetscSection    section;
346   PetscScalar    *values;
347   PetscInt        numFields, dim, spDim, totDim, numValues, pStart, pEnd, p, f, d, v, comp, h, maxHeight;
348   PetscErrorCode  ierr;
349 
350   PetscFunctionBegin;
351   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
352   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
353   ierr = PetscMalloc2(numFields, &sp, numFields, &numComp);CHKERRQ(ierr);
354   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
355   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
356   ierr = DMPlexGetMaxProjectionHeight(dm,&maxHeight);CHKERRQ(ierr);
357   if (maxHeight < 0 || maxHeight > dim) {SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"maximum projection height %d not in [0, %d)\n", maxHeight,dim);}
358   if (maxHeight > 0) {
359     ierr = PetscMalloc1(numFields,&cellsp);CHKERRQ(ierr);
360   }
361   else {
362     cellsp = sp;
363   }
364   for (h = 0; h <= maxHeight; h++) {
365     ierr = DMPlexGetHeightStratum(dm, h, &pStart, &pEnd);CHKERRQ(ierr);
366     if (pEnd <= pStart) continue;
367     totDim = 0;
368     for (f = 0; f < numFields; ++f) {
369       PetscObject  obj;
370       PetscClassId id;
371 
372       ierr = DMGetField(dm, f, &obj);CHKERRQ(ierr);
373       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
374       if (id == PETSCFE_CLASSID) {
375         PetscFE fe = (PetscFE) obj;
376 
377         ierr = PetscFEGetNumComponents(fe, &numComp[f]);CHKERRQ(ierr);
378         if (!h) {
379           ierr  = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr);
380           sp[f] = cellsp[f];
381           ierr  = PetscObjectReference((PetscObject) sp[f]);CHKERRQ(ierr);
382         }
383         else {
384           ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], h, &sp[f]);CHKERRQ(ierr);
385           if (!sp[f]) {
386             continue;
387           }
388         }
389       } else if (id == PETSCFV_CLASSID) {
390         PetscFV         fv = (PetscFV) obj;
391         PetscQuadrature q;
392 
393         if (h) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP, "Projection height > 0 not supported for finite volume");
394         ierr = PetscFVGetNumComponents(fv, &numComp[f]);CHKERRQ(ierr);
395         ierr = PetscFVGetQuadrature(fv, &q);CHKERRQ(ierr);
396         ierr = PetscDualSpaceCreate(PetscObjectComm((PetscObject) fv), &sp[f]);CHKERRQ(ierr);
397         ierr = PetscDualSpaceSetDM(sp[f], dm);CHKERRQ(ierr);
398         ierr = PetscDualSpaceSetType(sp[f], PETSCDUALSPACESIMPLE);CHKERRQ(ierr);
399         ierr = PetscDualSpaceSimpleSetDimension(sp[f], 1);CHKERRQ(ierr);
400         ierr = PetscDualSpaceSimpleSetFunctional(sp[f], 0, q);CHKERRQ(ierr);
401       } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
402       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
403       totDim += spDim*numComp[f];
404     }
405     ierr = DMPlexVecGetClosure(dm, section, localX, pStart, &numValues, NULL);CHKERRQ(ierr);
406     if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section point closure size %d != dual space dimension %d", numValues, totDim);
407     if (!totDim) continue;
408     ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
409     for (p = pStart; p < pEnd; ++p) {
410       PetscFECellGeom geom;
411 
412       ierr          = DMPlexComputeCellGeometryFEM(dm, p, NULL, &geom.v0, &geom.J, NULL, &geom.detJ);CHKERRQ(ierr);
413       geom.dim      = dim = h;
414       geom.dimEmbed = dimEmbed;
415       for (f = 0, v = 0; f < numFields; ++f) {
416         PetscFE      fe;
417         void * const ctx = ctxs ? ctxs[f] : NULL;
418 
419         if (!sp[f]) continue;
420         ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
421         for (d = 0; d < spDim; ++d) {
422           if (funcs[f]) {
423             ierr = PetscDualSpaceApply(sp[f], d, &geom, numComp[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr);
424           } else {
425             for (comp = 0; comp < numComp[f]; ++comp) values[v+comp] = 0.0;
426           }
427           v += numComp[f];
428         }
429       }
430       ierr = DMPlexVecSetClosure(dm, section, localX, p, values, mode);CHKERRQ(ierr);
431     }
432     ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
433     if (h || !maxHeight) {
434       for (f = 0; f < numFields; f++) {ierr = PetscDualSpaceDestroy(&sp[f]);CHKERRQ(ierr);}
435     }
436   }
437   ierr = PetscFree2(sp, numComp);CHKERRQ(ierr);
438   if (maxHeight > 0) {
439     for (f = 0; f < numFields; f++) {ierr = PetscDualSpaceDestroy(&cellsp[f]);CHKERRQ(ierr);}
440     ierr = PetscFree(cellsp);CHKERRQ(ierr);
441   }
442   PetscFunctionReturn(0);
443 }
444 
445 #undef __FUNCT__
446 #define __FUNCT__ "DMPlexProjectFieldLocal"
447 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)
448 {
449   DM                dmAux;
450   PetscDS           prob, probAux;
451   Vec               A;
452   PetscSection      section, sectionAux;
453   PetscScalar      *values, *u, *u_x, *a, *a_x;
454   PetscReal        *x, *v0, *J, *invJ, detJ, **basisField, **basisFieldDer, **basisFieldAux, **basisFieldDerAux;
455   PetscInt          Nf, dim, spDim, totDim, numValues, cStart, cEnd, c, f, d, v, comp, maxHeight;
456   PetscErrorCode    ierr;
457 
458   PetscFunctionBegin;
459   ierr = DMPlexGetMaxProjectionHeight(dm,&maxHeight);CHKERRQ(ierr);
460   if (maxHeight > 0) {SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Field projection for height > 0 not supported yet");}
461   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
462   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
463   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
464   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
465   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
466   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
467   ierr = PetscDSGetTabulation(prob, &basisField, &basisFieldDer);CHKERRQ(ierr);
468   ierr = PetscDSGetEvaluationArrays(prob, &u, NULL, &u_x);CHKERRQ(ierr);
469   ierr = PetscDSGetRefCoordArrays(prob, &x, NULL);CHKERRQ(ierr);
470   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
471   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
472   if (dmAux) {
473     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
474     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
475     ierr = PetscDSGetTabulation(prob, &basisFieldAux, &basisFieldDerAux);CHKERRQ(ierr);
476     ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);CHKERRQ(ierr);
477   }
478   ierr = DMPlexInsertBoundaryValues(dm, localU, 0.0, NULL, NULL, NULL);CHKERRQ(ierr);
479   ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr);
480   if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section cell closure size %d != dual space dimension %d", numValues, totDim);
481   ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
482   ierr = PetscMalloc3(dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
483   for (c = cStart; c < cEnd; ++c) {
484     PetscScalar *coefficients = NULL, *coefficientsAux = NULL;
485 
486     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
487     ierr = DMPlexVecGetClosure(dm, section, localU, c, NULL, &coefficients);CHKERRQ(ierr);
488     if (dmAux) {ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &coefficientsAux);CHKERRQ(ierr);}
489     for (f = 0, v = 0; f < Nf; ++f) {
490       PetscFE        fe;
491       PetscDualSpace sp;
492       PetscInt       Ncf;
493 
494       ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
495       ierr = PetscFEGetDualSpace(fe, &sp);CHKERRQ(ierr);
496       ierr = PetscFEGetNumComponents(fe, &Ncf);CHKERRQ(ierr);
497       ierr = PetscDualSpaceGetDimension(sp, &spDim);CHKERRQ(ierr);
498       for (d = 0; d < spDim; ++d) {
499         PetscQuadrature  quad;
500         const PetscReal *points, *weights;
501         PetscInt         numPoints, q;
502 
503         if (funcs[f]) {
504           ierr = PetscDualSpaceGetFunctional(sp, d, &quad);CHKERRQ(ierr);
505           ierr = PetscQuadratureGetData(quad, NULL, &numPoints, &points, &weights);CHKERRQ(ierr);
506           ierr = PetscFEGetTabulation(fe, numPoints, points, &basisField[f], &basisFieldDer[f], NULL);CHKERRQ(ierr);
507           for (q = 0; q < numPoints; ++q) {
508             CoordinatesRefToReal(dim, dim, v0, J, &points[q*dim], x);
509             ierr = EvaluateFieldJets(prob,    PETSC_FALSE, q, invJ, coefficients,    NULL, u, u_x, NULL);CHKERRQ(ierr);
510             ierr = EvaluateFieldJets(probAux, PETSC_FALSE, q, invJ, coefficientsAux, NULL, a, a_x, NULL);CHKERRQ(ierr);
511             (*funcs[f])(u, NULL, u_x, a, NULL, a_x, x, &values[v]);
512           }
513           ierr = PetscFERestoreTabulation(fe, numPoints, points, &basisField[f], &basisFieldDer[f], NULL);CHKERRQ(ierr);
514         } else {
515           for (comp = 0; comp < Ncf; ++comp) values[v+comp] = 0.0;
516         }
517         v += Ncf;
518       }
519     }
520     ierr = DMPlexVecRestoreClosure(dm, section, localU, c, NULL, &coefficients);CHKERRQ(ierr);
521     if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &coefficientsAux);CHKERRQ(ierr);}
522     ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr);
523   }
524   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
525   ierr = PetscFree3(v0,J,invJ);CHKERRQ(ierr);
526   PetscFunctionReturn(0);
527 }
528 
529 #undef __FUNCT__
530 #define __FUNCT__ "DMPlexInsertBoundaryValues_FEM_Internal"
531 static PetscErrorCode DMPlexInsertBoundaryValues_FEM_Internal(DM dm, PetscInt field, DMLabel label, PetscInt numids, const PetscInt ids[], void (*func)(const PetscReal[], PetscScalar *, void *), void *ctx, Vec locX)
532 {
533   void        (**funcs)(const PetscReal x[], PetscScalar *u, void *ctx);
534   void         **ctxs;
535   PetscInt       numFields;
536   PetscErrorCode ierr;
537 
538   PetscFunctionBegin;
539   ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr);
540   ierr = PetscCalloc2(numFields,&funcs,numFields,&ctxs);CHKERRQ(ierr);
541   funcs[field] = func;
542   ctxs[field]  = ctx;
543   ierr = DMPlexProjectFunctionLabelLocal(dm, label, numids, ids, funcs, ctxs, INSERT_BC_VALUES, locX);CHKERRQ(ierr);
544   ierr = PetscFree2(funcs,ctxs);CHKERRQ(ierr);
545   PetscFunctionReturn(0);
546 }
547 
548 #undef __FUNCT__
549 #define __FUNCT__ "DMPlexInsertBoundaryValues_FVM_Internal"
550 static PetscErrorCode DMPlexInsertBoundaryValues_FVM_Internal(DM dm, PetscReal time, Vec faceGeometry, Vec cellGeometry, Vec Grad,
551                                                               PetscInt field, DMLabel label, PetscInt numids, const PetscInt ids[], PetscErrorCode (*func)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*), void *ctx, Vec locX)
552 {
553   PetscFV            fv;
554   DM                 dmFace, dmCell, dmGrad;
555   const PetscScalar *facegeom, *cellgeom, *grad;
556   PetscScalar       *x, *fx;
557   PetscInt           dim, fStart, fEnd, pdim, i;
558   PetscErrorCode     ierr;
559 
560   PetscFunctionBegin;
561   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
562   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
563   ierr = DMGetField(dm, field, (PetscObject *) &fv);CHKERRQ(ierr);
564   ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr);
565   ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
566   ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr);
567   ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
568   if (Grad) {
569     ierr = VecGetDM(Grad, &dmGrad);CHKERRQ(ierr);
570     ierr = VecGetArrayRead(Grad, &grad);CHKERRQ(ierr);
571     ierr = PetscFVGetNumComponents(fv, &pdim);CHKERRQ(ierr);
572     ierr = DMGetWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr);
573   }
574   ierr = VecGetArray(locX, &x);CHKERRQ(ierr);
575   for (i = 0; i < numids; ++i) {
576     IS              faceIS;
577     const PetscInt *faces;
578     PetscInt        numFaces, f;
579 
580     ierr = DMLabelGetStratumIS(label, ids[i], &faceIS);CHKERRQ(ierr);
581     if (!faceIS) continue; /* No points with that id on this process */
582     ierr = ISGetLocalSize(faceIS, &numFaces);CHKERRQ(ierr);
583     ierr = ISGetIndices(faceIS, &faces);CHKERRQ(ierr);
584     for (f = 0; f < numFaces; ++f) {
585       const PetscInt         face = faces[f], *cells;
586       const PetscFVFaceGeom *fg;
587 
588       if ((face < fStart) || (face >= fEnd)) continue; /* Refinement adds non-faces to labels */
589       ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr);
590       ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr);
591       if (Grad) {
592         const PetscFVCellGeom *cg;
593         const PetscScalar     *cx, *cgrad;
594         PetscScalar           *xG;
595         PetscReal              dx[3];
596         PetscInt               d;
597 
598         ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cg);CHKERRQ(ierr);
599         ierr = DMPlexPointLocalRead(dm, cells[0], x, &cx);CHKERRQ(ierr);
600         ierr = DMPlexPointLocalRead(dmGrad, cells[0], grad, &cgrad);CHKERRQ(ierr);
601         ierr = DMPlexPointLocalRef(dm, cells[1], x, &xG);CHKERRQ(ierr);
602         DMPlex_WaxpyD_Internal(dim, -1, cg->centroid, fg->centroid, dx);
603         for (d = 0; d < pdim; ++d) fx[d] = cx[d] + DMPlex_DotD_Internal(dim, &cgrad[d*dim], dx);
604         ierr = (*func)(time, fg->centroid, fg->normal, fx, xG, ctx);CHKERRQ(ierr);
605       } else {
606         const PetscScalar *xI;
607         PetscScalar       *xG;
608 
609         ierr = DMPlexPointLocalRead(dm, cells[0], x, &xI);CHKERRQ(ierr);
610         ierr = DMPlexPointLocalRef(dm, cells[1], x, &xG);CHKERRQ(ierr);
611         ierr = (*func)(time, fg->centroid, fg->normal, xI, xG, ctx);CHKERRQ(ierr);
612       }
613     }
614     ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr);
615     ierr = ISDestroy(&faceIS);CHKERRQ(ierr);
616   }
617   ierr = VecRestoreArray(locX, &x);CHKERRQ(ierr);
618   if (Grad) {
619     ierr = DMRestoreWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr);
620     ierr = VecRestoreArrayRead(Grad, &grad);CHKERRQ(ierr);
621   }
622   ierr = VecRestoreArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
623   ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
624   PetscFunctionReturn(0);
625 }
626 
627 #undef __FUNCT__
628 #define __FUNCT__ "DMPlexInsertBoundaryValues"
629 PetscErrorCode DMPlexInsertBoundaryValues(DM dm, Vec locX, PetscReal time, Vec faceGeomFVM, Vec cellGeomFVM, Vec gradFVM)
630 {
631   PetscInt       numBd, b;
632   PetscErrorCode ierr;
633 
634   PetscFunctionBegin;
635   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
636   PetscValidHeaderSpecific(locX, VEC_CLASSID, 2);
637   if (faceGeomFVM) {PetscValidHeaderSpecific(faceGeomFVM, VEC_CLASSID, 4);}
638   if (cellGeomFVM) {PetscValidHeaderSpecific(cellGeomFVM, VEC_CLASSID, 5);}
639   if (gradFVM)     {PetscValidHeaderSpecific(gradFVM, VEC_CLASSID, 6);}
640   ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
641   for (b = 0; b < numBd; ++b) {
642     PetscBool       isEssential;
643     const char     *labelname;
644     DMLabel         label;
645     PetscInt        field;
646     PetscObject     obj;
647     PetscClassId    id;
648     void          (*func)();
649     PetscInt        numids;
650     const PetscInt *ids;
651     void           *ctx;
652 
653     ierr = DMPlexGetBoundary(dm, b, &isEssential, NULL, &labelname, &field, &func, &numids, &ids, &ctx);CHKERRQ(ierr);
654     if (!isEssential) continue;
655     ierr = DMPlexGetLabel(dm, labelname, &label);CHKERRQ(ierr);
656     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
657     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
658     if (id == PETSCFE_CLASSID) {
659       ierr = DMPlexInsertBoundaryValues_FEM_Internal(dm, field, label, numids, ids, (void (*)(const PetscReal[], PetscScalar *, void *)) func, ctx, locX);CHKERRQ(ierr);
660     } else if (id == PETSCFV_CLASSID) {
661       ierr = DMPlexInsertBoundaryValues_FVM_Internal(dm, time, faceGeomFVM, cellGeomFVM, gradFVM,
662                                                      field, label, numids, ids, (PetscErrorCode (*)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*)) func, ctx, locX);CHKERRQ(ierr);
663     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
664   }
665   PetscFunctionReturn(0);
666 }
667 
668 #undef __FUNCT__
669 #define __FUNCT__ "DMPlexComputeL2Diff"
670 /*@C
671   DMPlexComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.
672 
673   Input Parameters:
674 + dm    - The DM
675 . funcs - The functions to evaluate for each field component
676 . ctxs  - Optional array of contexts to pass to each function, or NULL.
677 - X     - The coefficient vector u_h
678 
679   Output Parameter:
680 . diff - The diff ||u - u_h||_2
681 
682   Level: developer
683 
684 .seealso: DMPlexProjectFunction(), DMPlexComputeL2FieldDiff(), DMPlexComputeL2GradientDiff()
685 @*/
686 PetscErrorCode DMPlexComputeL2Diff(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
687 {
688   const PetscInt   debug = 0;
689   PetscSection     section;
690   PetscQuadrature  quad;
691   Vec              localX;
692   PetscScalar     *funcVal, *interpolant;
693   PetscReal       *coords, *v0, *J, *invJ, detJ;
694   PetscReal        localDiff = 0.0;
695   const PetscReal *quadPoints, *quadWeights;
696   PetscInt         dim, numFields, numComponents = 0, numQuadPoints, cStart, cEnd, c, field, fieldOffset;
697   PetscErrorCode   ierr;
698 
699   PetscFunctionBegin;
700   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
701   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
702   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
703   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
704   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
705   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
706   for (field = 0; field < numFields; ++field) {
707     PetscObject  obj;
708     PetscClassId id;
709     PetscInt     Nc;
710 
711     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
712     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
713     if (id == PETSCFE_CLASSID) {
714       PetscFE fe = (PetscFE) obj;
715 
716       ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
717       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
718     } else if (id == PETSCFV_CLASSID) {
719       PetscFV fv = (PetscFV) obj;
720 
721       ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr);
722       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
723     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
724     numComponents += Nc;
725   }
726   ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
727   ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
728   ierr = PetscMalloc6(numComponents,&funcVal,numComponents,&interpolant,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
729   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
730   for (c = cStart; c < cEnd; ++c) {
731     PetscScalar *x = NULL;
732     PetscReal    elemDiff = 0.0;
733 
734     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
735     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
736     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
737 
738     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
739       PetscObject  obj;
740       PetscClassId id;
741       void * const ctx = ctxs ? ctxs[field] : NULL;
742       PetscInt     Nb, Nc, q, fc;
743 
744       ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
745       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
746       if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);}
747       else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;}
748       else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
749       if (debug) {
750         char title[1024];
751         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
752         ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr);
753       }
754       for (q = 0; q < numQuadPoints; ++q) {
755         CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], coords);
756         (*funcs[field])(coords, funcVal, ctx);
757         if (id == PETSCFE_CLASSID)      {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
758         else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
759         else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
760         for (fc = 0; fc < Nc; ++fc) {
761           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);}
762           elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ;
763         }
764       }
765       fieldOffset += Nb*Nc;
766     }
767     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
768     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
769     localDiff += elemDiff;
770   }
771   ierr  = PetscFree6(funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr);
772   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
773   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
774   *diff = PetscSqrtReal(*diff);
775   PetscFunctionReturn(0);
776 }
777 
778 #undef __FUNCT__
779 #define __FUNCT__ "DMPlexComputeL2GradientDiff"
780 /*@C
781   DMPlexComputeL2GradientDiff - This function computes the L_2 difference between the gradient of a function u and an FEM interpolant solution grad u_h.
782 
783   Input Parameters:
784 + dm    - The DM
785 . funcs - The gradient functions to evaluate for each field component
786 . ctxs  - Optional array of contexts to pass to each function, or NULL.
787 . X     - The coefficient vector u_h
788 - n     - The vector to project along
789 
790   Output Parameter:
791 . diff - The diff ||(grad u - grad u_h) . n||_2
792 
793   Level: developer
794 
795 .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff()
796 @*/
797 PetscErrorCode DMPlexComputeL2GradientDiff(DM dm, void (**funcs)(const PetscReal [], const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, const PetscReal n[], PetscReal *diff)
798 {
799   const PetscInt  debug = 0;
800   PetscSection    section;
801   PetscQuadrature quad;
802   Vec             localX;
803   PetscScalar    *funcVal, *interpolantVec;
804   PetscReal      *coords, *realSpaceDer, *v0, *J, *invJ, detJ;
805   PetscReal       localDiff = 0.0;
806   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp;
807   PetscErrorCode  ierr;
808 
809   PetscFunctionBegin;
810   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
811   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
812   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
813   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
814   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
815   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
816   for (field = 0; field < numFields; ++field) {
817     PetscFE  fe;
818     PetscInt Nc;
819 
820     ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
821     ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
822     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
823     numComponents += Nc;
824   }
825   /* ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */
826   ierr = PetscMalloc7(numComponents,&funcVal,dim,&coords,dim,&realSpaceDer,dim,&v0,dim*dim,&J,dim*dim,&invJ,dim,&interpolantVec);CHKERRQ(ierr);
827   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
828   for (c = cStart; c < cEnd; ++c) {
829     PetscScalar *x = NULL;
830     PetscReal    elemDiff = 0.0;
831 
832     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
833     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
834     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
835 
836     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
837       PetscFE          fe;
838       void * const     ctx = ctxs ? ctxs[field] : NULL;
839       const PetscReal *quadPoints, *quadWeights;
840       PetscReal       *basisDer;
841       PetscInt         numQuadPoints, Nb, Ncomp, q, d, e, fc, f, g;
842 
843       ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
844       ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
845       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
846       ierr = PetscFEGetNumComponents(fe, &Ncomp);CHKERRQ(ierr);
847       ierr = PetscFEGetDefaultTabulation(fe, NULL, &basisDer, NULL);CHKERRQ(ierr);
848       if (debug) {
849         char title[1024];
850         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
851         ierr = DMPrintCellVector(c, title, Nb*Ncomp, &x[fieldOffset]);CHKERRQ(ierr);
852       }
853       for (q = 0; q < numQuadPoints; ++q) {
854         for (d = 0; d < dim; d++) {
855           coords[d] = v0[d];
856           for (e = 0; e < dim; e++) {
857             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
858           }
859         }
860         (*funcs[field])(coords, n, funcVal, ctx);
861         for (fc = 0; fc < Ncomp; ++fc) {
862           PetscScalar interpolant = 0.0;
863 
864           for (d = 0; d < dim; ++d) interpolantVec[d] = 0.0;
865           for (f = 0; f < Nb; ++f) {
866             const PetscInt fidx = f*Ncomp+fc;
867 
868             for (d = 0; d < dim; ++d) {
869               realSpaceDer[d] = 0.0;
870               for (g = 0; g < dim; ++g) {
871                 realSpaceDer[d] += invJ[g*dim+d]*basisDer[(q*Nb*Ncomp+fidx)*dim+g];
872               }
873               interpolantVec[d] += x[fieldOffset+fidx]*realSpaceDer[d];
874             }
875           }
876           for (d = 0; d < dim; ++d) interpolant += interpolantVec[d]*n[d];
877           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);}
878           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ;
879         }
880       }
881       comp        += Ncomp;
882       fieldOffset += Nb*Ncomp;
883     }
884     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
885     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
886     localDiff += elemDiff;
887   }
888   ierr  = PetscFree7(funcVal,coords,realSpaceDer,v0,J,invJ,interpolantVec);CHKERRQ(ierr);
889   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
890   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
891   *diff = PetscSqrtReal(*diff);
892   PetscFunctionReturn(0);
893 }
894 
895 #undef __FUNCT__
896 #define __FUNCT__ "DMPlexComputeL2FieldDiff"
897 /*@C
898   DMPlexComputeL2FieldDiff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h, separated into field components.
899 
900   Input Parameters:
901 + dm    - The DM
902 . funcs - The functions to evaluate for each field component
903 . ctxs  - Optional array of contexts to pass to each function, or NULL.
904 - X     - The coefficient vector u_h
905 
906   Output Parameter:
907 . diff - The array of differences, ||u^f - u^f_h||_2
908 
909   Level: developer
910 
911 .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff(), DMPlexComputeL2GradientDiff()
912 @*/
913 PetscErrorCode DMPlexComputeL2FieldDiff(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal diff[])
914 {
915   const PetscInt   debug = 0;
916   PetscSection     section;
917   PetscQuadrature  quad;
918   Vec              localX;
919   PetscScalar     *funcVal, *interpolant;
920   PetscReal       *coords, *v0, *J, *invJ, detJ;
921   PetscReal       *localDiff;
922   const PetscReal *quadPoints, *quadWeights;
923   PetscInt         dim, numFields, numComponents = 0, numQuadPoints, cStart, cEnd, c, field, fieldOffset;
924   PetscErrorCode   ierr;
925 
926   PetscFunctionBegin;
927   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
928   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
929   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
930   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
931   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
932   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
933   for (field = 0; field < numFields; ++field) {
934     PetscObject  obj;
935     PetscClassId id;
936     PetscInt     Nc;
937 
938     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
939     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
940     if (id == PETSCFE_CLASSID) {
941       PetscFE fe = (PetscFE) obj;
942 
943       ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
944       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
945     } else if (id == PETSCFV_CLASSID) {
946       PetscFV fv = (PetscFV) obj;
947 
948       ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr);
949       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
950     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
951     numComponents += Nc;
952   }
953   ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
954   ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
955   ierr = PetscCalloc7(numFields,&localDiff,numComponents,&funcVal,numComponents,&interpolant,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
956   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
957   for (c = cStart; c < cEnd; ++c) {
958     PetscScalar *x = NULL;
959 
960     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
961     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
962     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
963 
964     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
965       PetscObject  obj;
966       PetscClassId id;
967       void * const ctx = ctxs ? ctxs[field] : NULL;
968       PetscInt     Nb, Nc, q, fc;
969 
970       PetscReal       elemDiff = 0.0;
971 
972       ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
973       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
974       if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);}
975       else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;}
976       else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
977       if (debug) {
978         char title[1024];
979         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
980         ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr);
981       }
982       for (q = 0; q < numQuadPoints; ++q) {
983         CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], coords);
984         (*funcs[field])(coords, funcVal, ctx);
985         if (id == PETSCFE_CLASSID)      {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
986         else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
987         else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
988         for (fc = 0; fc < Nc; ++fc) {
989           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);}
990           elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ;
991         }
992       }
993       fieldOffset += Nb*Nc;
994       localDiff[field] += elemDiff;
995     }
996     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
997   }
998   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
999   ierr = MPI_Allreduce(localDiff, diff, numFields, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
1000   for (field = 0; field < numFields; ++field) diff[field] = PetscSqrtReal(diff[field]);
1001   ierr = PetscFree7(localDiff,funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr);
1002   PetscFunctionReturn(0);
1003 }
1004 
1005 #undef __FUNCT__
1006 #define __FUNCT__ "DMPlexComputeIntegralFEM"
1007 /*@
1008   DMPlexComputeIntegralFEM - Form the local integral F from the local input X using pointwise functions specified by the user
1009 
1010   Input Parameters:
1011 + dm - The mesh
1012 . X  - Local input vector
1013 - user - The user context
1014 
1015   Output Parameter:
1016 . integral - Local integral for each field
1017 
1018   Level: developer
1019 
1020 .seealso: DMPlexComputeResidualFEM()
1021 @*/
1022 PetscErrorCode DMPlexComputeIntegralFEM(DM dm, Vec X, PetscReal *integral, void *user)
1023 {
1024   DM_Plex          *mesh  = (DM_Plex *) dm->data;
1025   DM                dmAux;
1026   Vec               localX, A;
1027   PetscDS           prob, probAux = NULL;
1028   PetscQuadrature   q;
1029   PetscSection      section, sectionAux;
1030   PetscFECellGeom  *cgeom;
1031   PetscScalar      *u, *a = NULL;
1032   PetscInt          dim, Nf, f, numCells, cStart, cEnd, c;
1033   PetscInt          totDim, totDimAux;
1034   PetscErrorCode    ierr;
1035 
1036   PetscFunctionBegin;
1037   /*ierr = PetscLogEventBegin(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);*/
1038   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
1039   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
1040   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
1041   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1042   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1043   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
1044   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
1045   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
1046   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1047   numCells = cEnd - cStart;
1048   for (f = 0; f < Nf; ++f) {integral[f]    = 0.0;}
1049   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
1050   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
1051   if (dmAux) {
1052     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
1053     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
1054     ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
1055   }
1056   ierr = DMPlexInsertBoundaryValues(dm, localX, 0.0, NULL, NULL, NULL);CHKERRQ(ierr);
1057   ierr = PetscMalloc2(numCells*totDim,&u,numCells,&cgeom);CHKERRQ(ierr);
1058   if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);}
1059   for (c = cStart; c < cEnd; ++c) {
1060     PetscScalar *x = NULL;
1061     PetscInt     i;
1062 
1063     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, cgeom[c].v0, cgeom[c].J, cgeom[c].invJ, &cgeom[c].detJ);CHKERRQ(ierr);
1064     if (cgeom[c].detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", cgeom[c].detJ, c);
1065     ierr = DMPlexVecGetClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr);
1066     for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i];
1067     ierr = DMPlexVecRestoreClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr);
1068     if (dmAux) {
1069       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
1070       for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i];
1071       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
1072     }
1073   }
1074   for (f = 0; f < Nf; ++f) {
1075     PetscFE  fe;
1076     PetscInt numQuadPoints, Nb;
1077     /* Conforming batches */
1078     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1079     /* Remainder */
1080     PetscInt Nr, offset;
1081 
1082     ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
1083     ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
1084     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
1085     ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
1086     ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
1087     blockSize = Nb*numQuadPoints;
1088     batchSize = numBlocks * blockSize;
1089     ierr =  PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1090     numChunks = numCells / (numBatches*batchSize);
1091     Ne        = numChunks*numBatches*batchSize;
1092     Nr        = numCells % (numBatches*batchSize);
1093     offset    = numCells - Nr;
1094     ierr = PetscFEIntegrate(fe, prob, f, Ne, cgeom, u, probAux, a, integral);CHKERRQ(ierr);
1095     ierr = PetscFEIntegrate(fe, prob, f, Nr, &cgeom[offset], &u[offset*totDim], probAux, &a[offset*totDimAux], integral);CHKERRQ(ierr);
1096   }
1097   ierr = PetscFree2(u,cgeom);CHKERRQ(ierr);
1098   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
1099   if (mesh->printFEM) {
1100     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Local integral:");CHKERRQ(ierr);
1101     for (f = 0; f < Nf; ++f) {ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), " %g", integral[f]);CHKERRQ(ierr);}
1102     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "\n");CHKERRQ(ierr);
1103   }
1104   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
1105   /* TODO: Allreduce for integral */
1106   /*ierr = PetscLogEventEnd(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);*/
1107   PetscFunctionReturn(0);
1108 }
1109 
1110 #undef __FUNCT__
1111 #define __FUNCT__ "DMPlexComputeInterpolatorFEM"
1112 /*@
1113   DMPlexComputeInterpolatorFEM - Form the local portion of the interpolation matrix I from the coarse DM to the uniformly refined DM.
1114 
1115   Input Parameters:
1116 + dmf  - The fine mesh
1117 . dmc  - The coarse mesh
1118 - user - The user context
1119 
1120   Output Parameter:
1121 . In  - The interpolation matrix
1122 
1123   Note:
1124   The first member of the user context must be an FEMContext.
1125 
1126   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
1127   like a GPU, or vectorize on a multicore machine.
1128 
1129   Level: developer
1130 
1131 .seealso: DMPlexComputeJacobianFEM()
1132 @*/
1133 PetscErrorCode DMPlexComputeInterpolatorFEM(DM dmc, DM dmf, Mat In, void *user)
1134 {
1135   DM_Plex          *mesh  = (DM_Plex *) dmc->data;
1136   const char       *name  = "Interpolator";
1137   PetscDS           prob;
1138   PetscFE          *feRef;
1139   PetscSection      fsection, fglobalSection;
1140   PetscSection      csection, cglobalSection;
1141   PetscScalar      *elemMat;
1142   PetscInt          dim, Nf, f, fieldI, fieldJ, offsetI, offsetJ, cStart, cEnd, c;
1143   PetscInt          cTotDim, rTotDim = 0;
1144   PetscErrorCode    ierr;
1145 
1146   PetscFunctionBegin;
1147   ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1148   ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr);
1149   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
1150   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
1151   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
1152   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
1153   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
1154   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
1155   ierr = DMGetDS(dmf, &prob);CHKERRQ(ierr);
1156   ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr);
1157   for (f = 0; f < Nf; ++f) {
1158     PetscFE  fe;
1159     PetscInt rNb, Nc;
1160 
1161     ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
1162     ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr);
1163     ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr);
1164     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1165     rTotDim += rNb*Nc;
1166   }
1167   ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr);
1168   ierr = PetscMalloc1(rTotDim*cTotDim,&elemMat);CHKERRQ(ierr);
1169   ierr = PetscMemzero(elemMat, rTotDim*cTotDim * sizeof(PetscScalar));CHKERRQ(ierr);
1170   for (fieldI = 0, offsetI = 0; fieldI < Nf; ++fieldI) {
1171     PetscDualSpace   Qref;
1172     PetscQuadrature  f;
1173     const PetscReal *qpoints, *qweights;
1174     PetscReal       *points;
1175     PetscInt         npoints = 0, Nc, Np, fpdim, i, k, p, d;
1176 
1177     /* Compose points from all dual basis functionals */
1178     ierr = PetscFEGetDualSpace(feRef[fieldI], &Qref);CHKERRQ(ierr);
1179     ierr = PetscFEGetNumComponents(feRef[fieldI], &Nc);CHKERRQ(ierr);
1180     ierr = PetscDualSpaceGetDimension(Qref, &fpdim);CHKERRQ(ierr);
1181     for (i = 0; i < fpdim; ++i) {
1182       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1183       ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, NULL);CHKERRQ(ierr);
1184       npoints += Np;
1185     }
1186     ierr = PetscMalloc1(npoints*dim,&points);CHKERRQ(ierr);
1187     for (i = 0, k = 0; i < fpdim; ++i) {
1188       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1189       ierr = PetscQuadratureGetData(f, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr);
1190       for (p = 0; p < Np; ++p, ++k) for (d = 0; d < dim; ++d) points[k*dim+d] = qpoints[p*dim+d];
1191     }
1192 
1193     for (fieldJ = 0, offsetJ = 0; fieldJ < Nf; ++fieldJ) {
1194       PetscFE    fe;
1195       PetscReal *B;
1196       PetscInt   NcJ, cpdim, j;
1197 
1198       /* Evaluate basis at points */
1199       ierr = PetscDSGetDiscretization(prob, fieldJ, (PetscObject *) &fe);CHKERRQ(ierr);
1200       ierr = PetscFEGetNumComponents(fe, &NcJ);CHKERRQ(ierr);
1201       ierr = PetscFEGetDimension(fe, &cpdim);CHKERRQ(ierr);
1202       /* For now, fields only interpolate themselves */
1203       if (fieldI == fieldJ) {
1204         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);
1205         ierr = PetscFEGetTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);
1206         for (i = 0, k = 0; i < fpdim; ++i) {
1207           ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1208           ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, &qweights);CHKERRQ(ierr);
1209           for (p = 0; p < Np; ++p, ++k) {
1210             for (j = 0; j < cpdim; ++j) {
1211               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];
1212             }
1213           }
1214         }
1215         ierr = PetscFERestoreTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr);
1216       }
1217       offsetJ += cpdim*NcJ;
1218     }
1219     offsetI += fpdim*Nc;
1220     ierr = PetscFree(points);CHKERRQ(ierr);
1221   }
1222   if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(0, name, rTotDim, cTotDim, elemMat);CHKERRQ(ierr);}
1223   /* Preallocate matrix */
1224   {
1225     PetscHashJK ht;
1226     PetscLayout rLayout;
1227     PetscInt   *dnz, *onz, *cellCIndices, *cellFIndices;
1228     PetscInt    locRows, rStart, rEnd, cell, r;
1229 
1230     ierr = MatGetLocalSize(In, &locRows, NULL);CHKERRQ(ierr);
1231     ierr = PetscLayoutCreate(PetscObjectComm((PetscObject) In), &rLayout);CHKERRQ(ierr);
1232     ierr = PetscLayoutSetLocalSize(rLayout, locRows);CHKERRQ(ierr);
1233     ierr = PetscLayoutSetBlockSize(rLayout, 1);CHKERRQ(ierr);
1234     ierr = PetscLayoutSetUp(rLayout);CHKERRQ(ierr);
1235     ierr = PetscLayoutGetRange(rLayout, &rStart, &rEnd);CHKERRQ(ierr);
1236     ierr = PetscLayoutDestroy(&rLayout);CHKERRQ(ierr);
1237     ierr = PetscCalloc4(locRows,&dnz,locRows,&onz,cTotDim,&cellCIndices,rTotDim,&cellFIndices);CHKERRQ(ierr);
1238     ierr = PetscHashJKCreate(&ht);CHKERRQ(ierr);
1239     for (cell = cStart; cell < cEnd; ++cell) {
1240       ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, cell, cellCIndices, cellFIndices);CHKERRQ(ierr);
1241       for (r = 0; r < rTotDim; ++r) {
1242         PetscHashJKKey  key;
1243         PetscHashJKIter missing, iter;
1244 
1245         key.j = cellFIndices[r];
1246         if (key.j < 0) continue;
1247         for (c = 0; c < cTotDim; ++c) {
1248           key.k = cellCIndices[c];
1249           if (key.k < 0) continue;
1250           ierr = PetscHashJKPut(ht, key, &missing, &iter);CHKERRQ(ierr);
1251           if (missing) {
1252             ierr = PetscHashJKSet(ht, iter, 1);CHKERRQ(ierr);
1253             if ((key.k >= rStart) && (key.k < rEnd)) ++dnz[key.j-rStart];
1254             else                                     ++onz[key.j-rStart];
1255           }
1256         }
1257       }
1258     }
1259     ierr = PetscHashJKDestroy(&ht);CHKERRQ(ierr);
1260     ierr = MatXAIJSetPreallocation(In, 1, dnz, onz, NULL, NULL);CHKERRQ(ierr);
1261     ierr = MatSetOption(In, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
1262     ierr = PetscFree4(dnz,onz,cellCIndices,cellFIndices);CHKERRQ(ierr);
1263   }
1264   /* Fill matrix */
1265   ierr = MatZeroEntries(In);CHKERRQ(ierr);
1266   for (c = cStart; c < cEnd; ++c) {
1267     ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES);CHKERRQ(ierr);
1268   }
1269   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);}
1270   ierr = PetscFree(feRef);CHKERRQ(ierr);
1271   ierr = PetscFree(elemMat);CHKERRQ(ierr);
1272   ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1273   ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1274   if (mesh->printFEM) {
1275     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
1276     ierr = MatChop(In, 1.0e-10);CHKERRQ(ierr);
1277     ierr = MatView(In, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1278   }
1279   ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1280   PetscFunctionReturn(0);
1281 }
1282 
1283 #undef __FUNCT__
1284 #define __FUNCT__ "DMPlexComputeInjectorFEM"
1285 PetscErrorCode DMPlexComputeInjectorFEM(DM dmc, DM dmf, VecScatter *sc, void *user)
1286 {
1287   PetscDS        prob;
1288   PetscFE       *feRef;
1289   Vec            fv, cv;
1290   IS             fis, cis;
1291   PetscSection   fsection, fglobalSection, csection, cglobalSection;
1292   PetscInt      *cmap, *cellCIndices, *cellFIndices, *cindices, *findices;
1293   PetscInt       cTotDim, fTotDim = 0, Nf, f, field, cStart, cEnd, c, dim, d, startC, offsetC, offsetF, m;
1294   PetscErrorCode ierr;
1295 
1296   PetscFunctionBegin;
1297   ierr = PetscLogEventBegin(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1298   ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr);
1299   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
1300   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
1301   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
1302   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
1303   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
1304   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
1305   ierr = DMGetDS(dmc, &prob);CHKERRQ(ierr);
1306   ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr);
1307   for (f = 0; f < Nf; ++f) {
1308     PetscFE  fe;
1309     PetscInt fNb, Nc;
1310 
1311     ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
1312     ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr);
1313     ierr = PetscFEGetDimension(feRef[f], &fNb);CHKERRQ(ierr);
1314     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1315     fTotDim += fNb*Nc;
1316   }
1317   ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr);
1318   ierr = PetscMalloc1(cTotDim,&cmap);CHKERRQ(ierr);
1319   for (field = 0, offsetC = 0, offsetF = 0; field < Nf; ++field) {
1320     PetscFE        feC;
1321     PetscDualSpace QF, QC;
1322     PetscInt       NcF, NcC, fpdim, cpdim;
1323 
1324     ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &feC);CHKERRQ(ierr);
1325     ierr = PetscFEGetNumComponents(feC, &NcC);CHKERRQ(ierr);
1326     ierr = PetscFEGetNumComponents(feRef[field], &NcF);CHKERRQ(ierr);
1327     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);
1328     ierr = PetscFEGetDualSpace(feRef[field], &QF);CHKERRQ(ierr);
1329     ierr = PetscDualSpaceGetDimension(QF, &fpdim);CHKERRQ(ierr);
1330     ierr = PetscFEGetDualSpace(feC, &QC);CHKERRQ(ierr);
1331     ierr = PetscDualSpaceGetDimension(QC, &cpdim);CHKERRQ(ierr);
1332     for (c = 0; c < cpdim; ++c) {
1333       PetscQuadrature  cfunc;
1334       const PetscReal *cqpoints;
1335       PetscInt         NpC;
1336 
1337       ierr = PetscDualSpaceGetFunctional(QC, c, &cfunc);CHKERRQ(ierr);
1338       ierr = PetscQuadratureGetData(cfunc, NULL, &NpC, &cqpoints, NULL);CHKERRQ(ierr);
1339       if (NpC != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not know how to do injection for moments");
1340       for (f = 0; f < fpdim; ++f) {
1341         PetscQuadrature  ffunc;
1342         const PetscReal *fqpoints;
1343         PetscReal        sum = 0.0;
1344         PetscInt         NpF, comp;
1345 
1346         ierr = PetscDualSpaceGetFunctional(QF, f, &ffunc);CHKERRQ(ierr);
1347         ierr = PetscQuadratureGetData(ffunc, NULL, &NpF, &fqpoints, NULL);CHKERRQ(ierr);
1348         if (NpC != NpF) continue;
1349         for (d = 0; d < dim; ++d) sum += PetscAbsReal(cqpoints[d] - fqpoints[d]);
1350         if (sum > 1.0e-9) continue;
1351         for (comp = 0; comp < NcC; ++comp) {
1352           cmap[(offsetC+c)*NcC+comp] = (offsetF+f)*NcF+comp;
1353         }
1354         break;
1355       }
1356     }
1357     offsetC += cpdim*NcC;
1358     offsetF += fpdim*NcF;
1359   }
1360   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);}
1361   ierr = PetscFree(feRef);CHKERRQ(ierr);
1362 
1363   ierr = DMGetGlobalVector(dmf, &fv);CHKERRQ(ierr);
1364   ierr = DMGetGlobalVector(dmc, &cv);CHKERRQ(ierr);
1365   ierr = VecGetOwnershipRange(cv, &startC, NULL);CHKERRQ(ierr);
1366   ierr = PetscSectionGetConstrainedStorageSize(cglobalSection, &m);CHKERRQ(ierr);
1367   ierr = PetscMalloc2(cTotDim,&cellCIndices,fTotDim,&cellFIndices);CHKERRQ(ierr);
1368   ierr = PetscMalloc1(m,&cindices);CHKERRQ(ierr);
1369   ierr = PetscMalloc1(m,&findices);CHKERRQ(ierr);
1370   for (d = 0; d < m; ++d) cindices[d] = findices[d] = -1;
1371   for (c = cStart; c < cEnd; ++c) {
1372     ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, c, cellCIndices, cellFIndices);CHKERRQ(ierr);
1373     for (d = 0; d < cTotDim; ++d) {
1374       if (cellCIndices[d] < 0) continue;
1375       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]]);
1376       cindices[cellCIndices[d]-startC] = cellCIndices[d];
1377       findices[cellCIndices[d]-startC] = cellFIndices[cmap[d]];
1378     }
1379   }
1380   ierr = PetscFree(cmap);CHKERRQ(ierr);
1381   ierr = PetscFree2(cellCIndices,cellFIndices);CHKERRQ(ierr);
1382 
1383   ierr = ISCreateGeneral(PETSC_COMM_SELF, m, cindices, PETSC_OWN_POINTER, &cis);CHKERRQ(ierr);
1384   ierr = ISCreateGeneral(PETSC_COMM_SELF, m, findices, PETSC_OWN_POINTER, &fis);CHKERRQ(ierr);
1385   ierr = VecScatterCreate(cv, cis, fv, fis, sc);CHKERRQ(ierr);
1386   ierr = ISDestroy(&cis);CHKERRQ(ierr);
1387   ierr = ISDestroy(&fis);CHKERRQ(ierr);
1388   ierr = DMRestoreGlobalVector(dmf, &fv);CHKERRQ(ierr);
1389   ierr = DMRestoreGlobalVector(dmc, &cv);CHKERRQ(ierr);
1390   ierr = PetscLogEventEnd(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1391   PetscFunctionReturn(0);
1392 }
1393