xref: /petsc/src/dm/impls/plex/plexfem.c (revision b92fff86961fc3c128ddb069ff024ef7e7ed3f53)
1 #include <petsc-private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2 
3 #include <petscfe.h>
4 
5 #undef __FUNCT__
6 #define __FUNCT__ "DMPlexGetScale"
7 PetscErrorCode DMPlexGetScale(DM dm, PetscUnit unit, PetscReal *scale)
8 {
9   DM_Plex *mesh = (DM_Plex*) dm->data;
10 
11   PetscFunctionBegin;
12   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
13   PetscValidPointer(scale, 3);
14   *scale = mesh->scale[unit];
15   PetscFunctionReturn(0);
16 }
17 
18 #undef __FUNCT__
19 #define __FUNCT__ "DMPlexSetScale"
20 PetscErrorCode DMPlexSetScale(DM dm, PetscUnit unit, PetscReal scale)
21 {
22   DM_Plex *mesh = (DM_Plex*) dm->data;
23 
24   PetscFunctionBegin;
25   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
26   mesh->scale[unit] = scale;
27   PetscFunctionReturn(0);
28 }
29 
30 PETSC_STATIC_INLINE PetscInt epsilon(PetscInt i, PetscInt j, PetscInt k)
31 {
32   switch (i) {
33   case 0:
34     switch (j) {
35     case 0: return 0;
36     case 1:
37       switch (k) {
38       case 0: return 0;
39       case 1: return 0;
40       case 2: return 1;
41       }
42     case 2:
43       switch (k) {
44       case 0: return 0;
45       case 1: return -1;
46       case 2: return 0;
47       }
48     }
49   case 1:
50     switch (j) {
51     case 0:
52       switch (k) {
53       case 0: return 0;
54       case 1: return 0;
55       case 2: return -1;
56       }
57     case 1: return 0;
58     case 2:
59       switch (k) {
60       case 0: return 1;
61       case 1: return 0;
62       case 2: return 0;
63       }
64     }
65   case 2:
66     switch (j) {
67     case 0:
68       switch (k) {
69       case 0: return 0;
70       case 1: return 1;
71       case 2: return 0;
72       }
73     case 1:
74       switch (k) {
75       case 0: return -1;
76       case 1: return 0;
77       case 2: return 0;
78       }
79     case 2: return 0;
80     }
81   }
82   return 0;
83 }
84 
85 #undef __FUNCT__
86 #define __FUNCT__ "DMPlexCreateRigidBody"
87 /*@C
88   DMPlexCreateRigidBody - create rigid body modes from coordinates
89 
90   Collective on DM
91 
92   Input Arguments:
93 + dm - the DM
94 . section - the local section associated with the rigid field, or NULL for the default section
95 - globalSection - the global section associated with the rigid field, or NULL for the default section
96 
97   Output Argument:
98 . sp - the null space
99 
100   Note: This is necessary to take account of Dirichlet conditions on the displacements
101 
102   Level: advanced
103 
104 .seealso: MatNullSpaceCreate()
105 @*/
106 PetscErrorCode DMPlexCreateRigidBody(DM dm, PetscSection section, PetscSection globalSection, MatNullSpace *sp)
107 {
108   MPI_Comm       comm;
109   Vec            coordinates, localMode, mode[6];
110   PetscSection   coordSection;
111   PetscScalar   *coords;
112   PetscInt       dim, vStart, vEnd, v, n, m, d, i, j;
113   PetscErrorCode ierr;
114 
115   PetscFunctionBegin;
116   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
117   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
118   if (dim == 1) {
119     ierr = MatNullSpaceCreate(comm, PETSC_TRUE, 0, NULL, sp);CHKERRQ(ierr);
120     PetscFunctionReturn(0);
121   }
122   if (!section)       {ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);}
123   if (!globalSection) {ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);}
124   ierr = PetscSectionGetConstrainedStorageSize(globalSection, &n);CHKERRQ(ierr);
125   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
126   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
127   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
128   m    = (dim*(dim+1))/2;
129   ierr = VecCreate(comm, &mode[0]);CHKERRQ(ierr);
130   ierr = VecSetSizes(mode[0], n, PETSC_DETERMINE);CHKERRQ(ierr);
131   ierr = VecSetUp(mode[0]);CHKERRQ(ierr);
132   for (i = 1; i < m; ++i) {ierr = VecDuplicate(mode[0], &mode[i]);CHKERRQ(ierr);}
133   /* Assume P1 */
134   ierr = DMGetLocalVector(dm, &localMode);CHKERRQ(ierr);
135   for (d = 0; d < dim; ++d) {
136     PetscScalar values[3] = {0.0, 0.0, 0.0};
137 
138     values[d] = 1.0;
139     ierr      = VecSet(localMode, 0.0);CHKERRQ(ierr);
140     for (v = vStart; v < vEnd; ++v) {
141       ierr = DMPlexVecSetClosure(dm, section, localMode, v, values, INSERT_VALUES);CHKERRQ(ierr);
142     }
143     ierr = DMLocalToGlobalBegin(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
144     ierr = DMLocalToGlobalEnd(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
145   }
146   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
147   for (d = dim; d < dim*(dim+1)/2; ++d) {
148     PetscInt i, j, k = dim > 2 ? d - dim : d;
149 
150     ierr = VecSet(localMode, 0.0);CHKERRQ(ierr);
151     for (v = vStart; v < vEnd; ++v) {
152       PetscScalar values[3] = {0.0, 0.0, 0.0};
153       PetscInt    off;
154 
155       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
156       for (i = 0; i < dim; ++i) {
157         for (j = 0; j < dim; ++j) {
158           values[j] += epsilon(i, j, k)*PetscRealPart(coords[off+i]);
159         }
160       }
161       ierr = DMPlexVecSetClosure(dm, section, localMode, v, values, INSERT_VALUES);CHKERRQ(ierr);
162     }
163     ierr = DMLocalToGlobalBegin(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
164     ierr = DMLocalToGlobalEnd(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
165   }
166   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
167   ierr = DMRestoreLocalVector(dm, &localMode);CHKERRQ(ierr);
168   for (i = 0; i < dim; ++i) {ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);}
169   /* Orthonormalize system */
170   for (i = dim; i < m; ++i) {
171     PetscScalar dots[6];
172 
173     ierr = VecMDot(mode[i], i, mode, dots);CHKERRQ(ierr);
174     for (j = 0; j < i; ++j) dots[j] *= -1.0;
175     ierr = VecMAXPY(mode[i], i, dots, mode);CHKERRQ(ierr);
176     ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);
177   }
178   ierr = MatNullSpaceCreate(comm, PETSC_FALSE, m, mode, sp);CHKERRQ(ierr);
179   for (i = 0; i< m; ++i) {ierr = VecDestroy(&mode[i]);CHKERRQ(ierr);}
180   PetscFunctionReturn(0);
181 }
182 
183 #undef __FUNCT__
184 #define __FUNCT__ "DMPlexProjectFunctionLabelLocal"
185 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)
186 {
187   PetscDualSpace *sp;
188   PetscSection    section;
189   PetscScalar    *values;
190   PetscReal      *v0, *J, detJ;
191   PetscInt        numFields, numComp, dim, spDim, totDim = 0, numValues, cStart, cEnd, f, d, v, i, comp;
192   PetscErrorCode  ierr;
193 
194   PetscFunctionBegin;
195   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
196   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
197   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
198   ierr = PetscMalloc3(numFields,&sp,dim,&v0,dim*dim,&J);CHKERRQ(ierr);
199   for (f = 0; f < numFields; ++f) {
200     ierr = PetscFEGetDualSpace(fe[f], &sp[f]);CHKERRQ(ierr);
201     ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr);
202     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
203     totDim += spDim*numComp;
204   }
205   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
206   ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr);
207   if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section cell closure size %d != dual space dimension %d", numValues, totDim);
208   ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
209   for (i = 0; i < numIds; ++i) {
210     IS              pointIS;
211     const PetscInt *points;
212     PetscInt        n, p;
213 
214     ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr);
215     ierr = ISGetLocalSize(pointIS, &n);CHKERRQ(ierr);
216     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
217     for (p = 0; p < n; ++p) {
218       const PetscInt    point = points[p];
219       PetscCellGeometry geom;
220 
221       if ((point < cStart) || (point >= cEnd)) continue;
222       ierr = DMPlexComputeCellGeometry(dm, point, v0, J, NULL, &detJ);CHKERRQ(ierr);
223       geom.v0   = v0;
224       geom.J    = J;
225       geom.detJ = &detJ;
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 = DMPlexVecSetClosure(dm, section, localX, 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 = PetscFree3(sp,v0,J);CHKERRQ(ierr);
246   PetscFunctionReturn(0);
247 }
248 
249 #undef __FUNCT__
250 #define __FUNCT__ "DMPlexProjectFunctionLocal"
251 PetscErrorCode DMPlexProjectFunctionLocal(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
252 {
253   PetscDualSpace *sp;
254   PetscSection    section;
255   PetscScalar    *values;
256   PetscReal      *v0, *J, detJ;
257   PetscInt        numFields, numComp, dim, spDim, totDim = 0, numValues, cStart, cEnd, c, f, d, v, comp;
258   PetscErrorCode  ierr;
259 
260   PetscFunctionBegin;
261   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
262   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
263   ierr = PetscMalloc1(numFields, &sp);CHKERRQ(ierr);
264   for (f = 0; f < numFields; ++f) {
265     ierr = PetscFEGetDualSpace(fe[f], &sp[f]);CHKERRQ(ierr);
266     ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr);
267     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
268     totDim += spDim*numComp;
269   }
270   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
271   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
272   ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr);
273   if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section cell closure size %d != dual space dimension %d", numValues, totDim);
274   ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
275   ierr = PetscMalloc2(dim,&v0,dim*dim,&J);CHKERRQ(ierr);
276   for (c = cStart; c < cEnd; ++c) {
277     PetscCellGeometry geom;
278 
279     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, NULL, &detJ);CHKERRQ(ierr);
280     geom.v0   = v0;
281     geom.J    = J;
282     geom.detJ = &detJ;
283     for (f = 0, v = 0; f < numFields; ++f) {
284       void * const ctx = ctxs ? ctxs[f] : NULL;
285       ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr);
286       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
287       for (d = 0; d < spDim; ++d) {
288         if (funcs[f]) {
289           ierr = PetscDualSpaceApply(sp[f], d, geom, numComp, funcs[f], ctx, &values[v]);CHKERRQ(ierr);
290         } else {
291           for (comp = 0; comp < numComp; ++comp) values[v+comp] = 0.0;
292         }
293         v += numComp;
294       }
295     }
296     ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr);
297   }
298   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
299   ierr = PetscFree2(v0,J);CHKERRQ(ierr);
300   ierr = PetscFree(sp);CHKERRQ(ierr);
301   PetscFunctionReturn(0);
302 }
303 
304 #undef __FUNCT__
305 #define __FUNCT__ "DMPlexProjectFunction"
306 /*@C
307   DMPlexProjectFunction - This projects the given function into the function space provided.
308 
309   Input Parameters:
310 + dm      - The DM
311 . fe      - The PetscFE associated with the field
312 . funcs   - The coordinate functions to evaluate, one per field
313 . ctxs    - Optional array of contexts to pass to each coordinate function.  ctxs itself may be null.
314 - mode    - The insertion mode for values
315 
316   Output Parameter:
317 . X - vector
318 
319   Level: developer
320 
321 .seealso: DMPlexComputeL2Diff()
322 @*/
323 PetscErrorCode DMPlexProjectFunction(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec X)
324 {
325   Vec            localX;
326   PetscErrorCode ierr;
327 
328   PetscFunctionBegin;
329   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
330   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
331   ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, ctxs, mode, localX);CHKERRQ(ierr);
332   ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr);
333   ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr);
334   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
335   PetscFunctionReturn(0);
336 }
337 
338 #undef __FUNCT__
339 #define __FUNCT__ "DMPlexInsertBoundaryValuesFEM"
340 PetscErrorCode DMPlexInsertBoundaryValuesFEM(DM dm, Vec localX)
341 {
342   void        (**funcs)(const PetscReal x[], PetscScalar *u, void *ctx);
343   void         **ctxs;
344   PetscFE       *fe;
345   PetscInt       numFields, f, numBd, b;
346   PetscErrorCode ierr;
347 
348   PetscFunctionBegin;
349   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
350   PetscValidHeaderSpecific(localX, VEC_CLASSID, 2);
351   ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr);
352   ierr = PetscMalloc3(numFields,&fe,numFields,&funcs,numFields,&ctxs);CHKERRQ(ierr);
353   for (f = 0; f < numFields; ++f) {ierr = DMGetField(dm, f, (PetscObject *) &fe[f]);CHKERRQ(ierr);}
354   /* OPT: Could attempt to do multiple BCs at once */
355   ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
356   for (b = 0; b < numBd; ++b) {
357     DMLabel         label;
358     const PetscInt *ids;
359     const char     *name;
360     PetscInt        numids, field;
361     PetscBool       isEssential;
362     void          (*func)();
363     void           *ctx;
364 
365     /* TODO: We need to set only the part indicated by the ids */
366     ierr = DMPlexGetBoundary(dm, b, &isEssential, &name, &field, &func, &numids, &ids, &ctx);CHKERRQ(ierr);
367     ierr = DMPlexGetLabel(dm, name, &label);CHKERRQ(ierr);
368     for (f = 0; f < numFields; ++f) {
369       funcs[f] = field == f ? (void (*)(const PetscReal[], PetscScalar *, void *)) func : NULL;
370       ctxs[f]  = field == f ? ctx : NULL;
371     }
372     ierr = DMPlexProjectFunctionLabelLocal(dm, label, numids, ids, fe, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
373   }
374   ierr = PetscFree3(fe,funcs,ctxs);CHKERRQ(ierr);
375   PetscFunctionReturn(0);
376 }
377 
378 /* Assuming dim == 3 */
379 typedef struct {
380   PetscScalar normal[3];   /* Area-scaled normals */
381   PetscScalar centroid[3]; /* Location of centroid (quadrature point) */
382   PetscScalar grad[2][3];  /* Face contribution to gradient in left and right cell */
383 } FaceGeom;
384 
385 #undef __FUNCT__
386 #define __FUNCT__ "DMPlexInsertBoundaryValuesFVM"
387 PetscErrorCode DMPlexInsertBoundaryValuesFVM(DM dm, PetscReal time, Vec locX)
388 {
389   DM                 dmFace;
390   Vec                faceGeometry;
391   DMLabel            label;
392   const PetscScalar *facegeom;
393   PetscScalar       *x;
394   PetscInt           numBd, b;
395   PetscErrorCode     ierr;
396 
397   PetscFunctionBegin;
398   /* TODO Pull this geometry calculation into the library */
399   ierr = PetscObjectQuery((PetscObject) dm, "FaceGeometry", (PetscObject *) &faceGeometry);CHKERRQ(ierr);
400   ierr = DMPlexGetLabel(dm, "Face Sets", &label);CHKERRQ(ierr);
401   ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
402   ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr);
403   ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
404   ierr = VecGetArray(locX, &x);CHKERRQ(ierr);
405   for (b = 0; b < numBd; ++b) {
406     PetscErrorCode (*func)(PetscReal,const PetscScalar*,const PetscScalar*,const PetscScalar*,PetscScalar*,void*);
407     const PetscInt  *ids;
408     PetscInt         numids, i;
409     void            *ctx;
410 
411     ierr = DMPlexGetBoundary(dm, b, NULL, NULL, NULL, (void (**)()) &func, &numids, &ids, &ctx);CHKERRQ(ierr);
412     for (i = 0; i < numids; ++i) {
413       IS              faceIS;
414       const PetscInt *faces;
415       PetscInt        numFaces, f;
416 
417       ierr = DMLabelGetStratumIS(label, ids[i], &faceIS);CHKERRQ(ierr);
418       if (!faceIS) continue; /* No points with that id on this process */
419       ierr = ISGetLocalSize(faceIS, &numFaces);CHKERRQ(ierr);
420       ierr = ISGetIndices(faceIS, &faces);CHKERRQ(ierr);
421       for (f = 0; f < numFaces; ++f) {
422         const PetscInt     face = faces[f], *cells;
423         const PetscScalar *xI;
424         PetscScalar       *xG;
425         const FaceGeom    *fg;
426 
427         ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr);
428         ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr);
429         ierr = DMPlexPointLocalRead(dm, cells[0], x, &xI);CHKERRQ(ierr);
430         ierr = DMPlexPointLocalRef(dm, cells[1], x, &xG);CHKERRQ(ierr);
431         ierr = (*func)(time, fg->centroid, fg->normal, xI, xG, ctx);CHKERRQ(ierr);
432       }
433       ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr);
434       ierr = ISDestroy(&faceIS);CHKERRQ(ierr);
435     }
436   }
437   ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
438   ierr = VecRestoreArray(locX, &x);CHKERRQ(ierr);
439   PetscFunctionReturn(0);
440 }
441 
442 #undef __FUNCT__
443 #define __FUNCT__ "DMPlexComputeL2Diff"
444 /*@C
445   DMPlexComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.
446 
447   Input Parameters:
448 + dm    - The DM
449 . fe    - The PetscFE object for each field
450 . funcs - The functions to evaluate for each field component
451 . ctxs  - Optional array of contexts to pass to each function, or NULL.
452 - X     - The coefficient vector u_h
453 
454   Output Parameter:
455 . diff - The diff ||u - u_h||_2
456 
457   Level: developer
458 
459 .seealso: DMPlexProjectFunction(), DMPlexComputeL2GradientDiff()
460 @*/
461 PetscErrorCode DMPlexComputeL2Diff(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
462 {
463   const PetscInt  debug = 0;
464   PetscSection    section;
465   PetscQuadrature quad;
466   Vec             localX;
467   PetscScalar    *funcVal;
468   PetscReal      *coords, *v0, *J, *invJ, detJ;
469   PetscReal       localDiff = 0.0;
470   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp;
471   PetscErrorCode  ierr;
472 
473   PetscFunctionBegin;
474   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
475   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
476   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
477   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
478   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
479   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
480   for (field = 0; field < numFields; ++field) {
481     PetscInt Nc;
482 
483     ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr);
484     numComponents += Nc;
485   }
486   ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
487   ierr = PetscMalloc5(numComponents,&funcVal,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
488   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
489   ierr = PetscFEGetQuadrature(fe[0], &quad);CHKERRQ(ierr);
490   for (c = cStart; c < cEnd; ++c) {
491     PetscScalar *x = NULL;
492     PetscReal    elemDiff = 0.0;
493 
494     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
495     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
496     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
497 
498     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
499       void * const     ctx = ctxs ? ctxs[field] : NULL;
500       const PetscReal *quadPoints, *quadWeights;
501       PetscReal       *basis;
502       PetscInt         numQuadPoints, numBasisFuncs, numBasisComps, q, d, e, fc, f;
503 
504       ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
505       ierr = PetscFEGetDimension(fe[field], &numBasisFuncs);CHKERRQ(ierr);
506       ierr = PetscFEGetNumComponents(fe[field], &numBasisComps);CHKERRQ(ierr);
507       ierr = PetscFEGetDefaultTabulation(fe[field], &basis, NULL, NULL);CHKERRQ(ierr);
508       if (debug) {
509         char title[1024];
510         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
511         ierr = DMPrintCellVector(c, title, numBasisFuncs*numBasisComps, &x[fieldOffset]);CHKERRQ(ierr);
512       }
513       for (q = 0; q < numQuadPoints; ++q) {
514         for (d = 0; d < dim; d++) {
515           coords[d] = v0[d];
516           for (e = 0; e < dim; e++) {
517             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
518           }
519         }
520         (*funcs[field])(coords, funcVal, ctx);
521         for (fc = 0; fc < numBasisComps; ++fc) {
522           PetscScalar interpolant = 0.0;
523 
524           for (f = 0; f < numBasisFuncs; ++f) {
525             const PetscInt fidx = f*numBasisComps+fc;
526             interpolant += x[fieldOffset+fidx]*basis[q*numBasisFuncs*numBasisComps+fidx];
527           }
528           if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "    elem %d field %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ);CHKERRQ(ierr);}
529           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ;
530         }
531       }
532       comp        += numBasisComps;
533       fieldOffset += numBasisFuncs*numBasisComps;
534     }
535     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
536     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
537     localDiff += elemDiff;
538   }
539   ierr  = PetscFree5(funcVal,coords,v0,J,invJ);CHKERRQ(ierr);
540   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
541   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
542   *diff = PetscSqrtReal(*diff);
543   PetscFunctionReturn(0);
544 }
545 
546 #undef __FUNCT__
547 #define __FUNCT__ "DMPlexComputeL2GradientDiff"
548 /*@C
549   DMPlexComputeL2GradientDiff - This function computes the L_2 difference between the gradient of a function u and an FEM interpolant solution grad u_h.
550 
551   Input Parameters:
552 + dm    - The DM
553 . fe    - The PetscFE object for each field
554 . funcs - The gradient functions to evaluate for each field component
555 . ctxs  - Optional array of contexts to pass to each function, or NULL.
556 . X     - The coefficient vector u_h
557 - n     - The vector to project along
558 
559   Output Parameter:
560 . diff - The diff ||(grad u - grad u_h) . n||_2
561 
562   Level: developer
563 
564 .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff()
565 @*/
566 PetscErrorCode DMPlexComputeL2GradientDiff(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, const PetscReal n[], PetscReal *diff)
567 {
568   const PetscInt  debug = 0;
569   PetscSection    section;
570   PetscQuadrature quad;
571   Vec             localX;
572   PetscScalar    *funcVal, *interpolantVec;
573   PetscReal      *coords, *realSpaceDer, *v0, *J, *invJ, detJ;
574   PetscReal       localDiff = 0.0;
575   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp;
576   PetscErrorCode  ierr;
577 
578   PetscFunctionBegin;
579   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
580   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
581   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
582   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
583   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
584   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
585   for (field = 0; field < numFields; ++field) {
586     PetscInt Nc;
587 
588     ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr);
589     numComponents += Nc;
590   }
591   /* ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */
592   ierr = PetscMalloc7(numComponents,&funcVal,dim,&coords,dim,&realSpaceDer,dim,&v0,dim*dim,&J,dim*dim,&invJ,dim,&interpolantVec);CHKERRQ(ierr);
593   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
594   ierr = PetscFEGetQuadrature(fe[0], &quad);CHKERRQ(ierr);
595   for (c = cStart; c < cEnd; ++c) {
596     PetscScalar *x = NULL;
597     PetscReal    elemDiff = 0.0;
598 
599     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
600     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
601     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
602 
603     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
604       void * const     ctx = ctxs ? ctxs[field] : NULL;
605       const PetscReal *quadPoints, *quadWeights;
606       PetscReal       *basisDer;
607       PetscInt         numQuadPoints, Nb, Ncomp, q, d, e, fc, f, g;
608 
609       ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
610       ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr);
611       ierr = PetscFEGetNumComponents(fe[field], &Ncomp);CHKERRQ(ierr);
612       ierr = PetscFEGetDefaultTabulation(fe[field], NULL, &basisDer, NULL);CHKERRQ(ierr);
613       if (debug) {
614         char title[1024];
615         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
616         ierr = DMPrintCellVector(c, title, Nb*Ncomp, &x[fieldOffset]);CHKERRQ(ierr);
617       }
618       for (q = 0; q < numQuadPoints; ++q) {
619         for (d = 0; d < dim; d++) {
620           coords[d] = v0[d];
621           for (e = 0; e < dim; e++) {
622             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
623           }
624         }
625         (*funcs[field])(coords, n, funcVal, ctx);
626         for (fc = 0; fc < Ncomp; ++fc) {
627           PetscScalar interpolant = 0.0;
628 
629           for (d = 0; d < dim; ++d) interpolantVec[d] = 0.0;
630           for (f = 0; f < Nb; ++f) {
631             const PetscInt fidx = f*Ncomp+fc;
632 
633             for (d = 0; d < dim; ++d) {
634               realSpaceDer[d] = 0.0;
635               for (g = 0; g < dim; ++g) {
636                 realSpaceDer[d] += invJ[g*dim+d]*basisDer[(q*Nb*Ncomp+fidx)*dim+g];
637               }
638               interpolantVec[d] += x[fieldOffset+fidx]*realSpaceDer[d];
639             }
640           }
641           for (d = 0; d < dim; ++d) interpolant += interpolantVec[d]*n[d];
642           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);}
643           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ;
644         }
645       }
646       comp        += Ncomp;
647       fieldOffset += Nb*Ncomp;
648     }
649     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
650     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
651     localDiff += elemDiff;
652   }
653   ierr  = PetscFree7(funcVal,coords,realSpaceDer,v0,J,invJ,interpolantVec);CHKERRQ(ierr);
654   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
655   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
656   *diff = PetscSqrtReal(*diff);
657   PetscFunctionReturn(0);
658 }
659 
660 #undef __FUNCT__
661 #define __FUNCT__ "DMPlexComputeResidualFEM"
662 /*@
663   DMPlexComputeResidualFEM - Form the local residual F from the local input X using pointwise functions specified by the user
664 
665   Input Parameters:
666 + dm - The mesh
667 . X  - Local input vector
668 - user - The user context
669 
670   Output Parameter:
671 . F  - Local output vector
672 
673   Note:
674   The first member of the user context must be an FEMContext.
675 
676   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
677   like a GPU, or vectorize on a multicore machine.
678 
679   Level: developer
680 
681 .seealso: DMPlexComputeJacobianActionFEM()
682 @*/
683 PetscErrorCode DMPlexComputeResidualFEM(DM dm, Vec X, Vec F, void *user)
684 {
685   DM_Plex          *mesh  = (DM_Plex *) dm->data;
686   PetscFEM         *fem   = (PetscFEM *) user;
687   PetscFE          *fe    = fem->fe;
688   PetscFE          *feAux = fem->feAux;
689   PetscFE          *feBd  = fem->feBd;
690   const char       *name  = "Residual";
691   DM                dmAux;
692   Vec               A;
693   PetscQuadrature   q;
694   PetscCellGeometry geom;
695   PetscSection      section, sectionAux;
696   PetscReal        *v0, *J, *invJ, *detJ;
697   PetscScalar      *elemVec, *u, *a = NULL;
698   PetscInt          dim, Nf, NfAux = 0, f, numCells, cStart, cEnd, c;
699   PetscInt          cellDof = 0, numComponents = 0;
700   PetscInt          cellDofAux = 0, numComponentsAux = 0;
701   PetscErrorCode    ierr;
702 
703   PetscFunctionBegin;
704   ierr = PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
705   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
706   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
707   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
708   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
709   numCells = cEnd - cStart;
710   for (f = 0; f < Nf; ++f) {
711     PetscInt Nb, Nc;
712 
713     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
714     ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr);
715     cellDof       += Nb*Nc;
716     numComponents += Nc;
717   }
718   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
719   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
720   if (dmAux) {
721     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
722     ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr);
723   }
724   for (f = 0; f < NfAux; ++f) {
725     PetscInt Nb, Nc;
726 
727     ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr);
728     ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr);
729     cellDofAux       += Nb*Nc;
730     numComponentsAux += Nc;
731   }
732   ierr = DMPlexInsertBoundaryValuesFEM(dm, X);CHKERRQ(ierr);
733   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
734   ierr = PetscMalloc6(numCells*cellDof,&u,numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*cellDof,&elemVec);CHKERRQ(ierr);
735   if (dmAux) {ierr = PetscMalloc1(numCells*cellDofAux, &a);CHKERRQ(ierr);}
736   for (c = cStart; c < cEnd; ++c) {
737     PetscScalar *x = NULL;
738     PetscInt     i;
739 
740     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
741     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
742     ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
743     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
744     ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
745     if (dmAux) {
746       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
747       for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x[i];
748       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
749     }
750   }
751   for (f = 0; f < Nf; ++f) {
752     void   (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f0Funcs[f];
753     void   (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f1Funcs[f];
754     PetscInt numQuadPoints, Nb;
755     /* Conforming batches */
756     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
757     /* Remainder */
758     PetscInt Nr, offset;
759 
760     ierr = PetscFEGetQuadrature(fe[f], &q);CHKERRQ(ierr);
761     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
762     ierr = PetscFEGetTileSizes(fe[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
763     ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
764     blockSize = Nb*numQuadPoints;
765     batchSize = numBlocks * blockSize;
766     ierr =  PetscFESetTileSizes(fe[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
767     numChunks = numCells / (numBatches*batchSize);
768     Ne        = numChunks*numBatches*batchSize;
769     Nr        = numCells % (numBatches*batchSize);
770     offset    = numCells - Nr;
771     geom.v0   = v0;
772     geom.J    = J;
773     geom.invJ = invJ;
774     geom.detJ = detJ;
775     ierr = PetscFEIntegrateResidual(fe[f], Ne, Nf, fe, f, geom, u, NfAux, feAux, a, f0, f1, elemVec);CHKERRQ(ierr);
776     geom.v0   = &v0[offset*dim];
777     geom.J    = &J[offset*dim*dim];
778     geom.invJ = &invJ[offset*dim*dim];
779     geom.detJ = &detJ[offset];
780     ierr = PetscFEIntegrateResidual(fe[f], Nr, Nf, fe, f, geom, &u[offset*cellDof], NfAux, feAux, &a[offset*cellDofAux], f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr);
781   }
782   for (c = cStart; c < cEnd; ++c) {
783     if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, name, cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);}
784     ierr = DMPlexVecSetClosure(dm, section, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr);
785   }
786   ierr = PetscFree6(u,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr);
787   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
788   if (feBd) {
789     DMLabel  depth;
790     PetscInt numBd, bd;
791 
792     for (f = 0, cellDof = 0, numComponents = 0; f < Nf; ++f) {
793       PetscInt Nb, Nc;
794 
795       ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr);
796       ierr = PetscFEGetNumComponents(feBd[f], &Nc);CHKERRQ(ierr);
797       cellDof       += Nb*Nc;
798       numComponents += Nc;
799     }
800     ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr);
801     ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
802     for (bd = 0; bd < numBd; ++bd) {
803       const char     *bdLabel;
804       DMLabel         label;
805       IS              pointIS;
806       const PetscInt *points;
807       const PetscInt *values;
808       PetscReal      *n;
809       PetscInt        field, numValues, numPoints, p, dep, numFaces;
810       PetscBool       isEssential;
811 
812       ierr = DMPlexGetBoundary(dm, bd, &isEssential, &bdLabel, &field, NULL, &numValues, &values, NULL);CHKERRQ(ierr);
813       if (numValues != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Bug me and I will fix this");
814       ierr = DMPlexGetLabel(dm, bdLabel, &label);CHKERRQ(ierr);
815       ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr);
816       ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr);
817       ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
818       for (p = 0, numFaces = 0; p < numPoints; ++p) {
819         ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
820         if (dep == dim-1) ++numFaces;
821       }
822       ierr = PetscMalloc7(numFaces*cellDof,&u,numFaces*dim,&v0,numFaces*dim,&n,numFaces*dim*dim,&J,numFaces*dim*dim,&invJ,numFaces,&detJ,numFaces*cellDof,&elemVec);CHKERRQ(ierr);
823       for (p = 0, f = 0; p < numPoints; ++p) {
824         const PetscInt point = points[p];
825         PetscScalar   *x     = NULL;
826         PetscInt       i;
827 
828         ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
829         if (dep != dim-1) continue;
830         ierr = DMPlexComputeCellGeometry(dm, point, &v0[f*dim], &J[f*dim*dim], &invJ[f*dim*dim], &detJ[f]);CHKERRQ(ierr);
831         ierr = DMPlexComputeCellGeometryFVM(dm, point, NULL, NULL, &n[f*dim]);
832         if (detJ[f] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[f], point);
833         ierr = DMPlexVecGetClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
834         for (i = 0; i < cellDof; ++i) u[f*cellDof+i] = x[i];
835         ierr = DMPlexVecRestoreClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
836         ++f;
837       }
838       for (f = 0; f < Nf; ++f) {
839         void   (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f0BdFuncs[f];
840         void   (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f1BdFuncs[f];
841         PetscInt numQuadPoints, Nb;
842         /* Conforming batches */
843         PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
844         /* Remainder */
845         PetscInt Nr, offset;
846 
847         ierr = PetscFEGetQuadrature(feBd[f], &q);CHKERRQ(ierr);
848         ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr);
849         ierr = PetscFEGetTileSizes(feBd[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
850         ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
851         blockSize = Nb*numQuadPoints;
852         batchSize = numBlocks * blockSize;
853         ierr =  PetscFESetTileSizes(feBd[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
854         numChunks = numFaces / (numBatches*batchSize);
855         Ne        = numChunks*numBatches*batchSize;
856         Nr        = numFaces % (numBatches*batchSize);
857         offset    = numFaces - Nr;
858         geom.v0   = v0;
859         geom.n    = n;
860         geom.J    = J;
861         geom.invJ = invJ;
862         geom.detJ = detJ;
863         ierr = PetscFEIntegrateBdResidual(feBd[f], Ne, Nf, feBd, f, geom, u, 0, NULL, NULL, f0, f1, elemVec);CHKERRQ(ierr);
864         geom.v0   = &v0[offset*dim];
865         geom.n    = &n[offset*dim];
866         geom.J    = &J[offset*dim*dim];
867         geom.invJ = &invJ[offset*dim*dim];
868         geom.detJ = &detJ[offset];
869         ierr = PetscFEIntegrateBdResidual(feBd[f], Nr, Nf, feBd, f, geom, &u[offset*cellDof], 0, NULL, NULL, f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr);
870       }
871       for (p = 0, f = 0; p < numPoints; ++p) {
872         const PetscInt point = points[p];
873 
874         ierr = DMLabelGetValue(depth, point, &dep);CHKERRQ(ierr);
875         if (dep != dim-1) continue;
876         if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "BdResidual", cellDof, &elemVec[f*cellDof]);CHKERRQ(ierr);}
877         ierr = DMPlexVecSetClosure(dm, NULL, F, point, &elemVec[f*cellDof], ADD_VALUES);CHKERRQ(ierr);
878         ++f;
879       }
880       ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
881       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
882       ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemVec);CHKERRQ(ierr);
883     }
884   }
885   if (mesh->printFEM) {ierr = DMPrintLocalVec(dm, name, mesh->printTol, F);CHKERRQ(ierr);}
886   ierr = PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
887   PetscFunctionReturn(0);
888 }
889 
890 #undef __FUNCT__
891 #define __FUNCT__ "DMPlexComputeIFunctionFEM"
892 /*@
893   DMPlexComputeIFunctionFEM - Form the local implicit function F from the local input X, X_t using pointwise functions specified by the user
894 
895   Input Parameters:
896 + dm - The mesh
897 . time - The current time
898 . X  - Local input vector
899 . X_t  - Time derivative of the local input vector
900 - user - The user context
901 
902   Output Parameter:
903 . F  - Local output vector
904 
905   Note:
906   The first member of the user context must be an FEMContext.
907 
908   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
909   like a GPU, or vectorize on a multicore machine.
910 
911   Level: developer
912 
913 .seealso: DMPlexComputeResidualFEM()
914 @*/
915 PetscErrorCode DMPlexComputeIFunctionFEM(DM dm, PetscReal time, Vec X, Vec X_t, Vec F, void *user)
916 {
917   DM_Plex          *mesh  = (DM_Plex *) dm->data;
918   PetscFEM         *fem   = (PetscFEM *) user;
919   PetscFE          *fe    = fem->fe;
920   PetscFE          *feAux = fem->feAux;
921   PetscFE          *feBd  = fem->feBd;
922   const char       *name  = "Residual";
923   DM                dmAux;
924   Vec               A;
925   PetscQuadrature   q;
926   PetscCellGeometry geom;
927   PetscSection      section, sectionAux;
928   PetscReal        *v0, *J, *invJ, *detJ;
929   PetscScalar      *elemVec, *u, *u_t, *a = NULL;
930   PetscInt          dim, Nf, NfAux = 0, f, numCells, cStart, cEnd, c;
931   PetscInt          cellDof = 0, numComponents = 0;
932   PetscInt          cellDofAux = 0, numComponentsAux = 0;
933   PetscErrorCode    ierr;
934 
935   PetscFunctionBegin;
936   ierr = PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
937   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
938   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
939   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
940   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
941   numCells = cEnd - cStart;
942   for (f = 0; f < Nf; ++f) {
943     PetscInt Nb, Nc;
944 
945     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
946     ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr);
947     cellDof       += Nb*Nc;
948     numComponents += Nc;
949   }
950   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
951   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
952   if (dmAux) {
953     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
954     ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr);
955   }
956   for (f = 0; f < NfAux; ++f) {
957     PetscInt Nb, Nc;
958 
959     ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr);
960     ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr);
961     cellDofAux       += Nb*Nc;
962     numComponentsAux += Nc;
963   }
964   ierr = DMPlexInsertBoundaryValuesFEM(dm, X);CHKERRQ(ierr);
965   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
966   ierr = PetscMalloc7(numCells*cellDof,&u,numCells*cellDof,&u_t,numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*cellDof,&elemVec);CHKERRQ(ierr);
967   if (dmAux) {ierr = PetscMalloc1(numCells*cellDofAux, &a);CHKERRQ(ierr);}
968   for (c = cStart; c < cEnd; ++c) {
969     PetscScalar *x = NULL, *x_t = NULL;
970     PetscInt     i;
971 
972     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
973     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
974     ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
975     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
976     ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
977     ierr = DMPlexVecGetClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr);
978     for (i = 0; i < cellDof; ++i) u_t[c*cellDof+i] = x_t[i];
979     ierr = DMPlexVecRestoreClosure(dm, section, X_t, c, NULL, &x_t);CHKERRQ(ierr);
980     if (dmAux) {
981       PetscScalar *x_a = NULL;
982       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x_a);CHKERRQ(ierr);
983       for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x_a[i];
984       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x_a);CHKERRQ(ierr);
985     }
986   }
987   for (f = 0; f < Nf; ++f) {
988     void   (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f0IFuncs[f];
989     void   (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f1IFuncs[f];
990     PetscInt numQuadPoints, Nb;
991     /* Conforming batches */
992     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
993     /* Remainder */
994     PetscInt Nr, offset;
995 
996     ierr = PetscFEGetQuadrature(fe[f], &q);CHKERRQ(ierr);
997     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
998     ierr = PetscFEGetTileSizes(fe[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
999     ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
1000     blockSize = Nb*numQuadPoints;
1001     batchSize = numBlocks * blockSize;
1002     ierr =  PetscFESetTileSizes(fe[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1003     numChunks = numCells / (numBatches*batchSize);
1004     Ne        = numChunks*numBatches*batchSize;
1005     Nr        = numCells % (numBatches*batchSize);
1006     offset    = numCells - Nr;
1007     geom.v0   = v0;
1008     geom.J    = J;
1009     geom.invJ = invJ;
1010     geom.detJ = detJ;
1011     ierr = PetscFEIntegrateIFunction(fe[f], Ne, Nf, fe, f, geom, u, u_t, NfAux, feAux, a, f0, f1, elemVec);CHKERRQ(ierr);
1012     geom.v0   = &v0[offset*dim];
1013     geom.J    = &J[offset*dim*dim];
1014     geom.invJ = &invJ[offset*dim*dim];
1015     geom.detJ = &detJ[offset];
1016     ierr = PetscFEIntegrateIFunction(fe[f], Nr, Nf, fe, f, geom, &u[offset*cellDof], &u_t[offset*cellDof], NfAux, feAux, &a[offset*cellDofAux], f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr);
1017   }
1018   for (c = cStart; c < cEnd; ++c) {
1019     if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, name, cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);}
1020     ierr = DMPlexVecSetClosure(dm, section, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr);
1021   }
1022   ierr = PetscFree7(u,u_t,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr);
1023   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
1024   if (feBd) {
1025     DMLabel         label, depth;
1026     IS              pointIS;
1027     const PetscInt *points;
1028     PetscInt        dep, numPoints, p, numFaces;
1029     PetscReal      *n;
1030 
1031     ierr = DMPlexGetLabel(dm, "boundary", &label);CHKERRQ(ierr);
1032     ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr);
1033     ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr);
1034     ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr);
1035     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
1036     for (f = 0, cellDof = 0, numComponents = 0; f < Nf; ++f) {
1037       PetscInt Nb, Nc;
1038 
1039       ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr);
1040       ierr = PetscFEGetNumComponents(feBd[f], &Nc);CHKERRQ(ierr);
1041       cellDof       += Nb*Nc;
1042       numComponents += Nc;
1043     }
1044     for (p = 0, numFaces = 0; p < numPoints; ++p) {
1045       ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
1046       if (dep == dim-1) ++numFaces;
1047     }
1048     ierr = PetscMalloc7(numFaces*cellDof,&u,numFaces*dim,&v0,numFaces*dim,&n,numFaces*dim*dim,&J,numFaces*dim*dim,&invJ,numFaces,&detJ,numFaces*cellDof,&elemVec);CHKERRQ(ierr);
1049     ierr = PetscMalloc1(numFaces*cellDof,&u_t);CHKERRQ(ierr);
1050     for (p = 0, f = 0; p < numPoints; ++p) {
1051       const PetscInt point = points[p];
1052       PetscScalar   *x     = NULL;
1053       PetscInt       i;
1054 
1055       ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
1056       if (dep != dim-1) continue;
1057       ierr = DMPlexComputeCellGeometry(dm, point, &v0[f*dim], &J[f*dim*dim], &invJ[f*dim*dim], &detJ[f]);CHKERRQ(ierr);
1058       ierr = DMPlexComputeCellGeometryFVM(dm, point, NULL, NULL, &n[f*dim]);
1059       if (detJ[f] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[f], point);
1060       ierr = DMPlexVecGetClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
1061       for (i = 0; i < cellDof; ++i) u[f*cellDof+i] = x[i];
1062       ierr = DMPlexVecRestoreClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
1063       ierr = DMPlexVecGetClosure(dm, section, X_t, point, NULL, &x);CHKERRQ(ierr);
1064       for (i = 0; i < cellDof; ++i) u_t[f*cellDof+i] = x[i];
1065       ierr = DMPlexVecRestoreClosure(dm, section, X_t, point, NULL, &x);CHKERRQ(ierr);
1066       ++f;
1067     }
1068     for (f = 0; f < Nf; ++f) {
1069       void   (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f0BdIFuncs[f];
1070       void   (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f1BdIFuncs[f];
1071       PetscInt numQuadPoints, Nb;
1072       /* Conforming batches */
1073       PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1074       /* Remainder */
1075       PetscInt Nr, offset;
1076 
1077       ierr = PetscFEGetQuadrature(feBd[f], &q);CHKERRQ(ierr);
1078       ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr);
1079       ierr = PetscFEGetTileSizes(feBd[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
1080       ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
1081       blockSize = Nb*numQuadPoints;
1082       batchSize = numBlocks * blockSize;
1083       ierr =  PetscFESetTileSizes(feBd[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1084       numChunks = numFaces / (numBatches*batchSize);
1085       Ne        = numChunks*numBatches*batchSize;
1086       Nr        = numFaces % (numBatches*batchSize);
1087       offset    = numFaces - Nr;
1088       geom.v0   = v0;
1089       geom.n    = n;
1090       geom.J    = J;
1091       geom.invJ = invJ;
1092       geom.detJ = detJ;
1093       ierr = PetscFEIntegrateBdIFunction(feBd[f], Ne, Nf, feBd, f, geom, u, u_t, 0, NULL, NULL, f0, f1, elemVec);CHKERRQ(ierr);
1094       geom.v0   = &v0[offset*dim];
1095       geom.n    = &n[offset*dim];
1096       geom.J    = &J[offset*dim*dim];
1097       geom.invJ = &invJ[offset*dim*dim];
1098       geom.detJ = &detJ[offset];
1099       ierr = PetscFEIntegrateBdIFunction(feBd[f], Nr, Nf, feBd, f, geom, &u[offset*cellDof], &u_t[offset*cellDof], 0, NULL, NULL, f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr);
1100     }
1101     for (p = 0, f = 0; p < numPoints; ++p) {
1102       const PetscInt point = points[p];
1103 
1104       ierr = DMLabelGetValue(depth, point, &dep);CHKERRQ(ierr);
1105       if (dep != dim-1) continue;
1106       if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "BdResidual", cellDof, &elemVec[f*cellDof]);CHKERRQ(ierr);}
1107       ierr = DMPlexVecSetClosure(dm, NULL, F, point, &elemVec[f*cellDof], ADD_VALUES);CHKERRQ(ierr);
1108       ++f;
1109     }
1110     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
1111     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
1112     ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemVec);CHKERRQ(ierr);
1113     ierr = PetscFree(u_t);CHKERRQ(ierr);
1114   }
1115   if (mesh->printFEM) {ierr = DMPrintLocalVec(dm, name, mesh->printTol, F);CHKERRQ(ierr);}
1116   ierr = PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
1117   PetscFunctionReturn(0);
1118 }
1119 
1120 #undef __FUNCT__
1121 #define __FUNCT__ "DMPlexComputeJacobianActionFEM"
1122 /*@C
1123   DMPlexComputeJacobianActionFEM - Form the local action of Jacobian J(u) on the local input X using pointwise functions specified by the user
1124 
1125   Input Parameters:
1126 + dm - The mesh
1127 . J  - The Jacobian shell matrix
1128 . X  - Local input vector
1129 - user - The user context
1130 
1131   Output Parameter:
1132 . F  - Local output vector
1133 
1134   Note:
1135   The first member of the user context must be an FEMContext.
1136 
1137   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
1138   like a GPU, or vectorize on a multicore machine.
1139 
1140   Level: developer
1141 
1142 .seealso: DMPlexComputeResidualFEM()
1143 @*/
1144 PetscErrorCode DMPlexComputeJacobianActionFEM(DM dm, Mat Jac, Vec X, Vec F, void *user)
1145 {
1146   DM_Plex          *mesh = (DM_Plex *) dm->data;
1147   PetscFEM         *fem  = (PetscFEM *) user;
1148   PetscFE          *fe   = fem->fe;
1149   PetscQuadrature   quad;
1150   PetscCellGeometry geom;
1151   PetscSection      section;
1152   JacActionCtx     *jctx;
1153   PetscReal        *v0, *J, *invJ, *detJ;
1154   PetscScalar      *elemVec, *u, *a;
1155   PetscInt          dim, numFields, field, numCells, cStart, cEnd, c;
1156   PetscInt          cellDof = 0;
1157   PetscErrorCode    ierr;
1158 
1159   PetscFunctionBegin;
1160   /* ierr = PetscLogEventBegin(DMPLEX_JacobianActionFEM,dm,0,0,0);CHKERRQ(ierr); */
1161   ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr);
1162   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
1163   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1164   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
1165   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1166   numCells = cEnd - cStart;
1167   for (field = 0; field < numFields; ++field) {
1168     PetscInt Nb, Nc;
1169 
1170     ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr);
1171     ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr);
1172     cellDof += Nb*Nc;
1173   }
1174   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
1175   ierr = PetscMalloc7(numCells*cellDof,&u,numCells*cellDof,&a,numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*cellDof,&elemVec);CHKERRQ(ierr);
1176   for (c = cStart; c < cEnd; ++c) {
1177     PetscScalar *x = NULL;
1178     PetscInt     i;
1179 
1180     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
1181     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
1182     ierr = DMPlexVecGetClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr);
1183     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
1184     ierr = DMPlexVecRestoreClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr);
1185     ierr = DMPlexVecGetClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr);
1186     for (i = 0; i < cellDof; ++i) a[c*cellDof+i] = x[i];
1187     ierr = DMPlexVecRestoreClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr);
1188   }
1189   for (field = 0; field < numFields; ++field) {
1190     PetscInt numQuadPoints, Nb;
1191     /* Conforming batches */
1192     PetscInt numBlocks  = 1;
1193     PetscInt numBatches = 1;
1194     PetscInt numChunks, Ne, blockSize, batchSize;
1195     /* Remainder */
1196     PetscInt Nr, offset;
1197 
1198     ierr = PetscFEGetQuadrature(fe[field], &quad);CHKERRQ(ierr);
1199     ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr);
1200     ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
1201     blockSize = Nb*numQuadPoints;
1202     batchSize = numBlocks * blockSize;
1203     numChunks = numCells / (numBatches*batchSize);
1204     Ne        = numChunks*numBatches*batchSize;
1205     Nr        = numCells % (numBatches*batchSize);
1206     offset    = numCells - Nr;
1207     geom.v0   = v0;
1208     geom.J    = J;
1209     geom.invJ = invJ;
1210     geom.detJ = detJ;
1211     ierr = PetscFEIntegrateJacobianAction(fe[field], Ne, numFields, fe, field, geom, u, a, fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, elemVec);CHKERRQ(ierr);
1212     geom.v0   = &v0[offset*dim];
1213     geom.J    = &J[offset*dim*dim];
1214     geom.invJ = &invJ[offset*dim*dim];
1215     geom.detJ = &detJ[offset];
1216     ierr = PetscFEIntegrateJacobianAction(fe[field], Nr, numFields, fe, field, geom, &u[offset*cellDof], &a[offset*cellDof],
1217                                           fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, &elemVec[offset*cellDof]);CHKERRQ(ierr);
1218   }
1219   for (c = cStart; c < cEnd; ++c) {
1220     if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, "Jacobian Action", cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);}
1221     ierr = DMPlexVecSetClosure(dm, NULL, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr);
1222   }
1223   ierr = PetscFree7(u,a,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr);
1224   if (mesh->printFEM) {
1225     PetscMPIInt rank, numProcs;
1226     PetscInt    p;
1227 
1228     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
1229     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm), &numProcs);CHKERRQ(ierr);
1230     ierr = PetscPrintf(PetscObjectComm((PetscObject)dm), "Jacobian Action:\n");CHKERRQ(ierr);
1231     for (p = 0; p < numProcs; ++p) {
1232       if (p == rank) {ierr = VecView(F, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);}
1233       ierr = PetscBarrier((PetscObject) dm);CHKERRQ(ierr);
1234     }
1235   }
1236   /* ierr = PetscLogEventEnd(DMPLEX_JacobianActionFEM,dm,0,0,0);CHKERRQ(ierr); */
1237   PetscFunctionReturn(0);
1238 }
1239 
1240 #undef __FUNCT__
1241 #define __FUNCT__ "DMPlexComputeJacobianFEM"
1242 /*@
1243   DMPlexComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user.
1244 
1245   Input Parameters:
1246 + dm - The mesh
1247 . X  - Local input vector
1248 - user - The user context
1249 
1250   Output Parameter:
1251 . Jac  - Jacobian matrix
1252 
1253   Note:
1254   The first member of the user context must be an FEMContext.
1255 
1256   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
1257   like a GPU, or vectorize on a multicore machine.
1258 
1259   Level: developer
1260 
1261 .seealso: FormFunctionLocal()
1262 @*/
1263 PetscErrorCode DMPlexComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP,void *user)
1264 {
1265   DM_Plex          *mesh  = (DM_Plex *) dm->data;
1266   PetscFEM         *fem   = (PetscFEM *) user;
1267   PetscFE          *fe    = fem->fe;
1268   PetscFE          *feAux = fem->feAux;
1269   PetscFE          *feBd  = fem->feBd;
1270   const char       *name  = "Jacobian";
1271   DM                dmAux;
1272   Vec               A;
1273   PetscQuadrature   quad;
1274   PetscCellGeometry geom;
1275   PetscSection      section, globalSection, sectionAux;
1276   PetscReal        *v0, *J, *invJ, *detJ;
1277   PetscScalar      *elemMat, *u, *a;
1278   PetscInt          dim, Nf, NfAux = 0, f, fieldI, fieldJ, numCells, cStart, cEnd, c;
1279   PetscInt          cellDof = 0, numComponents = 0;
1280   PetscInt          cellDofAux = 0, numComponentsAux = 0;
1281   PetscBool         isShell;
1282   PetscErrorCode    ierr;
1283 
1284   PetscFunctionBegin;
1285   ierr = PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
1286   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
1287   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1288   ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);
1289   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
1290   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1291   numCells = cEnd - cStart;
1292   for (f = 0; f < Nf; ++f) {
1293     PetscInt Nb, Nc;
1294 
1295     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
1296     ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr);
1297     cellDof       += Nb*Nc;
1298     numComponents += Nc;
1299   }
1300   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
1301   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
1302   if (dmAux) {
1303     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
1304     ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr);
1305   }
1306   for (f = 0; f < NfAux; ++f) {
1307     PetscInt Nb, Nc;
1308 
1309     ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr);
1310     ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr);
1311     cellDofAux       += Nb*Nc;
1312     numComponentsAux += Nc;
1313   }
1314   ierr = DMPlexInsertBoundaryValuesFEM(dm, X);CHKERRQ(ierr);
1315   ierr = MatZeroEntries(JacP);CHKERRQ(ierr);
1316   ierr = PetscMalloc6(numCells*cellDof,&u,numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*cellDof*cellDof,&elemMat);CHKERRQ(ierr);
1317   if (dmAux) {ierr = PetscMalloc1(numCells*cellDofAux, &a);CHKERRQ(ierr);}
1318   for (c = cStart; c < cEnd; ++c) {
1319     PetscScalar *x = NULL;
1320     PetscInt     i;
1321 
1322     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
1323     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
1324     ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
1325     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
1326     ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
1327     if (dmAux) {
1328       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
1329       for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x[i];
1330       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
1331     }
1332   }
1333   ierr = PetscMemzero(elemMat, numCells*cellDof*cellDof * sizeof(PetscScalar));CHKERRQ(ierr);
1334   for (fieldI = 0; fieldI < Nf; ++fieldI) {
1335     PetscInt numQuadPoints, Nb;
1336     /* Conforming batches */
1337     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1338     /* Remainder */
1339     PetscInt Nr, offset;
1340 
1341     ierr = PetscFEGetQuadrature(fe[fieldI], &quad);CHKERRQ(ierr);
1342     ierr = PetscFEGetDimension(fe[fieldI], &Nb);CHKERRQ(ierr);
1343     ierr = PetscFEGetTileSizes(fe[fieldI], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
1344     ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
1345     blockSize = Nb*numQuadPoints;
1346     batchSize = numBlocks * blockSize;
1347     ierr = PetscFESetTileSizes(fe[fieldI], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1348     numChunks = numCells / (numBatches*batchSize);
1349     Ne        = numChunks*numBatches*batchSize;
1350     Nr        = numCells % (numBatches*batchSize);
1351     offset    = numCells - Nr;
1352     for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
1353       void   (*g0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g0Funcs[fieldI*Nf+fieldJ];
1354       void   (*g1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g1Funcs[fieldI*Nf+fieldJ];
1355       void   (*g2)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g2Funcs[fieldI*Nf+fieldJ];
1356       void   (*g3)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g3Funcs[fieldI*Nf+fieldJ];
1357 
1358       geom.v0   = v0;
1359       geom.J    = J;
1360       geom.invJ = invJ;
1361       geom.detJ = detJ;
1362       ierr = PetscFEIntegrateJacobian(fe[fieldI], Ne, Nf, fe, fieldI, fieldJ, geom, u, NfAux, feAux, a, g0, g1, g2, g3, elemMat);CHKERRQ(ierr);
1363       geom.v0   = &v0[offset*dim];
1364       geom.J    = &J[offset*dim*dim];
1365       geom.invJ = &invJ[offset*dim*dim];
1366       geom.detJ = &detJ[offset];
1367       ierr = PetscFEIntegrateJacobian(fe[fieldI], Nr, Nf, fe, fieldI, fieldJ, geom, &u[offset*cellDof], NfAux, feAux, &a[offset*cellDofAux], g0, g1, g2, g3, &elemMat[offset*cellDof*cellDof]);CHKERRQ(ierr);
1368     }
1369   }
1370   for (c = cStart; c < cEnd; ++c) {
1371     if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, name, cellDof, cellDof, &elemMat[c*cellDof*cellDof]);CHKERRQ(ierr);}
1372     ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, c, &elemMat[c*cellDof*cellDof], ADD_VALUES);CHKERRQ(ierr);
1373   }
1374   ierr = PetscFree6(u,v0,J,invJ,detJ,elemMat);CHKERRQ(ierr);
1375   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
1376   if (feBd) {
1377     DMLabel  depth;
1378     PetscInt numBd, bd;
1379 
1380     for (f = 0, cellDof = 0, numComponents = 0; f < Nf; ++f) {
1381       PetscInt Nb, Nc;
1382 
1383       ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr);
1384       ierr = PetscFEGetNumComponents(feBd[f], &Nc);CHKERRQ(ierr);
1385       cellDof       += Nb*Nc;
1386       numComponents += Nc;
1387     }
1388     ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr);
1389     ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
1390     for (bd = 0; bd < numBd; ++bd) {
1391       const char     *bdLabel;
1392       DMLabel         label;
1393       IS              pointIS;
1394       const PetscInt *points;
1395       const PetscInt *values;
1396       PetscReal      *n;
1397       PetscInt        field, numValues, numPoints, p, dep, numFaces;
1398       PetscBool       isEssential;
1399 
1400       ierr = DMPlexGetBoundary(dm, bd, &isEssential, &bdLabel, &field, NULL, &numValues, &values, NULL);CHKERRQ(ierr);
1401       if (numValues != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Bug me and I will fix this");
1402       ierr = DMPlexGetLabel(dm, bdLabel, &label);CHKERRQ(ierr);
1403       ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr);
1404       ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr);
1405       ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
1406       for (p = 0, numFaces = 0; p < numPoints; ++p) {
1407         ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
1408         if (dep == dim-1) ++numFaces;
1409       }
1410       ierr = PetscMalloc7(numFaces*cellDof,&u,numFaces*dim,&v0,numFaces*dim,&n,numFaces*dim*dim,&J,numFaces*dim*dim,&invJ,numFaces,&detJ,numFaces*cellDof*cellDof,&elemMat);CHKERRQ(ierr);
1411       for (p = 0, f = 0; p < numPoints; ++p) {
1412         const PetscInt point = points[p];
1413         PetscScalar   *x     = NULL;
1414         PetscInt       i;
1415 
1416         ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
1417         if (dep != dim-1) continue;
1418         ierr = DMPlexComputeCellGeometry(dm, point, &v0[f*dim], &J[f*dim*dim], &invJ[f*dim*dim], &detJ[f]);CHKERRQ(ierr);
1419         ierr = DMPlexComputeCellGeometryFVM(dm, point, NULL, NULL, &n[f*dim]);
1420         if (detJ[f] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[f], point);
1421         ierr = DMPlexVecGetClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
1422         for (i = 0; i < cellDof; ++i) u[f*cellDof+i] = x[i];
1423         ierr = DMPlexVecRestoreClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
1424         ++f;
1425       }
1426       ierr = PetscMemzero(elemMat, numFaces*cellDof*cellDof * sizeof(PetscScalar));CHKERRQ(ierr);
1427       for (fieldI = 0; fieldI < Nf; ++fieldI) {
1428         PetscInt numQuadPoints, Nb;
1429         /* Conforming batches */
1430         PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1431         /* Remainder */
1432         PetscInt Nr, offset;
1433 
1434         ierr = PetscFEGetQuadrature(feBd[fieldI], &quad);CHKERRQ(ierr);
1435         ierr = PetscFEGetDimension(feBd[fieldI], &Nb);CHKERRQ(ierr);
1436         ierr = PetscFEGetTileSizes(feBd[fieldI], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
1437         ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
1438         blockSize = Nb*numQuadPoints;
1439         batchSize = numBlocks * blockSize;
1440         ierr =  PetscFESetTileSizes(feBd[fieldI], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1441         numChunks = numFaces / (numBatches*batchSize);
1442         Ne        = numChunks*numBatches*batchSize;
1443         Nr        = numFaces % (numBatches*batchSize);
1444         offset    = numFaces - Nr;
1445         for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
1446           void   (*g0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->g0BdFuncs[fieldI*Nf+fieldJ];
1447           void   (*g1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->g1BdFuncs[fieldI*Nf+fieldJ];
1448           void   (*g2)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->g2BdFuncs[fieldI*Nf+fieldJ];
1449           void   (*g3)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->g3BdFuncs[fieldI*Nf+fieldJ];
1450 
1451           geom.v0   = v0;
1452           geom.n    = n;
1453           geom.J    = J;
1454           geom.invJ = invJ;
1455           geom.detJ = detJ;
1456           ierr = PetscFEIntegrateBdJacobian(feBd[fieldI], Ne, Nf, feBd, fieldI, fieldJ, geom, u, 0, NULL, NULL, g0, g1, g2, g3, elemMat);CHKERRQ(ierr);
1457           geom.v0   = &v0[offset*dim];
1458           geom.n    = &n[offset*dim];
1459           geom.J    = &J[offset*dim*dim];
1460           geom.invJ = &invJ[offset*dim*dim];
1461           geom.detJ = &detJ[offset];
1462           ierr = PetscFEIntegrateBdJacobian(feBd[fieldI], Nr, Nf, feBd, fieldI, fieldJ, geom, &u[offset*cellDof], 0, NULL, NULL, g0, g1, g2, g3, &elemMat[offset*cellDof*cellDof]);CHKERRQ(ierr);
1463         }
1464       }
1465       for (p = 0, f = 0; p < numPoints; ++p) {
1466         const PetscInt point = points[p];
1467 
1468         ierr = DMLabelGetValue(depth, point, &dep);CHKERRQ(ierr);
1469         if (dep != dim-1) continue;
1470         if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(point, "BdJacobian", cellDof, cellDof, &elemMat[f*cellDof*cellDof]);CHKERRQ(ierr);}
1471         ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, point, &elemMat[f*cellDof*cellDof], ADD_VALUES);CHKERRQ(ierr);
1472         ++f;
1473       }
1474       ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
1475       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
1476       ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemMat);CHKERRQ(ierr);
1477     }
1478   }
1479   ierr = MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1480   ierr = MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1481   if (mesh->printFEM) {
1482     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
1483     ierr = MatChop(JacP, 1.0e-10);CHKERRQ(ierr);
1484     ierr = MatView(JacP, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1485   }
1486   ierr = PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
1487   ierr = PetscObjectTypeCompare((PetscObject) Jac, MATSHELL, &isShell);CHKERRQ(ierr);
1488   if (isShell) {
1489     JacActionCtx *jctx;
1490 
1491     ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr);
1492     ierr = VecCopy(X, jctx->u);CHKERRQ(ierr);
1493   }
1494   PetscFunctionReturn(0);
1495 }
1496 
1497 #undef __FUNCT__
1498 #define __FUNCT__ "DMPlexComputeInterpolatorFEM"
1499 /*@
1500   DMPlexComputeInterpolatorFEM - Form the local portion of the interpolation matrix I from the coarse DM to the uniformly refined DM.
1501 
1502   Input Parameters:
1503 + dmf  - The fine mesh
1504 . dmc  - The coarse mesh
1505 - user - The user context
1506 
1507   Output Parameter:
1508 . In  - The interpolation matrix
1509 
1510   Note:
1511   The first member of the user context must be an FEMContext.
1512 
1513   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
1514   like a GPU, or vectorize on a multicore machine.
1515 
1516   Level: developer
1517 
1518 .seealso: DMPlexComputeJacobianFEM()
1519 @*/
1520 PetscErrorCode DMPlexComputeInterpolatorFEM(DM dmc, DM dmf, Mat In, void *user)
1521 {
1522   DM_Plex          *mesh  = (DM_Plex *) dmc->data;
1523   PetscFEM         *fem   = (PetscFEM *) user;
1524   PetscFE          *fe    = fem->fe;
1525   const char       *name  = "Interpolator";
1526   PetscFE          *feRef;
1527   PetscSection      fsection, fglobalSection;
1528   PetscSection      csection, cglobalSection;
1529   PetscScalar      *elemMat;
1530   PetscInt          dim, Nf, f, fieldI, fieldJ, offsetI, offsetJ, cStart, cEnd, c;
1531   PetscInt          rCellDof = 0, cCellDof = 0;
1532   PetscErrorCode    ierr;
1533 
1534   PetscFunctionBegin;
1535 #if 0
1536   ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1537 #endif
1538   ierr = DMPlexGetDimension(dmf, &dim);CHKERRQ(ierr);
1539   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
1540   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
1541   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
1542   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
1543   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
1544   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
1545   ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr);
1546   for (f = 0; f < Nf; ++f) {
1547     PetscInt rNb, cNb, Nc;
1548 
1549     ierr = PetscFERefine(fe[f], &feRef[f]);CHKERRQ(ierr);
1550     ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr);
1551     ierr = PetscFEGetDimension(fe[f], &cNb);CHKERRQ(ierr);
1552     ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr);
1553     rCellDof += rNb*Nc;
1554     cCellDof += cNb*Nc;
1555   }
1556   ierr = MatZeroEntries(In);CHKERRQ(ierr);
1557   ierr = PetscMalloc1(rCellDof*cCellDof,&elemMat);CHKERRQ(ierr);
1558   ierr = PetscMemzero(elemMat, rCellDof*cCellDof * sizeof(PetscScalar));CHKERRQ(ierr);
1559   for (fieldI = 0, offsetI = 0; fieldI < Nf; ++fieldI) {
1560     PetscDualSpace   Qref;
1561     PetscQuadrature  f;
1562     const PetscReal *qpoints, *qweights;
1563     PetscReal       *points;
1564     PetscInt         npoints = 0, Nc, Np, fpdim, i, k, p, d;
1565 
1566     /* Compose points from all dual basis functionals */
1567     ierr = PetscFEGetNumComponents(fe[fieldI], &Nc);CHKERRQ(ierr);
1568     ierr = PetscFEGetDualSpace(feRef[fieldI], &Qref);CHKERRQ(ierr);
1569     ierr = PetscDualSpaceGetDimension(Qref, &fpdim);CHKERRQ(ierr);
1570     for (i = 0; i < fpdim; ++i) {
1571       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1572       ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, NULL);CHKERRQ(ierr);
1573       npoints += Np;
1574     }
1575     ierr = PetscMalloc1(npoints*dim,&points);CHKERRQ(ierr);
1576     for (i = 0, k = 0; i < fpdim; ++i) {
1577       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1578       ierr = PetscQuadratureGetData(f, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr);
1579       for (p = 0; p < Np; ++p, ++k) for (d = 0; d < dim; ++d) points[k*dim+d] = qpoints[p*dim+d];
1580     }
1581 
1582     for (fieldJ = 0, offsetJ = 0; fieldJ < Nf; ++fieldJ) {
1583       PetscReal *B;
1584       PetscInt   NcJ, cpdim, j;
1585 
1586       /* Evaluate basis at points */
1587       ierr = PetscFEGetNumComponents(fe[fieldJ], &NcJ);CHKERRQ(ierr);
1588       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);
1589       ierr = PetscFEGetDimension(fe[fieldJ], &cpdim);CHKERRQ(ierr);
1590       /* For now, fields only interpolate themselves */
1591       if (fieldI == fieldJ) {
1592         ierr = PetscFEGetTabulation(fe[fieldJ], npoints, points, &B, NULL, NULL);CHKERRQ(ierr);
1593         for (i = 0, k = 0; i < fpdim; ++i) {
1594           ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1595           ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, &qweights);CHKERRQ(ierr);
1596           for (p = 0; p < Np; ++p, ++k) {
1597             for (j = 0; j < cpdim; ++j) {
1598               for (c = 0; c < Nc; ++c) elemMat[(offsetI + i*Nc + c)*cCellDof + offsetJ + j*NcJ + c] += B[k*cpdim*NcJ+j*Nc+c]*qweights[p];
1599             }
1600           }
1601         }
1602         ierr = PetscFERestoreTabulation(fe[fieldJ], npoints, points, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr);
1603       }
1604       offsetJ += cpdim*NcJ;
1605     }
1606     offsetI += fpdim*Nc;
1607     ierr = PetscFree(points);CHKERRQ(ierr);
1608   }
1609   if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(0, name, rCellDof, cCellDof, elemMat);CHKERRQ(ierr);}
1610   for (c = cStart; c < cEnd; ++c) {
1611     ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES);CHKERRQ(ierr);
1612   }
1613   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);}
1614   ierr = PetscFree(feRef);CHKERRQ(ierr);
1615   ierr = PetscFree(elemMat);CHKERRQ(ierr);
1616   ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1617   ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1618   if (mesh->printFEM) {
1619     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
1620     ierr = MatChop(In, 1.0e-10);CHKERRQ(ierr);
1621     ierr = MatView(In, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1622   }
1623 #if 0
1624   ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1625 #endif
1626   PetscFunctionReturn(0);
1627 }
1628 
1629 #undef __FUNCT__
1630 #define __FUNCT__ "DMPlexAddBoundary"
1631 /* The ids can be overridden by the command line option -bc_<boundary name> */
1632 PetscErrorCode DMPlexAddBoundary(DM dm, PetscBool isEssential, const char name[], PetscInt field, void (*bcFunc)(), PetscInt numids, const PetscInt *ids, void *ctx)
1633 {
1634   DM_Plex       *mesh = (DM_Plex *) dm->data;
1635   DMBoundary     b;
1636   PetscErrorCode ierr;
1637 
1638   PetscFunctionBegin;
1639   ierr = PetscNew(&b);CHKERRQ(ierr);
1640   ierr = PetscStrallocpy(name, (char **) &b->name);CHKERRQ(ierr);
1641   ierr = PetscMalloc1(numids, &b->ids);CHKERRQ(ierr);
1642   ierr = PetscMemcpy(b->ids, ids, numids*sizeof(PetscInt));CHKERRQ(ierr);
1643   b->essential   = isEssential;
1644   b->field       = field;
1645   b->func        = bcFunc;
1646   b->numids      = numids;
1647   b->ctx         = ctx;
1648   b->next        = mesh->boundary;
1649   mesh->boundary = b;
1650   PetscFunctionReturn(0);
1651 }
1652 
1653 #undef __FUNCT__
1654 #define __FUNCT__ "DMPlexGetNumBoundary"
1655 PetscErrorCode DMPlexGetNumBoundary(DM dm, PetscInt *numBd)
1656 {
1657   DM_Plex   *mesh = (DM_Plex *) dm->data;
1658   DMBoundary b    = mesh->boundary;
1659 
1660   PetscFunctionBegin;
1661   *numBd = 0;
1662   while (b) {++(*numBd); b = b->next;}
1663   PetscFunctionReturn(0);
1664 }
1665 
1666 #undef __FUNCT__
1667 #define __FUNCT__ "DMPlexGetBoundary"
1668 PetscErrorCode DMPlexGetBoundary(DM dm, PetscInt bd, PetscBool *isEssential, const char **name, PetscInt *field, void (**func)(), PetscInt *numids, const PetscInt **ids, void **ctx)
1669 {
1670   DM_Plex   *mesh = (DM_Plex *) dm->data;
1671   DMBoundary b    = mesh->boundary;
1672   PetscInt   n    = 0;
1673 
1674   PetscFunctionBegin;
1675   while (b) {
1676     if (n == bd) break;
1677     b = b->next;
1678     ++n;
1679   }
1680   if (n != bd) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Boundary %d is not in [0, %d)", bd, n);
1681   if (isEssential) {
1682     PetscValidPointer(isEssential, 3);
1683     *isEssential = b->essential;
1684   }
1685   if (name) {
1686     PetscValidPointer(name, 4);
1687     *name = b->name;
1688   }
1689   if (field) {
1690     PetscValidPointer(field, 5);
1691     *field = b->field;
1692   }
1693   if (func) {
1694     PetscValidPointer(func, 6);
1695     *func = b->func;
1696   }
1697   if (numids) {
1698     PetscValidPointer(numids, 7);
1699     *numids = b->numids;
1700   }
1701   if (ids) {
1702     PetscValidPointer(ids, 8);
1703     *ids = b->ids;
1704   }
1705   if (ctx) {
1706     PetscValidPointer(ctx, 9);
1707     *ctx = b->ctx;
1708   }
1709   PetscFunctionReturn(0);
1710 }
1711