xref: /petsc/src/dm/impls/plex/plexfem.c (revision bbce034c4e250f2d9613c607136e917b6de06000)
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__ "DMPlexProjectFunctionLabelLocal"
186 PetscErrorCode DMPlexProjectFunctionLabelLocal(DM dm, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
187 {
188   PetscDualSpace *sp;
189   PetscSection    section;
190   PetscScalar    *values;
191   PetscBool      *fieldActive;
192   PetscInt        numFields, numComp, dim, spDim, totDim = 0, numValues, cStart, cEnd, f, d, v, i, comp;
193   PetscErrorCode  ierr;
194 
195   PetscFunctionBegin;
196   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
197   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
198   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
199   ierr = PetscMalloc1(numFields,&sp);CHKERRQ(ierr);
200   for (f = 0; f < numFields; ++f) {
201     ierr = PetscFEGetDualSpace(fe[f], &sp[f]);CHKERRQ(ierr);
202     ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr);
203     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
204     totDim += spDim*numComp;
205   }
206   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
207   ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr);
208   if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section cell closure size %d != dual space dimension %d", numValues, totDim);
209   ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
210   ierr = DMGetWorkArray(dm, numFields, PETSC_BOOL, &fieldActive);CHKERRQ(ierr);
211   for (f = 0; f < numFields; ++f) fieldActive[f] = funcs[f] ? PETSC_TRUE : PETSC_FALSE;
212   for (i = 0; i < numIds; ++i) {
213     IS              pointIS;
214     const PetscInt *points;
215     PetscInt        n, p;
216 
217     ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr);
218     ierr = ISGetLocalSize(pointIS, &n);CHKERRQ(ierr);
219     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
220     for (p = 0; p < n; ++p) {
221       const PetscInt  point = points[p];
222       PetscFECellGeom geom;
223 
224       if ((point < cStart) || (point >= cEnd)) continue;
225       ierr = DMPlexComputeCellGeometryFEM(dm, point, NULL, geom.v0, geom.J, NULL, &geom.detJ);CHKERRQ(ierr);
226       for (f = 0, v = 0; f < numFields; ++f) {
227         void * const ctx = ctxs ? ctxs[f] : NULL;
228         ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr);
229         ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
230         for (d = 0; d < spDim; ++d) {
231           if (funcs[f]) {
232             ierr = PetscDualSpaceApply(sp[f], d, &geom, numComp, funcs[f], ctx, &values[v]);CHKERRQ(ierr);
233           } else {
234             for (comp = 0; comp < numComp; ++comp) values[v+comp] = 0.0;
235           }
236           v += numComp;
237         }
238       }
239       ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, point, values, mode);CHKERRQ(ierr);
240     }
241     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
242     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
243   }
244   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
245   ierr = DMRestoreWorkArray(dm, numFields, PETSC_BOOL, &fieldActive);CHKERRQ(ierr);
246   ierr = PetscFree(sp);CHKERRQ(ierr);
247   PetscFunctionReturn(0);
248 }
249 
250 #undef __FUNCT__
251 #define __FUNCT__ "DMPlexProjectFunctionLocal"
252 PetscErrorCode DMPlexProjectFunctionLocal(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
253 {
254   PetscDualSpace *sp;
255   PetscInt       *numComp;
256   PetscSection    section;
257   PetscScalar    *values;
258   PetscInt        numFields, dim, spDim, totDim = 0, numValues, cStart, cEnd, c, f, d, v, comp;
259   PetscErrorCode  ierr;
260 
261   PetscFunctionBegin;
262   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
263   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
264   ierr = PetscMalloc2(numFields, &sp, numFields, &numComp);CHKERRQ(ierr);
265   for (f = 0; f < numFields; ++f) {
266     PetscObject  obj;
267     PetscClassId id;
268 
269     ierr = DMGetField(dm, f, &obj);CHKERRQ(ierr);
270     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
271     if (id == PETSCFE_CLASSID) {
272       PetscFE fe = (PetscFE) obj;
273 
274       ierr = PetscFEGetNumComponents(fe, &numComp[f]);CHKERRQ(ierr);
275       ierr = PetscFEGetDualSpace(fe, &sp[f]);CHKERRQ(ierr);
276       ierr = PetscObjectReference((PetscObject) sp[f]);CHKERRQ(ierr);
277     } else if (id == PETSCFV_CLASSID) {
278       PetscFV         fv = (PetscFV) obj;
279       PetscQuadrature q;
280 
281       ierr = PetscFVGetNumComponents(fv, &numComp[f]);CHKERRQ(ierr);
282       ierr = PetscFVGetQuadrature(fv, &q);CHKERRQ(ierr);
283       ierr = PetscDualSpaceCreate(PetscObjectComm((PetscObject) fv), &sp[f]);CHKERRQ(ierr);
284       ierr = PetscDualSpaceSetDM(sp[f], dm);CHKERRQ(ierr);
285       ierr = PetscDualSpaceSetType(sp[f], PETSCDUALSPACESIMPLE);CHKERRQ(ierr);
286       ierr = PetscDualSpaceSimpleSetDimension(sp[f], 1);CHKERRQ(ierr);
287       ierr = PetscDualSpaceSimpleSetFunctional(sp[f], 0, q);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 = DMGetDimension(dm, &dim);CHKERRQ(ierr);
293   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
294   ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr);
295   if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section cell closure size %d != dual space dimension %d", numValues, totDim);
296   ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
297   for (c = cStart; c < cEnd; ++c) {
298     PetscFECellGeom geom;
299 
300     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, geom.v0, geom.J, NULL, &geom.detJ);CHKERRQ(ierr);
301     for (f = 0, v = 0; f < numFields; ++f) {
302       void *const ctx = ctxs ? ctxs[f] : NULL;
303 
304       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
305       for (d = 0; d < spDim; ++d) {
306         if (funcs[f]) {
307           ierr = PetscDualSpaceApply(sp[f], d, &geom, numComp[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr);
308         } else {
309           for (comp = 0; comp < numComp[f]; ++comp) values[v+comp] = 0.0;
310         }
311         v += numComp[f];
312       }
313     }
314     ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr);
315   }
316   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
317   for (f = 0; f < numFields; ++f) {ierr = PetscDualSpaceDestroy(&sp[f]);CHKERRQ(ierr);}
318   ierr = PetscFree2(sp, numComp);CHKERRQ(ierr);
319   PetscFunctionReturn(0);
320 }
321 
322 #undef __FUNCT__
323 #define __FUNCT__ "DMPlexProjectFunction"
324 /*@C
325   DMPlexProjectFunction - This projects the given function into the function space provided.
326 
327   Input Parameters:
328 + dm      - The DM
329 . funcs   - The coordinate functions to evaluate, one per field
330 . ctxs    - Optional array of contexts to pass to each coordinate function.  ctxs itself may be null.
331 - mode    - The insertion mode for values
332 
333   Output Parameter:
334 . X - vector
335 
336   Level: developer
337 
338 .seealso: DMPlexComputeL2Diff()
339 @*/
340 PetscErrorCode DMPlexProjectFunction(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec X)
341 {
342   Vec            localX;
343   PetscErrorCode ierr;
344 
345   PetscFunctionBegin;
346   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
347   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
348   ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, mode, localX);CHKERRQ(ierr);
349   ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr);
350   ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr);
351   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
352   PetscFunctionReturn(0);
353 }
354 
355 #undef __FUNCT__
356 #define __FUNCT__ "DMPlexProjectFieldLocal"
357 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)
358 {
359   DM                dmAux;
360   PetscDS           prob, probAux;
361   Vec               A;
362   PetscSection      section, sectionAux;
363   PetscScalar      *values, *u, *u_x, *a, *a_x;
364   PetscReal        *x, *v0, *J, *invJ, detJ, **basisField, **basisFieldDer, **basisFieldAux, **basisFieldDerAux;
365   PetscInt          Nf, dim, spDim, totDim, numValues, cStart, cEnd, c, f, d, v, comp;
366   PetscErrorCode    ierr;
367 
368   PetscFunctionBegin;
369   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
370   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
371   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
372   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
373   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
374   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
375   ierr = PetscDSGetTabulation(prob, &basisField, &basisFieldDer);CHKERRQ(ierr);
376   ierr = PetscDSGetEvaluationArrays(prob, &u, NULL, &u_x);CHKERRQ(ierr);
377   ierr = PetscDSGetRefCoordArrays(prob, &x, NULL);CHKERRQ(ierr);
378   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
379   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
380   if (dmAux) {
381     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
382     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
383     ierr = PetscDSGetTabulation(prob, &basisFieldAux, &basisFieldDerAux);CHKERRQ(ierr);
384     ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);CHKERRQ(ierr);
385   }
386   ierr = DMPlexInsertBoundaryValuesFEM(dm, localU);CHKERRQ(ierr);
387   ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr);
388   if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section cell closure size %d != dual space dimension %d", numValues, totDim);
389   ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
390   ierr = PetscMalloc3(dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
391   for (c = cStart; c < cEnd; ++c) {
392     PetscScalar *coefficients = NULL, *coefficientsAux = NULL;
393 
394     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
395     ierr = DMPlexVecGetClosure(dm, section, localU, c, NULL, &coefficients);CHKERRQ(ierr);
396     if (dmAux) {ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &coefficientsAux);CHKERRQ(ierr);}
397     for (f = 0, v = 0; f < Nf; ++f) {
398       PetscFE        fe;
399       PetscDualSpace sp;
400       PetscInt       Ncf;
401 
402       ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
403       ierr = PetscFEGetDualSpace(fe, &sp);CHKERRQ(ierr);
404       ierr = PetscFEGetNumComponents(fe, &Ncf);CHKERRQ(ierr);
405       ierr = PetscDualSpaceGetDimension(sp, &spDim);CHKERRQ(ierr);
406       for (d = 0; d < spDim; ++d) {
407         PetscQuadrature  quad;
408         const PetscReal *points, *weights;
409         PetscInt         numPoints, q;
410 
411         if (funcs[f]) {
412           ierr = PetscDualSpaceGetFunctional(sp, d, &quad);CHKERRQ(ierr);
413           ierr = PetscQuadratureGetData(quad, NULL, &numPoints, &points, &weights);CHKERRQ(ierr);
414           ierr = PetscFEGetTabulation(fe, numPoints, points, &basisField[f], &basisFieldDer[f], NULL);CHKERRQ(ierr);
415           for (q = 0; q < numPoints; ++q) {
416             CoordinatesRefToReal(dim, dim, v0, J, &points[q*dim], x);
417             ierr = EvaluateFieldJets(prob,    PETSC_FALSE, q, invJ, coefficients,    NULL, u, u_x, NULL);CHKERRQ(ierr);
418             ierr = EvaluateFieldJets(probAux, PETSC_FALSE, q, invJ, coefficientsAux, NULL, a, a_x, NULL);CHKERRQ(ierr);
419             (*funcs[f])(u, NULL, u_x, a, NULL, a_x, x, &values[v]);
420           }
421           ierr = PetscFERestoreTabulation(fe, numPoints, points, &basisField[f], &basisFieldDer[f], NULL);CHKERRQ(ierr);
422         } else {
423           for (comp = 0; comp < Ncf; ++comp) values[v+comp] = 0.0;
424         }
425         v += Ncf;
426       }
427     }
428     ierr = DMPlexVecRestoreClosure(dm, section, localU, c, NULL, &coefficients);CHKERRQ(ierr);
429     if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &coefficientsAux);CHKERRQ(ierr);}
430     ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr);
431   }
432   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
433   ierr = PetscFree3(v0,J,invJ);CHKERRQ(ierr);
434   PetscFunctionReturn(0);
435 }
436 
437 #undef __FUNCT__
438 #define __FUNCT__ "DMPlexProjectField"
439 /*@C
440   DMPlexProjectField - This projects the given function of the fields into the function space provided.
441 
442   Input Parameters:
443 + dm      - The DM
444 . U       - The input field vector
445 . funcs   - The functions to evaluate, one per field
446 - mode    - The insertion mode for values
447 
448   Output Parameter:
449 . X       - The output vector
450 
451   Level: developer
452 
453 .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff()
454 @*/
455 PetscErrorCode DMPlexProjectField(DM dm, Vec U, void (**funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal [], PetscScalar []), InsertMode mode, Vec X)
456 {
457   Vec            localX, localU;
458   PetscErrorCode ierr;
459 
460   PetscFunctionBegin;
461   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
462   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
463   ierr = DMGetLocalVector(dm, &localU);CHKERRQ(ierr);
464   ierr = DMGlobalToLocalBegin(dm, U, INSERT_VALUES, localU);CHKERRQ(ierr);
465   ierr = DMGlobalToLocalEnd(dm, U, INSERT_VALUES, localU);CHKERRQ(ierr);
466   ierr = DMPlexProjectFieldLocal(dm, localU, funcs, mode, localX);CHKERRQ(ierr);
467   ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr);
468   ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr);
469   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
470   ierr = DMRestoreLocalVector(dm, &localU);CHKERRQ(ierr);
471   PetscFunctionReturn(0);
472 }
473 
474 #undef __FUNCT__
475 #define __FUNCT__ "DMPlexInsertBoundaryValuesFEM"
476 PetscErrorCode DMPlexInsertBoundaryValuesFEM(DM dm, Vec localX)
477 {
478   void        (**funcs)(const PetscReal x[], PetscScalar *u, void *ctx);
479   void         **ctxs;
480   PetscFE       *fe;
481   PetscInt       numFields, f, numBd, b;
482   PetscErrorCode ierr;
483 
484   PetscFunctionBegin;
485   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
486   PetscValidHeaderSpecific(localX, VEC_CLASSID, 2);
487   ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr);
488   ierr = PetscMalloc3(numFields,&fe,numFields,&funcs,numFields,&ctxs);CHKERRQ(ierr);
489   for (f = 0; f < numFields; ++f) {ierr = DMGetField(dm, f, (PetscObject *) &fe[f]);CHKERRQ(ierr);}
490   /* OPT: Could attempt to do multiple BCs at once */
491   ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
492   for (b = 0; b < numBd; ++b) {
493     DMLabel         label;
494     const PetscInt *ids;
495     const char     *labelname;
496     PetscInt        numids, field;
497     PetscBool       isEssential;
498     void          (*func)();
499     void           *ctx;
500 
501     ierr = DMPlexGetBoundary(dm, b, &isEssential, NULL, &labelname, &field, &func, &numids, &ids, &ctx);CHKERRQ(ierr);
502     ierr = DMPlexGetLabel(dm, labelname, &label);CHKERRQ(ierr);
503     for (f = 0; f < numFields; ++f) {
504       funcs[f] = field == f ? (void (*)(const PetscReal[], PetscScalar *, void *)) func : NULL;
505       ctxs[f]  = field == f ? ctx : NULL;
506     }
507     ierr = DMPlexProjectFunctionLabelLocal(dm, label, numids, ids, fe, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
508   }
509   ierr = PetscFree3(fe,funcs,ctxs);CHKERRQ(ierr);
510   PetscFunctionReturn(0);
511 }
512 
513 #undef __FUNCT__
514 #define __FUNCT__ "DMPlexComputeL2Diff"
515 /*@C
516   DMPlexComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.
517 
518   Input Parameters:
519 + dm    - The DM
520 . funcs - The functions to evaluate for each field component
521 . ctxs  - Optional array of contexts to pass to each function, or NULL.
522 - X     - The coefficient vector u_h
523 
524   Output Parameter:
525 . diff - The diff ||u - u_h||_2
526 
527   Level: developer
528 
529 .seealso: DMPlexProjectFunction(), DMPlexComputeL2FieldDiff(), DMPlexComputeL2GradientDiff()
530 @*/
531 PetscErrorCode DMPlexComputeL2Diff(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
532 {
533   const PetscInt   debug = 0;
534   PetscSection     section;
535   PetscQuadrature  quad;
536   Vec              localX;
537   PetscScalar     *funcVal, *interpolant;
538   PetscReal       *coords, *v0, *J, *invJ, detJ;
539   PetscReal        localDiff = 0.0;
540   const PetscReal *quadPoints, *quadWeights;
541   PetscInt         dim, numFields, numComponents = 0, numQuadPoints, cStart, cEnd, c, field, fieldOffset;
542   PetscErrorCode   ierr;
543 
544   PetscFunctionBegin;
545   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
546   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
547   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
548   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
549   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
550   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
551   for (field = 0; field < numFields; ++field) {
552     PetscObject  obj;
553     PetscClassId id;
554     PetscInt     Nc;
555 
556     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
557     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
558     if (id == PETSCFE_CLASSID) {
559       PetscFE fe = (PetscFE) obj;
560 
561       ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
562       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
563     } else if (id == PETSCFV_CLASSID) {
564       PetscFV fv = (PetscFV) obj;
565 
566       ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr);
567       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
568     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
569     numComponents += Nc;
570   }
571   ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
572   ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
573   ierr = PetscMalloc6(numComponents,&funcVal,numComponents,&interpolant,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
574   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
575   for (c = cStart; c < cEnd; ++c) {
576     PetscScalar *x = NULL;
577     PetscReal    elemDiff = 0.0;
578 
579     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
580     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
581     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
582 
583     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
584       PetscObject  obj;
585       PetscClassId id;
586       void * const ctx = ctxs ? ctxs[field] : NULL;
587       PetscInt     Nb, Nc, q, fc;
588 
589       ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
590       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
591       if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);}
592       else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;}
593       else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
594       if (debug) {
595         char title[1024];
596         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
597         ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr);
598       }
599       for (q = 0; q < numQuadPoints; ++q) {
600         CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], coords);
601         (*funcs[field])(coords, funcVal, ctx);
602         if (id == PETSCFE_CLASSID)      {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
603         else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
604         else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
605         for (fc = 0; fc < Nc; ++fc) {
606           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);}
607           elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ;
608         }
609       }
610       fieldOffset += Nb*Nc;
611     }
612     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
613     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
614     localDiff += elemDiff;
615   }
616   ierr  = PetscFree6(funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr);
617   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
618   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
619   *diff = PetscSqrtReal(*diff);
620   PetscFunctionReturn(0);
621 }
622 
623 #undef __FUNCT__
624 #define __FUNCT__ "DMPlexComputeL2GradientDiff"
625 /*@C
626   DMPlexComputeL2GradientDiff - This function computes the L_2 difference between the gradient of a function u and an FEM interpolant solution grad u_h.
627 
628   Input Parameters:
629 + dm    - The DM
630 . funcs - The gradient functions to evaluate for each field component
631 . ctxs  - Optional array of contexts to pass to each function, or NULL.
632 . X     - The coefficient vector u_h
633 - n     - The vector to project along
634 
635   Output Parameter:
636 . diff - The diff ||(grad u - grad u_h) . n||_2
637 
638   Level: developer
639 
640 .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff()
641 @*/
642 PetscErrorCode DMPlexComputeL2GradientDiff(DM dm, void (**funcs)(const PetscReal [], const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, const PetscReal n[], PetscReal *diff)
643 {
644   const PetscInt  debug = 0;
645   PetscSection    section;
646   PetscQuadrature quad;
647   Vec             localX;
648   PetscScalar    *funcVal, *interpolantVec;
649   PetscReal      *coords, *realSpaceDer, *v0, *J, *invJ, detJ;
650   PetscReal       localDiff = 0.0;
651   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp;
652   PetscErrorCode  ierr;
653 
654   PetscFunctionBegin;
655   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
656   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
657   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
658   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
659   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
660   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
661   for (field = 0; field < numFields; ++field) {
662     PetscFE  fe;
663     PetscInt Nc;
664 
665     ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
666     ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
667     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
668     numComponents += Nc;
669   }
670   /* ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */
671   ierr = PetscMalloc7(numComponents,&funcVal,dim,&coords,dim,&realSpaceDer,dim,&v0,dim*dim,&J,dim*dim,&invJ,dim,&interpolantVec);CHKERRQ(ierr);
672   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
673   for (c = cStart; c < cEnd; ++c) {
674     PetscScalar *x = NULL;
675     PetscReal    elemDiff = 0.0;
676 
677     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
678     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
679     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
680 
681     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
682       PetscFE          fe;
683       void * const     ctx = ctxs ? ctxs[field] : NULL;
684       const PetscReal *quadPoints, *quadWeights;
685       PetscReal       *basisDer;
686       PetscInt         numQuadPoints, Nb, Ncomp, q, d, e, fc, f, g;
687 
688       ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
689       ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
690       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
691       ierr = PetscFEGetNumComponents(fe, &Ncomp);CHKERRQ(ierr);
692       ierr = PetscFEGetDefaultTabulation(fe, NULL, &basisDer, NULL);CHKERRQ(ierr);
693       if (debug) {
694         char title[1024];
695         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
696         ierr = DMPrintCellVector(c, title, Nb*Ncomp, &x[fieldOffset]);CHKERRQ(ierr);
697       }
698       for (q = 0; q < numQuadPoints; ++q) {
699         for (d = 0; d < dim; d++) {
700           coords[d] = v0[d];
701           for (e = 0; e < dim; e++) {
702             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
703           }
704         }
705         (*funcs[field])(coords, n, funcVal, ctx);
706         for (fc = 0; fc < Ncomp; ++fc) {
707           PetscScalar interpolant = 0.0;
708 
709           for (d = 0; d < dim; ++d) interpolantVec[d] = 0.0;
710           for (f = 0; f < Nb; ++f) {
711             const PetscInt fidx = f*Ncomp+fc;
712 
713             for (d = 0; d < dim; ++d) {
714               realSpaceDer[d] = 0.0;
715               for (g = 0; g < dim; ++g) {
716                 realSpaceDer[d] += invJ[g*dim+d]*basisDer[(q*Nb*Ncomp+fidx)*dim+g];
717               }
718               interpolantVec[d] += x[fieldOffset+fidx]*realSpaceDer[d];
719             }
720           }
721           for (d = 0; d < dim; ++d) interpolant += interpolantVec[d]*n[d];
722           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);}
723           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ;
724         }
725       }
726       comp        += Ncomp;
727       fieldOffset += Nb*Ncomp;
728     }
729     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
730     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
731     localDiff += elemDiff;
732   }
733   ierr  = PetscFree7(funcVal,coords,realSpaceDer,v0,J,invJ,interpolantVec);CHKERRQ(ierr);
734   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
735   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
736   *diff = PetscSqrtReal(*diff);
737   PetscFunctionReturn(0);
738 }
739 
740 #undef __FUNCT__
741 #define __FUNCT__ "DMPlexComputeL2FieldDiff"
742 /*@C
743   DMPlexComputeL2FieldDiff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h, separated into field components.
744 
745   Input Parameters:
746 + dm    - The DM
747 . funcs - The functions to evaluate for each field component
748 . ctxs  - Optional array of contexts to pass to each function, or NULL.
749 - X     - The coefficient vector u_h
750 
751   Output Parameter:
752 . diff - The array of differences, ||u^f - u^f_h||_2
753 
754   Level: developer
755 
756 .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff(), DMPlexComputeL2GradientDiff()
757 @*/
758 PetscErrorCode DMPlexComputeL2FieldDiff(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal diff[])
759 {
760   const PetscInt   debug = 0;
761   PetscSection     section;
762   PetscQuadrature  quad;
763   Vec              localX;
764   PetscScalar     *funcVal, *interpolant;
765   PetscReal       *coords, *v0, *J, *invJ, detJ;
766   PetscReal       *localDiff;
767   const PetscReal *quadPoints, *quadWeights;
768   PetscInt         dim, numFields, numComponents = 0, numQuadPoints, cStart, cEnd, c, field, fieldOffset;
769   PetscErrorCode   ierr;
770 
771   PetscFunctionBegin;
772   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
773   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
774   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
775   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
776   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
777   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
778   for (field = 0; field < numFields; ++field) {
779     PetscObject  obj;
780     PetscClassId id;
781     PetscInt     Nc;
782 
783     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
784     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
785     if (id == PETSCFE_CLASSID) {
786       PetscFE fe = (PetscFE) obj;
787 
788       ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
789       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
790     } else if (id == PETSCFV_CLASSID) {
791       PetscFV fv = (PetscFV) obj;
792 
793       ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr);
794       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
795     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
796     numComponents += Nc;
797   }
798   ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
799   ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
800   ierr = PetscCalloc7(numFields,&localDiff,numComponents,&funcVal,numComponents,&interpolant,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
801   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
802   for (c = cStart; c < cEnd; ++c) {
803     PetscScalar *x = NULL;
804 
805     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
806     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
807     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
808 
809     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
810       PetscObject  obj;
811       PetscClassId id;
812       void * const ctx = ctxs ? ctxs[field] : NULL;
813       PetscInt     Nb, Nc, q, fc;
814 
815       PetscReal       elemDiff = 0.0;
816 
817       ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
818       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
819       if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);}
820       else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;}
821       else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
822       if (debug) {
823         char title[1024];
824         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
825         ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr);
826       }
827       for (q = 0; q < numQuadPoints; ++q) {
828         CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], coords);
829         (*funcs[field])(coords, funcVal, ctx);
830         if (id == PETSCFE_CLASSID)      {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
831         else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
832         else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
833         for (fc = 0; fc < Nc; ++fc) {
834           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);}
835           elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ;
836         }
837       }
838       fieldOffset += Nb*Nc;
839       localDiff[field] += elemDiff;
840     }
841     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
842   }
843   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
844   ierr = MPI_Allreduce(localDiff, diff, numFields, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
845   for (field = 0; field < numFields; ++field) diff[field] = PetscSqrtReal(diff[field]);
846   ierr = PetscFree7(localDiff,funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr);
847   PetscFunctionReturn(0);
848 }
849 
850 #undef __FUNCT__
851 #define __FUNCT__ "DMPlexComputeIntegralFEM"
852 /*@
853   DMPlexComputeIntegralFEM - Form the local integral F from the local input X using pointwise functions specified by the user
854 
855   Input Parameters:
856 + dm - The mesh
857 . X  - Local input vector
858 - user - The user context
859 
860   Output Parameter:
861 . integral - Local integral for each field
862 
863   Level: developer
864 
865 .seealso: DMPlexComputeResidualFEM()
866 @*/
867 PetscErrorCode DMPlexComputeIntegralFEM(DM dm, Vec X, PetscReal *integral, void *user)
868 {
869   DM_Plex          *mesh  = (DM_Plex *) dm->data;
870   DM                dmAux;
871   Vec               localX, A;
872   PetscDS           prob, probAux = NULL;
873   PetscQuadrature   q;
874   PetscSection      section, sectionAux;
875   PetscFECellGeom  *cgeom;
876   PetscScalar      *u, *a = NULL;
877   PetscInt          dim, Nf, f, numCells, cStart, cEnd, c;
878   PetscInt          totDim, totDimAux;
879   PetscErrorCode    ierr;
880 
881   PetscFunctionBegin;
882   /*ierr = PetscLogEventBegin(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);*/
883   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
884   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
885   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
886   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
887   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
888   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
889   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
890   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
891   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
892   numCells = cEnd - cStart;
893   for (f = 0; f < Nf; ++f) {integral[f]    = 0.0;}
894   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
895   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
896   if (dmAux) {
897     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
898     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
899     ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
900   }
901   ierr = DMPlexInsertBoundaryValuesFEM(dm, localX);CHKERRQ(ierr);
902   ierr = PetscMalloc2(numCells*totDim,&u,numCells,&cgeom);CHKERRQ(ierr);
903   if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);}
904   for (c = cStart; c < cEnd; ++c) {
905     PetscScalar *x = NULL;
906     PetscInt     i;
907 
908     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, cgeom[c].v0, cgeom[c].J, cgeom[c].invJ, &cgeom[c].detJ);CHKERRQ(ierr);
909     if (cgeom[c].detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", cgeom[c].detJ, c);
910     ierr = DMPlexVecGetClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr);
911     for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i];
912     ierr = DMPlexVecRestoreClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr);
913     if (dmAux) {
914       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
915       for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i];
916       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
917     }
918   }
919   for (f = 0; f < Nf; ++f) {
920     PetscFE  fe;
921     PetscInt numQuadPoints, Nb;
922     /* Conforming batches */
923     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
924     /* Remainder */
925     PetscInt Nr, offset;
926 
927     ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
928     ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
929     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
930     ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
931     ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
932     blockSize = Nb*numQuadPoints;
933     batchSize = numBlocks * blockSize;
934     ierr =  PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
935     numChunks = numCells / (numBatches*batchSize);
936     Ne        = numChunks*numBatches*batchSize;
937     Nr        = numCells % (numBatches*batchSize);
938     offset    = numCells - Nr;
939     ierr = PetscFEIntegrate(fe, prob, f, Ne, cgeom, u, probAux, a, integral);CHKERRQ(ierr);
940     ierr = PetscFEIntegrate(fe, prob, f, Nr, &cgeom[offset], &u[offset*totDim], probAux, &a[offset*totDimAux], integral);CHKERRQ(ierr);
941   }
942   ierr = PetscFree2(u,cgeom);CHKERRQ(ierr);
943   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
944   if (mesh->printFEM) {
945     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Local integral:");CHKERRQ(ierr);
946     for (f = 0; f < Nf; ++f) {ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), " %g", integral[f]);CHKERRQ(ierr);}
947     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "\n");CHKERRQ(ierr);
948   }
949   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
950   /* TODO: Allreduce for integral */
951   /*ierr = PetscLogEventEnd(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);*/
952   PetscFunctionReturn(0);
953 }
954 
955 #undef __FUNCT__
956 #define __FUNCT__ "DMPlexComputeInterpolatorFEM"
957 /*@
958   DMPlexComputeInterpolatorFEM - Form the local portion of the interpolation matrix I from the coarse DM to the uniformly refined DM.
959 
960   Input Parameters:
961 + dmf  - The fine mesh
962 . dmc  - The coarse mesh
963 - user - The user context
964 
965   Output Parameter:
966 . In  - The interpolation matrix
967 
968   Note:
969   The first member of the user context must be an FEMContext.
970 
971   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
972   like a GPU, or vectorize on a multicore machine.
973 
974   Level: developer
975 
976 .seealso: DMPlexComputeJacobianFEM()
977 @*/
978 PetscErrorCode DMPlexComputeInterpolatorFEM(DM dmc, DM dmf, Mat In, void *user)
979 {
980   DM_Plex          *mesh  = (DM_Plex *) dmc->data;
981   const char       *name  = "Interpolator";
982   PetscDS           prob;
983   PetscFE          *feRef;
984   PetscSection      fsection, fglobalSection;
985   PetscSection      csection, cglobalSection;
986   PetscScalar      *elemMat;
987   PetscInt          dim, Nf, f, fieldI, fieldJ, offsetI, offsetJ, cStart, cEnd, c;
988   PetscInt          cTotDim, rTotDim = 0;
989   PetscErrorCode    ierr;
990 
991   PetscFunctionBegin;
992   ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
993   ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr);
994   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
995   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
996   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
997   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
998   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
999   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
1000   ierr = DMGetDS(dmf, &prob);CHKERRQ(ierr);
1001   ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr);
1002   for (f = 0; f < Nf; ++f) {
1003     PetscFE  fe;
1004     PetscInt rNb, Nc;
1005 
1006     ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
1007     ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr);
1008     ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr);
1009     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1010     rTotDim += rNb*Nc;
1011   }
1012   ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr);
1013   ierr = PetscMalloc1(rTotDim*cTotDim,&elemMat);CHKERRQ(ierr);
1014   ierr = PetscMemzero(elemMat, rTotDim*cTotDim * sizeof(PetscScalar));CHKERRQ(ierr);
1015   for (fieldI = 0, offsetI = 0; fieldI < Nf; ++fieldI) {
1016     PetscDualSpace   Qref;
1017     PetscQuadrature  f;
1018     const PetscReal *qpoints, *qweights;
1019     PetscReal       *points;
1020     PetscInt         npoints = 0, Nc, Np, fpdim, i, k, p, d;
1021 
1022     /* Compose points from all dual basis functionals */
1023     ierr = PetscFEGetDualSpace(feRef[fieldI], &Qref);CHKERRQ(ierr);
1024     ierr = PetscFEGetNumComponents(feRef[fieldI], &Nc);CHKERRQ(ierr);
1025     ierr = PetscDualSpaceGetDimension(Qref, &fpdim);CHKERRQ(ierr);
1026     for (i = 0; i < fpdim; ++i) {
1027       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1028       ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, NULL);CHKERRQ(ierr);
1029       npoints += Np;
1030     }
1031     ierr = PetscMalloc1(npoints*dim,&points);CHKERRQ(ierr);
1032     for (i = 0, k = 0; i < fpdim; ++i) {
1033       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1034       ierr = PetscQuadratureGetData(f, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr);
1035       for (p = 0; p < Np; ++p, ++k) for (d = 0; d < dim; ++d) points[k*dim+d] = qpoints[p*dim+d];
1036     }
1037 
1038     for (fieldJ = 0, offsetJ = 0; fieldJ < Nf; ++fieldJ) {
1039       PetscFE    fe;
1040       PetscReal *B;
1041       PetscInt   NcJ, cpdim, j;
1042 
1043       /* Evaluate basis at points */
1044       ierr = PetscDSGetDiscretization(prob, fieldJ, (PetscObject *) &fe);CHKERRQ(ierr);
1045       ierr = PetscFEGetNumComponents(fe, &NcJ);CHKERRQ(ierr);
1046       ierr = PetscFEGetDimension(fe, &cpdim);CHKERRQ(ierr);
1047       /* For now, fields only interpolate themselves */
1048       if (fieldI == fieldJ) {
1049         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);
1050         ierr = PetscFEGetTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);
1051         for (i = 0, k = 0; i < fpdim; ++i) {
1052           ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1053           ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, &qweights);CHKERRQ(ierr);
1054           for (p = 0; p < Np; ++p, ++k) {
1055             for (j = 0; j < cpdim; ++j) {
1056               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];
1057             }
1058           }
1059         }
1060         ierr = PetscFERestoreTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr);
1061       }
1062       offsetJ += cpdim*NcJ;
1063     }
1064     offsetI += fpdim*Nc;
1065     ierr = PetscFree(points);CHKERRQ(ierr);
1066   }
1067   if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(0, name, rTotDim, cTotDim, elemMat);CHKERRQ(ierr);}
1068   /* Preallocate matrix */
1069   {
1070     PetscHashJK ht;
1071     PetscLayout rLayout;
1072     PetscInt   *dnz, *onz, *cellCIndices, *cellFIndices;
1073     PetscInt    locRows, rStart, rEnd, cell, r;
1074 
1075     ierr = MatGetLocalSize(In, &locRows, NULL);CHKERRQ(ierr);
1076     ierr = PetscLayoutCreate(PetscObjectComm((PetscObject) In), &rLayout);CHKERRQ(ierr);
1077     ierr = PetscLayoutSetLocalSize(rLayout, locRows);CHKERRQ(ierr);
1078     ierr = PetscLayoutSetBlockSize(rLayout, 1);CHKERRQ(ierr);
1079     ierr = PetscLayoutSetUp(rLayout);CHKERRQ(ierr);
1080     ierr = PetscLayoutGetRange(rLayout, &rStart, &rEnd);CHKERRQ(ierr);
1081     ierr = PetscLayoutDestroy(&rLayout);CHKERRQ(ierr);
1082     ierr = PetscCalloc4(locRows,&dnz,locRows,&onz,cTotDim,&cellCIndices,rTotDim,&cellFIndices);CHKERRQ(ierr);
1083     ierr = PetscHashJKCreate(&ht);CHKERRQ(ierr);
1084     for (cell = cStart; cell < cEnd; ++cell) {
1085       ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, cell, cellCIndices, cellFIndices);CHKERRQ(ierr);
1086       for (r = 0; r < rTotDim; ++r) {
1087         PetscHashJKKey  key;
1088         PetscHashJKIter missing, iter;
1089 
1090         key.j = cellFIndices[r];
1091         if (key.j < 0) continue;
1092         for (c = 0; c < cTotDim; ++c) {
1093           key.k = cellCIndices[c];
1094           if (key.k < 0) continue;
1095           ierr = PetscHashJKPut(ht, key, &missing, &iter);CHKERRQ(ierr);
1096           if (missing) {
1097             ierr = PetscHashJKSet(ht, iter, 1);CHKERRQ(ierr);
1098             if ((key.k >= rStart) && (key.k < rEnd)) ++dnz[key.j-rStart];
1099             else                                     ++onz[key.j-rStart];
1100           }
1101         }
1102       }
1103     }
1104     ierr = PetscHashJKDestroy(&ht);CHKERRQ(ierr);
1105     ierr = MatXAIJSetPreallocation(In, 1, dnz, onz, NULL, NULL);CHKERRQ(ierr);
1106     ierr = MatSetOption(In, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
1107     ierr = PetscFree4(dnz,onz,cellCIndices,cellFIndices);CHKERRQ(ierr);
1108   }
1109   /* Fill matrix */
1110   ierr = MatZeroEntries(In);CHKERRQ(ierr);
1111   for (c = cStart; c < cEnd; ++c) {
1112     ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES);CHKERRQ(ierr);
1113   }
1114   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);}
1115   ierr = PetscFree(feRef);CHKERRQ(ierr);
1116   ierr = PetscFree(elemMat);CHKERRQ(ierr);
1117   ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1118   ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1119   if (mesh->printFEM) {
1120     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
1121     ierr = MatChop(In, 1.0e-10);CHKERRQ(ierr);
1122     ierr = MatView(In, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1123   }
1124   ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1125   PetscFunctionReturn(0);
1126 }
1127 
1128 #undef __FUNCT__
1129 #define __FUNCT__ "DMPlexComputeInjectorFEM"
1130 PetscErrorCode DMPlexComputeInjectorFEM(DM dmc, DM dmf, VecScatter *sc, void *user)
1131 {
1132   PetscDS        prob;
1133   PetscFE       *feRef;
1134   Vec            fv, cv;
1135   IS             fis, cis;
1136   PetscSection   fsection, fglobalSection, csection, cglobalSection;
1137   PetscInt      *cmap, *cellCIndices, *cellFIndices, *cindices, *findices;
1138   PetscInt       cTotDim, fTotDim = 0, Nf, f, field, cStart, cEnd, c, dim, d, startC, offsetC, offsetF, m;
1139   PetscErrorCode ierr;
1140 
1141   PetscFunctionBegin;
1142   ierr = PetscLogEventBegin(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1143   ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr);
1144   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
1145   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
1146   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
1147   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
1148   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
1149   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
1150   ierr = DMGetDS(dmc, &prob);CHKERRQ(ierr);
1151   ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr);
1152   for (f = 0; f < Nf; ++f) {
1153     PetscFE  fe;
1154     PetscInt fNb, Nc;
1155 
1156     ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
1157     ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr);
1158     ierr = PetscFEGetDimension(feRef[f], &fNb);CHKERRQ(ierr);
1159     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1160     fTotDim += fNb*Nc;
1161   }
1162   ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr);
1163   ierr = PetscMalloc1(cTotDim,&cmap);CHKERRQ(ierr);
1164   for (field = 0, offsetC = 0, offsetF = 0; field < Nf; ++field) {
1165     PetscFE        feC;
1166     PetscDualSpace QF, QC;
1167     PetscInt       NcF, NcC, fpdim, cpdim;
1168 
1169     ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &feC);CHKERRQ(ierr);
1170     ierr = PetscFEGetNumComponents(feC, &NcC);CHKERRQ(ierr);
1171     ierr = PetscFEGetNumComponents(feRef[field], &NcF);CHKERRQ(ierr);
1172     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);
1173     ierr = PetscFEGetDualSpace(feRef[field], &QF);CHKERRQ(ierr);
1174     ierr = PetscDualSpaceGetDimension(QF, &fpdim);CHKERRQ(ierr);
1175     ierr = PetscFEGetDualSpace(feC, &QC);CHKERRQ(ierr);
1176     ierr = PetscDualSpaceGetDimension(QC, &cpdim);CHKERRQ(ierr);
1177     for (c = 0; c < cpdim; ++c) {
1178       PetscQuadrature  cfunc;
1179       const PetscReal *cqpoints;
1180       PetscInt         NpC;
1181 
1182       ierr = PetscDualSpaceGetFunctional(QC, c, &cfunc);CHKERRQ(ierr);
1183       ierr = PetscQuadratureGetData(cfunc, NULL, &NpC, &cqpoints, NULL);CHKERRQ(ierr);
1184       if (NpC != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not know how to do injection for moments");
1185       for (f = 0; f < fpdim; ++f) {
1186         PetscQuadrature  ffunc;
1187         const PetscReal *fqpoints;
1188         PetscReal        sum = 0.0;
1189         PetscInt         NpF, comp;
1190 
1191         ierr = PetscDualSpaceGetFunctional(QF, f, &ffunc);CHKERRQ(ierr);
1192         ierr = PetscQuadratureGetData(ffunc, NULL, &NpF, &fqpoints, NULL);CHKERRQ(ierr);
1193         if (NpC != NpF) continue;
1194         for (d = 0; d < dim; ++d) sum += PetscAbsReal(cqpoints[d] - fqpoints[d]);
1195         if (sum > 1.0e-9) continue;
1196         for (comp = 0; comp < NcC; ++comp) {
1197           cmap[(offsetC+c)*NcC+comp] = (offsetF+f)*NcF+comp;
1198         }
1199         break;
1200       }
1201     }
1202     offsetC += cpdim*NcC;
1203     offsetF += fpdim*NcF;
1204   }
1205   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);}
1206   ierr = PetscFree(feRef);CHKERRQ(ierr);
1207 
1208   ierr = DMGetGlobalVector(dmf, &fv);CHKERRQ(ierr);
1209   ierr = DMGetGlobalVector(dmc, &cv);CHKERRQ(ierr);
1210   ierr = VecGetOwnershipRange(cv, &startC, NULL);CHKERRQ(ierr);
1211   ierr = PetscSectionGetConstrainedStorageSize(cglobalSection, &m);CHKERRQ(ierr);
1212   ierr = PetscMalloc2(cTotDim,&cellCIndices,fTotDim,&cellFIndices);CHKERRQ(ierr);
1213   ierr = PetscMalloc1(m,&cindices);CHKERRQ(ierr);
1214   ierr = PetscMalloc1(m,&findices);CHKERRQ(ierr);
1215   for (d = 0; d < m; ++d) cindices[d] = findices[d] = -1;
1216   for (c = cStart; c < cEnd; ++c) {
1217     ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, c, cellCIndices, cellFIndices);CHKERRQ(ierr);
1218     for (d = 0; d < cTotDim; ++d) {
1219       if (cellCIndices[d] < 0) continue;
1220       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]]);
1221       cindices[cellCIndices[d]-startC] = cellCIndices[d];
1222       findices[cellCIndices[d]-startC] = cellFIndices[cmap[d]];
1223     }
1224   }
1225   ierr = PetscFree(cmap);CHKERRQ(ierr);
1226   ierr = PetscFree2(cellCIndices,cellFIndices);CHKERRQ(ierr);
1227 
1228   ierr = ISCreateGeneral(PETSC_COMM_SELF, m, cindices, PETSC_OWN_POINTER, &cis);CHKERRQ(ierr);
1229   ierr = ISCreateGeneral(PETSC_COMM_SELF, m, findices, PETSC_OWN_POINTER, &fis);CHKERRQ(ierr);
1230   ierr = VecScatterCreate(cv, cis, fv, fis, sc);CHKERRQ(ierr);
1231   ierr = ISDestroy(&cis);CHKERRQ(ierr);
1232   ierr = ISDestroy(&fis);CHKERRQ(ierr);
1233   ierr = DMRestoreGlobalVector(dmf, &fv);CHKERRQ(ierr);
1234   ierr = DMRestoreGlobalVector(dmc, &cv);CHKERRQ(ierr);
1235   ierr = PetscLogEventEnd(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1236   PetscFunctionReturn(0);
1237 }
1238