xref: /petsc/src/dm/impls/plex/plexfem.c (revision 7ba4506d7d545948e894dedc80e47cb5d5de7ba1)
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__ "DMPlexProjectFunctionLocal"
185 PetscErrorCode DMPlexProjectFunctionLocal(DM dm, 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, c, f, d, v, comp;
192   PetscErrorCode  ierr;
193 
194   PetscFunctionBegin;
195   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
196   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
197   ierr = PetscMalloc1(numFields, &sp);CHKERRQ(ierr);
198   for (f = 0; f < numFields; ++f) {
199     ierr = PetscFEGetDualSpace(fe[f], &sp[f]);CHKERRQ(ierr);
200     ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr);
201     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
202     totDim += spDim*numComp;
203   }
204   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
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   ierr = PetscMalloc2(dim,&v0,dim*dim,&J);CHKERRQ(ierr);
210   for (c = cStart; c < cEnd; ++c) {
211     PetscCellGeometry geom;
212 
213     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, NULL, &detJ);CHKERRQ(ierr);
214     geom.v0   = v0;
215     geom.J    = J;
216     geom.detJ = &detJ;
217     for (f = 0, v = 0; f < numFields; ++f) {
218       void * const ctx = ctxs ? ctxs[f] : NULL;
219       ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr);
220       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
221       for (d = 0; d < spDim; ++d) {
222         if (funcs[f]) {
223           ierr = PetscDualSpaceApply(sp[f], d, geom, numComp, funcs[f], ctx, &values[v]);CHKERRQ(ierr);
224         } else {
225           for (comp = 0; comp < numComp; ++comp) values[v+comp] = 0.0;
226         }
227         v += numComp;
228       }
229     }
230     ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr);
231   }
232   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
233   ierr = PetscFree2(v0,J);CHKERRQ(ierr);
234   ierr = PetscFree(sp);CHKERRQ(ierr);
235   PetscFunctionReturn(0);
236 }
237 
238 #undef __FUNCT__
239 #define __FUNCT__ "DMPlexProjectFunction"
240 /*@C
241   DMPlexProjectFunction - This projects the given function into the function space provided.
242 
243   Input Parameters:
244 + dm      - The DM
245 . fe      - The PetscFE associated with the field
246 . funcs   - The coordinate functions to evaluate, one per field
247 . ctxs    - Optional array of contexts to pass to each coordinate function.  ctxs itself may be null.
248 - mode    - The insertion mode for values
249 
250   Output Parameter:
251 . X - vector
252 
253   Level: developer
254 
255 .seealso: DMPlexComputeL2Diff()
256 @*/
257 PetscErrorCode DMPlexProjectFunction(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec X)
258 {
259   Vec            localX;
260   PetscErrorCode ierr;
261 
262   PetscFunctionBegin;
263   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
264   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
265   ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, ctxs, mode, localX);CHKERRQ(ierr);
266   ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr);
267   ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr);
268   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
269   PetscFunctionReturn(0);
270 }
271 
272 
273 #undef __FUNCT__
274 #define __FUNCT__ "DMPlexInsertBoundaryValues"
275 PetscErrorCode DMPlexInsertBoundaryValues(DM dm, Vec localX)
276 {
277   void        (**funcs)(const PetscReal x[], PetscScalar *u, void *ctx);
278   void         **ctxs;
279   PetscFE       *fe;
280   PetscInt       numFields, f, numBd, b;
281   PetscErrorCode ierr;
282 
283   PetscFunctionBegin;
284   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
285   PetscValidHeaderSpecific(localX, VEC_CLASSID, 2);
286   ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr);
287   ierr = PetscMalloc3(numFields,&fe,numFields,&funcs,numFields,&ctxs);CHKERRQ(ierr);
288   for (f = 0; f < numFields; ++f) {ierr = DMGetField(dm, f, (PetscObject *) &fe[f]);CHKERRQ(ierr);}
289   /* OPT: Could attempt to do multiple BCs at once */
290   ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
291   for (b = 0; b < numBd; ++b) {
292     const PetscInt *ids;
293     PetscInt        numids, field;
294     PetscBool       isEssential;
295     void          (*func)();
296     void           *ctx;
297 
298     /* TODO: We need to set only the part indicated by the ids */
299     ierr = DMPlexGetBoundary(dm, b, &isEssential, NULL, &field, &func, &numids, &ids, &ctx);CHKERRQ(ierr);
300     for (f = 0; f < numFields; ++f) {
301       funcs[f] = field == f ? (void (*)(const PetscReal[], PetscScalar *, void *)) func : NULL;
302       ctxs[f]  = field == f ? ctx : NULL;
303     }
304     ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
305   }
306   ierr = PetscFree3(fe,funcs,ctxs);CHKERRQ(ierr);
307   PetscFunctionReturn(0);
308 }
309 
310 #undef __FUNCT__
311 #define __FUNCT__ "DMPlexComputeL2Diff"
312 /*@C
313   DMPlexComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.
314 
315   Input Parameters:
316 + dm    - The DM
317 . fe    - The PetscFE object for each field
318 . funcs - The functions to evaluate for each field component
319 . ctxs  - Optional array of contexts to pass to each function, or NULL.
320 - X     - The coefficient vector u_h
321 
322   Output Parameter:
323 . diff - The diff ||u - u_h||_2
324 
325   Level: developer
326 
327 .seealso: DMPlexProjectFunction(), DMPlexComputeL2GradientDiff()
328 @*/
329 PetscErrorCode DMPlexComputeL2Diff(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
330 {
331   const PetscInt  debug = 0;
332   PetscSection    section;
333   PetscQuadrature quad;
334   Vec             localX;
335   PetscScalar    *funcVal;
336   PetscReal      *coords, *v0, *J, *invJ, detJ;
337   PetscReal       localDiff = 0.0;
338   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp;
339   PetscErrorCode  ierr;
340 
341   PetscFunctionBegin;
342   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
343   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
344   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
345   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
346   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
347   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
348   for (field = 0; field < numFields; ++field) {
349     PetscInt Nc;
350 
351     ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr);
352     numComponents += Nc;
353   }
354   ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
355   ierr = PetscMalloc5(numComponents,&funcVal,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
356   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
357   ierr = PetscFEGetQuadrature(fe[0], &quad);CHKERRQ(ierr);
358   for (c = cStart; c < cEnd; ++c) {
359     PetscScalar *x = NULL;
360     PetscReal    elemDiff = 0.0;
361 
362     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
363     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
364     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
365 
366     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
367       void * const     ctx = ctxs ? ctxs[field] : NULL;
368       const PetscReal *quadPoints, *quadWeights;
369       PetscReal       *basis;
370       PetscInt         numQuadPoints, numBasisFuncs, numBasisComps, q, d, e, fc, f;
371 
372       ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
373       ierr = PetscFEGetDimension(fe[field], &numBasisFuncs);CHKERRQ(ierr);
374       ierr = PetscFEGetNumComponents(fe[field], &numBasisComps);CHKERRQ(ierr);
375       ierr = PetscFEGetDefaultTabulation(fe[field], &basis, NULL, NULL);CHKERRQ(ierr);
376       if (debug) {
377         char title[1024];
378         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
379         ierr = DMPrintCellVector(c, title, numBasisFuncs*numBasisComps, &x[fieldOffset]);CHKERRQ(ierr);
380       }
381       for (q = 0; q < numQuadPoints; ++q) {
382         for (d = 0; d < dim; d++) {
383           coords[d] = v0[d];
384           for (e = 0; e < dim; e++) {
385             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
386           }
387         }
388         (*funcs[field])(coords, funcVal, ctx);
389         for (fc = 0; fc < numBasisComps; ++fc) {
390           PetscScalar interpolant = 0.0;
391 
392           for (f = 0; f < numBasisFuncs; ++f) {
393             const PetscInt fidx = f*numBasisComps+fc;
394             interpolant += x[fieldOffset+fidx]*basis[q*numBasisFuncs*numBasisComps+fidx];
395           }
396           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);}
397           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ;
398         }
399       }
400       comp        += numBasisComps;
401       fieldOffset += numBasisFuncs*numBasisComps;
402     }
403     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
404     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
405     localDiff += elemDiff;
406   }
407   ierr  = PetscFree5(funcVal,coords,v0,J,invJ);CHKERRQ(ierr);
408   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
409   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
410   *diff = PetscSqrtReal(*diff);
411   PetscFunctionReturn(0);
412 }
413 
414 #undef __FUNCT__
415 #define __FUNCT__ "DMPlexComputeL2GradientDiff"
416 /*@C
417   DMPlexComputeL2GradientDiff - This function computes the L_2 difference between the gradient of a function u and an FEM interpolant solution grad u_h.
418 
419   Input Parameters:
420 + dm    - The DM
421 . fe    - The PetscFE object for each field
422 . funcs - The gradient functions to evaluate for each field component
423 . ctxs  - Optional array of contexts to pass to each function, or NULL.
424 . X     - The coefficient vector u_h
425 - n     - The vector to project along
426 
427   Output Parameter:
428 . diff - The diff ||(grad u - grad u_h) . n||_2
429 
430   Level: developer
431 
432 .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff()
433 @*/
434 PetscErrorCode DMPlexComputeL2GradientDiff(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, const PetscReal n[], PetscReal *diff)
435 {
436   const PetscInt  debug = 0;
437   PetscSection    section;
438   PetscQuadrature quad;
439   Vec             localX;
440   PetscScalar    *funcVal, *interpolantVec;
441   PetscReal      *coords, *realSpaceDer, *v0, *J, *invJ, detJ;
442   PetscReal       localDiff = 0.0;
443   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp;
444   PetscErrorCode  ierr;
445 
446   PetscFunctionBegin;
447   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
448   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
449   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
450   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
451   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
452   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
453   for (field = 0; field < numFields; ++field) {
454     PetscInt Nc;
455 
456     ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr);
457     numComponents += Nc;
458   }
459   /* ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */
460   ierr = PetscMalloc7(numComponents,&funcVal,dim,&coords,dim,&realSpaceDer,dim,&v0,dim*dim,&J,dim*dim,&invJ,dim,&interpolantVec);CHKERRQ(ierr);
461   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
462   ierr = PetscFEGetQuadrature(fe[0], &quad);CHKERRQ(ierr);
463   for (c = cStart; c < cEnd; ++c) {
464     PetscScalar *x = NULL;
465     PetscReal    elemDiff = 0.0;
466 
467     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
468     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
469     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
470 
471     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
472       void * const     ctx = ctxs ? ctxs[field] : NULL;
473       const PetscReal *quadPoints, *quadWeights;
474       PetscReal       *basisDer;
475       PetscInt         numQuadPoints, Nb, Ncomp, q, d, e, fc, f, g;
476 
477       ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
478       ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr);
479       ierr = PetscFEGetNumComponents(fe[field], &Ncomp);CHKERRQ(ierr);
480       ierr = PetscFEGetDefaultTabulation(fe[field], NULL, &basisDer, NULL);CHKERRQ(ierr);
481       if (debug) {
482         char title[1024];
483         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
484         ierr = DMPrintCellVector(c, title, Nb*Ncomp, &x[fieldOffset]);CHKERRQ(ierr);
485       }
486       for (q = 0; q < numQuadPoints; ++q) {
487         for (d = 0; d < dim; d++) {
488           coords[d] = v0[d];
489           for (e = 0; e < dim; e++) {
490             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
491           }
492         }
493         (*funcs[field])(coords, n, funcVal, ctx);
494         for (fc = 0; fc < Ncomp; ++fc) {
495           PetscScalar interpolant = 0.0;
496 
497           for (d = 0; d < dim; ++d) interpolantVec[d] = 0.0;
498           for (f = 0; f < Nb; ++f) {
499             const PetscInt fidx = f*Ncomp+fc;
500 
501             for (d = 0; d < dim; ++d) {
502               realSpaceDer[d] = 0.0;
503               for (g = 0; g < dim; ++g) {
504                 realSpaceDer[d] += invJ[g*dim+d]*basisDer[(q*Nb*Ncomp+fidx)*dim+g];
505               }
506               interpolantVec[d] += x[fieldOffset+fidx]*realSpaceDer[d];
507             }
508           }
509           for (d = 0; d < dim; ++d) interpolant += interpolantVec[d]*n[d];
510           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);}
511           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ;
512         }
513       }
514       comp        += Ncomp;
515       fieldOffset += Nb*Ncomp;
516     }
517     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
518     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
519     localDiff += elemDiff;
520   }
521   ierr  = PetscFree7(funcVal,coords,realSpaceDer,v0,J,invJ,interpolantVec);CHKERRQ(ierr);
522   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
523   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
524   *diff = PetscSqrtReal(*diff);
525   PetscFunctionReturn(0);
526 }
527 
528 #undef __FUNCT__
529 #define __FUNCT__ "DMPlexComputeResidualFEM"
530 /*@
531   DMPlexComputeResidualFEM - Form the local residual F from the local input X using pointwise functions specified by the user
532 
533   Input Parameters:
534 + dm - The mesh
535 . X  - Local input vector
536 - user - The user context
537 
538   Output Parameter:
539 . F  - Local output vector
540 
541   Note:
542   The first member of the user context must be an FEMContext.
543 
544   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
545   like a GPU, or vectorize on a multicore machine.
546 
547   Level: developer
548 
549 .seealso: DMPlexComputeJacobianActionFEM()
550 @*/
551 PetscErrorCode DMPlexComputeResidualFEM(DM dm, Vec X, Vec F, void *user)
552 {
553   DM_Plex          *mesh  = (DM_Plex *) dm->data;
554   PetscFEM         *fem   = (PetscFEM *) user;
555   PetscFE          *fe    = fem->fe;
556   PetscFE          *feAux = fem->feAux;
557   PetscFE          *feBd  = fem->feBd;
558   const char       *name  = "Residual";
559   DM                dmAux;
560   Vec               A;
561   PetscQuadrature   q;
562   PetscCellGeometry geom;
563   PetscSection      section, sectionAux;
564   PetscReal        *v0, *J, *invJ, *detJ;
565   PetscScalar      *elemVec, *u, *a = NULL;
566   PetscInt          dim, Nf, NfAux = 0, f, numCells, cStart, cEnd, c;
567   PetscInt          cellDof = 0, numComponents = 0;
568   PetscInt          cellDofAux = 0, numComponentsAux = 0;
569   PetscErrorCode    ierr;
570 
571   PetscFunctionBegin;
572   ierr = PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
573   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
574   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
575   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
576   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
577   numCells = cEnd - cStart;
578   for (f = 0; f < Nf; ++f) {
579     PetscInt Nb, Nc;
580 
581     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
582     ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr);
583     cellDof       += Nb*Nc;
584     numComponents += Nc;
585   }
586   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
587   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
588   if (dmAux) {
589     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
590     ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr);
591   }
592   for (f = 0; f < NfAux; ++f) {
593     PetscInt Nb, Nc;
594 
595     ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr);
596     ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr);
597     cellDofAux       += Nb*Nc;
598     numComponentsAux += Nc;
599   }
600   ierr = DMPlexProjectFunctionLocal(dm, fe, fem->bcFuncs, fem->bcCtxs, INSERT_BC_VALUES, X);CHKERRQ(ierr);
601   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
602   ierr = PetscMalloc6(numCells*cellDof,&u,numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*cellDof,&elemVec);CHKERRQ(ierr);
603   if (dmAux) {ierr = PetscMalloc1(numCells*cellDofAux, &a);CHKERRQ(ierr);}
604   for (c = cStart; c < cEnd; ++c) {
605     PetscScalar *x = NULL;
606     PetscInt     i;
607 
608     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
609     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
610     ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
611     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
612     ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
613     if (dmAux) {
614       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
615       for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x[i];
616       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
617     }
618   }
619   for (f = 0; f < Nf; ++f) {
620     void   (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f0Funcs[f];
621     void   (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f1Funcs[f];
622     PetscInt numQuadPoints, Nb;
623     /* Conforming batches */
624     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
625     /* Remainder */
626     PetscInt Nr, offset;
627 
628     ierr = PetscFEGetQuadrature(fe[f], &q);CHKERRQ(ierr);
629     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
630     ierr = PetscFEGetTileSizes(fe[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
631     ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
632     blockSize = Nb*numQuadPoints;
633     batchSize = numBlocks * blockSize;
634     ierr =  PetscFESetTileSizes(fe[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
635     numChunks = numCells / (numBatches*batchSize);
636     Ne        = numChunks*numBatches*batchSize;
637     Nr        = numCells % (numBatches*batchSize);
638     offset    = numCells - Nr;
639     geom.v0   = v0;
640     geom.J    = J;
641     geom.invJ = invJ;
642     geom.detJ = detJ;
643     ierr = PetscFEIntegrateResidual(fe[f], Ne, Nf, fe, f, geom, u, NfAux, feAux, a, f0, f1, elemVec);CHKERRQ(ierr);
644     geom.v0   = &v0[offset*dim];
645     geom.J    = &J[offset*dim*dim];
646     geom.invJ = &invJ[offset*dim*dim];
647     geom.detJ = &detJ[offset];
648     ierr = PetscFEIntegrateResidual(fe[f], Nr, Nf, fe, f, geom, &u[offset*cellDof], NfAux, feAux, &a[offset*cellDofAux], f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr);
649   }
650   for (c = cStart; c < cEnd; ++c) {
651     if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, name, cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);}
652     ierr = DMPlexVecSetClosure(dm, section, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr);
653   }
654   ierr = PetscFree6(u,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr);
655   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
656   if (feBd) {
657     DMLabel         label, depth;
658     IS              pointIS;
659     const PetscInt *points;
660     PetscInt        dep, numPoints, p, numFaces;
661     PetscReal      *n;
662 
663     ierr = DMPlexGetLabel(dm, "boundary", &label);CHKERRQ(ierr);
664     ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr);
665     ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr);
666     ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr);
667     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
668     for (f = 0, cellDof = 0, numComponents = 0; f < Nf; ++f) {
669       PetscInt Nb, Nc;
670 
671       ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr);
672       ierr = PetscFEGetNumComponents(feBd[f], &Nc);CHKERRQ(ierr);
673       cellDof       += Nb*Nc;
674       numComponents += Nc;
675     }
676     for (p = 0, numFaces = 0; p < numPoints; ++p) {
677       ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
678       if (dep == dim-1) ++numFaces;
679     }
680     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);
681     for (p = 0, f = 0; p < numPoints; ++p) {
682       const PetscInt point = points[p];
683       PetscScalar   *x     = NULL;
684       PetscInt       i;
685 
686       ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr);
687       if (dep != dim-1) continue;
688       ierr = DMPlexComputeCellGeometry(dm, point, &v0[f*dim], &J[f*dim*dim], &invJ[f*dim*dim], &detJ[f]);CHKERRQ(ierr);
689       ierr = DMPlexComputeCellGeometryFVM(dm, point, NULL, NULL, &n[f*dim]);
690       if (detJ[f] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[f], point);
691       ierr = DMPlexVecGetClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
692       for (i = 0; i < cellDof; ++i) u[f*cellDof+i] = x[i];
693       ierr = DMPlexVecRestoreClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
694       ++f;
695     }
696     for (f = 0; f < Nf; ++f) {
697       void   (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f0BdFuncs[f];
698       void   (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f1BdFuncs[f];
699       PetscInt numQuadPoints, Nb;
700       /* Conforming batches */
701       PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
702       /* Remainder */
703       PetscInt Nr, offset;
704 
705       ierr = PetscFEGetQuadrature(feBd[f], &q);CHKERRQ(ierr);
706       ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr);
707       ierr = PetscFEGetTileSizes(feBd[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
708       ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
709       blockSize = Nb*numQuadPoints;
710       batchSize = numBlocks * blockSize;
711       ierr =  PetscFESetTileSizes(feBd[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
712       numChunks = numFaces / (numBatches*batchSize);
713       Ne        = numChunks*numBatches*batchSize;
714       Nr        = numFaces % (numBatches*batchSize);
715       offset    = numFaces - Nr;
716       geom.v0   = v0;
717       geom.n    = n;
718       geom.J    = J;
719       geom.invJ = invJ;
720       geom.detJ = detJ;
721       ierr = PetscFEIntegrateBdResidual(feBd[f], Ne, Nf, feBd, f, geom, u, 0, NULL, NULL, f0, f1, elemVec);CHKERRQ(ierr);
722       geom.v0   = &v0[offset*dim];
723       geom.n    = &n[offset*dim];
724       geom.J    = &J[offset*dim*dim];
725       geom.invJ = &invJ[offset*dim*dim];
726       geom.detJ = &detJ[offset];
727       ierr = PetscFEIntegrateBdResidual(feBd[f], Nr, Nf, feBd, f, geom, &u[offset*cellDof], 0, NULL, NULL, f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr);
728     }
729     for (p = 0, f = 0; p < numPoints; ++p) {
730       const PetscInt point = points[p];
731 
732       ierr = DMLabelGetValue(depth, point, &dep);CHKERRQ(ierr);
733       if (dep != dim-1) continue;
734       if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "BdResidual", cellDof, &elemVec[f*cellDof]);CHKERRQ(ierr);}
735       ierr = DMPlexVecSetClosure(dm, NULL, F, point, &elemVec[f*cellDof], ADD_VALUES);CHKERRQ(ierr);
736       ++f;
737     }
738     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
739     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
740     ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemVec);CHKERRQ(ierr);
741   }
742   if (mesh->printFEM) {ierr = DMPrintLocalVec(dm, name, mesh->printTol, F);CHKERRQ(ierr);}
743   ierr = PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
744   PetscFunctionReturn(0);
745 }
746 
747 #undef __FUNCT__
748 #define __FUNCT__ "DMPlexComputeJacobianActionFEM"
749 /*@C
750   DMPlexComputeJacobianActionFEM - Form the local action of Jacobian J(u) on the local input X using pointwise functions specified by the user
751 
752   Input Parameters:
753 + dm - The mesh
754 . J  - The Jacobian shell matrix
755 . X  - Local input vector
756 - user - The user context
757 
758   Output Parameter:
759 . F  - Local output vector
760 
761   Note:
762   The first member of the user context must be an FEMContext.
763 
764   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
765   like a GPU, or vectorize on a multicore machine.
766 
767   Level: developer
768 
769 .seealso: DMPlexComputeResidualFEM()
770 @*/
771 PetscErrorCode DMPlexComputeJacobianActionFEM(DM dm, Mat Jac, Vec X, Vec F, void *user)
772 {
773   DM_Plex          *mesh = (DM_Plex *) dm->data;
774   PetscFEM         *fem  = (PetscFEM *) user;
775   PetscFE          *fe   = fem->fe;
776   PetscQuadrature   quad;
777   PetscCellGeometry geom;
778   PetscSection      section;
779   JacActionCtx     *jctx;
780   PetscReal        *v0, *J, *invJ, *detJ;
781   PetscScalar      *elemVec, *u, *a;
782   PetscInt          dim, numFields, field, numCells, cStart, cEnd, c;
783   PetscInt          cellDof = 0;
784   PetscErrorCode    ierr;
785 
786   PetscFunctionBegin;
787   /* ierr = PetscLogEventBegin(DMPLEX_JacobianActionFEM,dm,0,0,0);CHKERRQ(ierr); */
788   ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr);
789   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
790   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
791   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
792   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
793   numCells = cEnd - cStart;
794   for (field = 0; field < numFields; ++field) {
795     PetscInt Nb, Nc;
796 
797     ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr);
798     ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr);
799     cellDof += Nb*Nc;
800   }
801   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
802   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);
803   for (c = cStart; c < cEnd; ++c) {
804     PetscScalar *x = NULL;
805     PetscInt     i;
806 
807     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
808     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
809     ierr = DMPlexVecGetClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr);
810     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
811     ierr = DMPlexVecRestoreClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr);
812     ierr = DMPlexVecGetClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr);
813     for (i = 0; i < cellDof; ++i) a[c*cellDof+i] = x[i];
814     ierr = DMPlexVecRestoreClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr);
815   }
816   for (field = 0; field < numFields; ++field) {
817     PetscInt numQuadPoints, Nb;
818     /* Conforming batches */
819     PetscInt numBlocks  = 1;
820     PetscInt numBatches = 1;
821     PetscInt numChunks, Ne, blockSize, batchSize;
822     /* Remainder */
823     PetscInt Nr, offset;
824 
825     ierr = PetscFEGetQuadrature(fe[field], &quad);CHKERRQ(ierr);
826     ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr);
827     ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
828     blockSize = Nb*numQuadPoints;
829     batchSize = numBlocks * blockSize;
830     numChunks = numCells / (numBatches*batchSize);
831     Ne        = numChunks*numBatches*batchSize;
832     Nr        = numCells % (numBatches*batchSize);
833     offset    = numCells - Nr;
834     geom.v0   = v0;
835     geom.J    = J;
836     geom.invJ = invJ;
837     geom.detJ = detJ;
838     ierr = PetscFEIntegrateJacobianAction(fe[field], Ne, numFields, fe, field, geom, u, a, fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, elemVec);CHKERRQ(ierr);
839     geom.v0   = &v0[offset*dim];
840     geom.J    = &J[offset*dim*dim];
841     geom.invJ = &invJ[offset*dim*dim];
842     geom.detJ = &detJ[offset];
843     ierr = PetscFEIntegrateJacobianAction(fe[field], Nr, numFields, fe, field, geom, &u[offset*cellDof], &a[offset*cellDof],
844                                           fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, &elemVec[offset*cellDof]);CHKERRQ(ierr);
845   }
846   for (c = cStart; c < cEnd; ++c) {
847     if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, "Jacobian Action", cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);}
848     ierr = DMPlexVecSetClosure(dm, NULL, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr);
849   }
850   ierr = PetscFree7(u,a,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr);
851   if (mesh->printFEM) {
852     PetscMPIInt rank, numProcs;
853     PetscInt    p;
854 
855     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
856     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm), &numProcs);CHKERRQ(ierr);
857     ierr = PetscPrintf(PetscObjectComm((PetscObject)dm), "Jacobian Action:\n");CHKERRQ(ierr);
858     for (p = 0; p < numProcs; ++p) {
859       if (p == rank) {ierr = VecView(F, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);}
860       ierr = PetscBarrier((PetscObject) dm);CHKERRQ(ierr);
861     }
862   }
863   /* ierr = PetscLogEventEnd(DMPLEX_JacobianActionFEM,dm,0,0,0);CHKERRQ(ierr); */
864   PetscFunctionReturn(0);
865 }
866 
867 #undef __FUNCT__
868 #define __FUNCT__ "DMPlexComputeJacobianFEM"
869 /*@
870   DMPlexComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user.
871 
872   Input Parameters:
873 + dm - The mesh
874 . X  - Local input vector
875 - user - The user context
876 
877   Output Parameter:
878 . Jac  - Jacobian matrix
879 
880   Note:
881   The first member of the user context must be an FEMContext.
882 
883   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
884   like a GPU, or vectorize on a multicore machine.
885 
886   Level: developer
887 
888 .seealso: FormFunctionLocal()
889 @*/
890 PetscErrorCode DMPlexComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP,void *user)
891 {
892   DM_Plex          *mesh  = (DM_Plex *) dm->data;
893   PetscFEM         *fem   = (PetscFEM *) user;
894   PetscFE          *fe    = fem->fe;
895   PetscFE          *feAux = fem->feAux;
896   const char       *name  = "Jacobian";
897   DM                dmAux;
898   Vec               A;
899   PetscQuadrature   quad;
900   PetscCellGeometry geom;
901   PetscSection      section, globalSection, sectionAux;
902   PetscReal        *v0, *J, *invJ, *detJ;
903   PetscScalar      *elemMat, *u, *a;
904   PetscInt          dim, Nf, NfAux = 0, f, fieldI, fieldJ, numCells, cStart, cEnd, c;
905   PetscInt          cellDof = 0, numComponents = 0;
906   PetscInt          cellDofAux = 0, numComponentsAux = 0;
907   PetscBool         isShell;
908   PetscErrorCode    ierr;
909 
910   PetscFunctionBegin;
911   ierr = PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
912   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
913   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
914   ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);
915   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
916   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
917   numCells = cEnd - cStart;
918   for (f = 0; f < Nf; ++f) {
919     PetscInt Nb, Nc;
920 
921     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
922     ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr);
923     cellDof       += Nb*Nc;
924     numComponents += Nc;
925   }
926   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
927   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
928   if (dmAux) {
929     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
930     ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr);
931   }
932   for (f = 0; f < NfAux; ++f) {
933     PetscInt Nb, Nc;
934 
935     ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr);
936     ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr);
937     cellDofAux       += Nb*Nc;
938     numComponentsAux += Nc;
939   }
940   ierr = DMPlexProjectFunctionLocal(dm, fe, fem->bcFuncs, fem->bcCtxs, INSERT_BC_VALUES, X);CHKERRQ(ierr);
941   ierr = MatZeroEntries(JacP);CHKERRQ(ierr);
942   ierr = PetscMalloc6(numCells*cellDof,&u,numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*cellDof*cellDof,&elemMat);CHKERRQ(ierr);
943   if (dmAux) {ierr = PetscMalloc1(numCells*cellDofAux, &a);CHKERRQ(ierr);}
944   for (c = cStart; c < cEnd; ++c) {
945     PetscScalar *x = NULL;
946     PetscInt     i;
947 
948     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
949     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
950     ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
951     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
952     ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
953     if (dmAux) {
954       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
955       for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x[i];
956       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
957     }
958   }
959   ierr = PetscMemzero(elemMat, numCells*cellDof*cellDof * sizeof(PetscScalar));CHKERRQ(ierr);
960   for (fieldI = 0; fieldI < Nf; ++fieldI) {
961     PetscInt numQuadPoints, Nb;
962     /* Conforming batches */
963     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
964     /* Remainder */
965     PetscInt Nr, offset;
966 
967     ierr = PetscFEGetQuadrature(fe[fieldI], &quad);CHKERRQ(ierr);
968     ierr = PetscFEGetDimension(fe[fieldI], &Nb);CHKERRQ(ierr);
969     ierr = PetscFEGetTileSizes(fe[fieldI], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
970     ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
971     blockSize = Nb*numQuadPoints;
972     batchSize = numBlocks * blockSize;
973     ierr = PetscFESetTileSizes(fe[fieldI], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
974     numChunks = numCells / (numBatches*batchSize);
975     Ne        = numChunks*numBatches*batchSize;
976     Nr        = numCells % (numBatches*batchSize);
977     offset    = numCells - Nr;
978     for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
979       void   (*g0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g0Funcs[fieldI*Nf+fieldJ];
980       void   (*g1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g1Funcs[fieldI*Nf+fieldJ];
981       void   (*g2)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g2Funcs[fieldI*Nf+fieldJ];
982       void   (*g3)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g3Funcs[fieldI*Nf+fieldJ];
983 
984       geom.v0   = v0;
985       geom.J    = J;
986       geom.invJ = invJ;
987       geom.detJ = detJ;
988       ierr = PetscFEIntegrateJacobian(fe[fieldI], Ne, Nf, fe, fieldI, fieldJ, geom, u, NfAux, feAux, a, g0, g1, g2, g3, elemMat);CHKERRQ(ierr);
989       geom.v0   = &v0[offset*dim];
990       geom.J    = &J[offset*dim*dim];
991       geom.invJ = &invJ[offset*dim*dim];
992       geom.detJ = &detJ[offset];
993       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);
994     }
995   }
996   for (c = cStart; c < cEnd; ++c) {
997     if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, name, cellDof, cellDof, &elemMat[c*cellDof*cellDof]);CHKERRQ(ierr);}
998     ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, c, &elemMat[c*cellDof*cellDof], ADD_VALUES);CHKERRQ(ierr);
999   }
1000   ierr = PetscFree6(u,v0,J,invJ,detJ,elemMat);CHKERRQ(ierr);
1001   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
1002   ierr = MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1003   ierr = MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1004   if (mesh->printFEM) {
1005     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
1006     ierr = MatChop(JacP, 1.0e-10);CHKERRQ(ierr);
1007     ierr = MatView(JacP, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1008   }
1009   ierr = PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
1010   ierr = PetscObjectTypeCompare((PetscObject) Jac, MATSHELL, &isShell);CHKERRQ(ierr);
1011   if (isShell) {
1012     JacActionCtx *jctx;
1013 
1014     ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr);
1015     ierr = VecCopy(X, jctx->u);CHKERRQ(ierr);
1016   }
1017   PetscFunctionReturn(0);
1018 }
1019 
1020 #if 0
1021 
1022 static void g0_identity_1d_static(const PetscScalar u[], const PetscScalar gradU[], const PetscScalar a[], const PetscScalar gradA[], const PetscReal x[], PetscScalar g3[])
1023 {
1024   g3[0] = 1.0;
1025 }
1026 
1027 static void g0_identity_2d_static(const PetscScalar u[], const PetscScalar gradU[], const PetscScalar a[], const PetscScalar gradA[], const PetscReal x[], PetscScalar g3[])
1028 {
1029   g3[0] = g3[3] = 1.0;
1030 }
1031 
1032 static void g0_identity_3d_static(const PetscScalar u[], const PetscScalar gradU[], const PetscScalar a[], const PetscScalar gradA[], const PetscReal x[], PetscScalar g3[])
1033 {
1034   g3[0] = g3[4] = g3[8] = 1.0;
1035 }
1036 
1037 #undef __FUNCT__
1038 #define __FUNCT__ "DMPlexComputeInterpolatorFEMBroken"
1039 PetscErrorCode DMPlexComputeInterpolatorFEMBroken(DM dmc, DM dmf, Mat I, void *user)
1040 {
1041   DM_Plex          *mesh  = (DM_Plex *) dmc->data;
1042   PetscFEM         *fem   = (PetscFEM *) user;
1043   PetscFE          *fe    = fem->fe;
1044   const char       *name  = "Interpolator";
1045   PetscFE          *feRef;
1046   PetscQuadrature   quad, quadOld;
1047   PetscCellGeometry geom;
1048   PetscSection      fsection, fglobalSection;
1049   PetscSection      csection, cglobalSection;
1050   PetscReal        *v0, *J, *invJ, *detJ;
1051   PetscScalar      *elemMat;
1052   PetscInt          dim, Nf, f, fieldI, fieldJ, numCells, cStart, cEnd, c;
1053   PetscInt          rCellDof = 0, cCellDof = 0, numComponents = 0;
1054   PetscErrorCode    ierr;
1055 
1056   PetscFunctionBegin;
1057 #if 0
1058   ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1059 #endif
1060   ierr = DMPlexGetDimension(dmf, &dim);CHKERRQ(ierr);
1061   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
1062   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
1063   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
1064   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
1065   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
1066   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
1067   numCells = cEnd - cStart;
1068   ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr);
1069   for (fieldI = 0; fieldI < Nf; ++fieldI) {
1070     ierr = PetscFERefine(fe[fieldI], &feRef[fieldI]);CHKERRQ(ierr);
1071   }
1072   for (f = 0; f < Nf; ++f) {
1073     PetscInt rNb, cNb, Nc;
1074 
1075     ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr);
1076     ierr = PetscFEGetDimension(fe[f], &cNb);CHKERRQ(ierr);
1077     ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr);
1078     numComponents += Nc;
1079     rCellDof += rNb*Nc;
1080     cCellDof += cNb*Nc;
1081   }
1082   ierr = MatZeroEntries(I);CHKERRQ(ierr);
1083   ierr = PetscMalloc5(numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*rCellDof*cCellDof,&elemMat);CHKERRQ(ierr);
1084   for (c = cStart; c < cEnd; ++c) {
1085     ierr = DMPlexComputeCellGeometry(dmc, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
1086     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
1087   }
1088   ierr = PetscMemzero(elemMat, numCells*rCellDof*cCellDof * sizeof(PetscScalar));CHKERRQ(ierr);
1089   for (fieldI = 0; fieldI < Nf; ++fieldI) {
1090     PetscInt numQuadPoints, Nb;
1091     /* Conforming batches */
1092     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1093     /* Remainder */
1094     PetscInt Nr, offset;
1095 
1096     /* Make new fine FE which refines the ref cell and the quadrature rule */
1097     ierr = PetscFEGetQuadrature(feRef[fieldI], &quad);CHKERRQ(ierr);
1098     ierr = PetscFEGetDimension(feRef[fieldI], &Nb);CHKERRQ(ierr);
1099     ierr = PetscFEGetTileSizes(feRef[fieldI], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
1100     ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
1101     blockSize = Nb*numQuadPoints;
1102     batchSize = numBlocks * blockSize;
1103     ierr = PetscFESetTileSizes(feRef[fieldI], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1104     numChunks = numCells / (numBatches*batchSize);
1105     Ne        = numChunks*numBatches*batchSize;
1106     Nr        = numCells % (numBatches*batchSize);
1107     offset    = numCells - Nr;
1108 
1109     for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
1110       /* void   (*g0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g0Funcs[fieldI*Nf+fieldJ]; */
1111       void   (*g0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = g0_identity_2d_static;
1112 
1113       /* Replace quadrature in coarse FE with refined quadrature */
1114       ierr = PetscFEGetQuadrature(fe[fieldJ], &quadOld);CHKERRQ(ierr);
1115       ierr = PetscObjectReference((PetscObject) quadOld);CHKERRQ(ierr);
1116       ierr = PetscFESetQuadrature(fe[fieldJ], quad);CHKERRQ(ierr);
1117       geom.v0   = v0;
1118       geom.J    = J;
1119       geom.invJ = invJ;
1120       geom.detJ = detJ;
1121       ierr = PetscFEIntegrateInterpolator_Basic(feRef[fieldI], Ne, Nf, feRef, fieldI, fe, fieldJ, geom, g0, elemMat);CHKERRQ(ierr);
1122       geom.v0   = &v0[offset*dim];
1123       geom.J    = &J[offset*dim*dim];
1124       geom.invJ = &invJ[offset*dim*dim];
1125       geom.detJ = &detJ[offset];
1126       ierr = PetscFEIntegrateInterpolator_Basic(feRef[fieldI], Nr, Nf, feRef, fieldI, fe, fieldJ, geom, g0, &elemMat[offset*rCellDof*cCellDof]);CHKERRQ(ierr);
1127       ierr = PetscFESetQuadrature(fe[fieldJ], quadOld);CHKERRQ(ierr);
1128     }
1129   }
1130   for (c = cStart; c < cEnd; ++c) {
1131     if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, name, rCellDof, cCellDof, &elemMat[c*rCellDof*cCellDof]);CHKERRQ(ierr);}
1132     ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, I, c, &elemMat[c*rCellDof*cCellDof], ADD_VALUES);CHKERRQ(ierr);
1133   }
1134   ierr = PetscFree5(v0,J,invJ,detJ,elemMat);CHKERRQ(ierr);
1135   ierr = PetscFree(feRef);CHKERRQ(ierr);
1136   ierr = MatAssemblyBegin(I, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1137   ierr = MatAssemblyEnd(I, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1138   if (mesh->printFEM) {
1139     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
1140     ierr = MatChop(I, 1.0e-10);CHKERRQ(ierr);
1141     ierr = MatView(I, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1142   }
1143 #if 0
1144   ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1145 #endif
1146   PetscFunctionReturn(0);
1147 }
1148 #endif
1149 
1150 #undef __FUNCT__
1151 #define __FUNCT__ "DMPlexComputeInterpolatorFEM"
1152 /*@
1153   DMPlexComputeInterpolatorFEM - Form the local portion of the interpolation matrix I from the coarse DM to the uniformly refined DM.
1154 
1155   Input Parameters:
1156 + dmf  - The fine mesh
1157 . dmc  - The coarse mesh
1158 - user - The user context
1159 
1160   Output Parameter:
1161 . In  - The interpolation matrix
1162 
1163   Note:
1164   The first member of the user context must be an FEMContext.
1165 
1166   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
1167   like a GPU, or vectorize on a multicore machine.
1168 
1169   Level: developer
1170 
1171 .seealso: DMPlexComputeJacobianFEM()
1172 @*/
1173 PetscErrorCode DMPlexComputeInterpolatorFEM(DM dmc, DM dmf, Mat In, void *user)
1174 {
1175   DM_Plex          *mesh  = (DM_Plex *) dmc->data;
1176   PetscFEM         *fem   = (PetscFEM *) user;
1177   PetscFE          *fe    = fem->fe;
1178   const char       *name  = "Interpolator";
1179   PetscFE          *feRef;
1180   PetscSection      fsection, fglobalSection;
1181   PetscSection      csection, cglobalSection;
1182   PetscScalar      *elemMat;
1183   PetscInt          dim, Nf, f, fieldI, fieldJ, offsetI, offsetJ, cStart, cEnd, c;
1184   PetscInt          rCellDof = 0, cCellDof = 0;
1185   PetscErrorCode    ierr;
1186 
1187   PetscFunctionBegin;
1188 #if 0
1189   ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1190 #endif
1191   ierr = DMPlexGetDimension(dmf, &dim);CHKERRQ(ierr);
1192   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
1193   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
1194   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
1195   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
1196   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
1197   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
1198   ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr);
1199   for (f = 0; f < Nf; ++f) {
1200     PetscInt rNb, cNb, Nc;
1201 
1202     ierr = PetscFERefine(fe[f], &feRef[f]);CHKERRQ(ierr);
1203     ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr);
1204     ierr = PetscFEGetDimension(fe[f], &cNb);CHKERRQ(ierr);
1205     ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr);
1206     rCellDof += rNb*Nc;
1207     cCellDof += cNb*Nc;
1208   }
1209   ierr = MatZeroEntries(In);CHKERRQ(ierr);
1210   ierr = PetscMalloc1(rCellDof*cCellDof,&elemMat);CHKERRQ(ierr);
1211   ierr = PetscMemzero(elemMat, rCellDof*cCellDof * sizeof(PetscScalar));CHKERRQ(ierr);
1212   for (fieldI = 0, offsetI = 0; fieldI < Nf; ++fieldI) {
1213     PetscDualSpace   Qref;
1214     PetscQuadrature  f;
1215     const PetscReal *qpoints, *qweights;
1216     PetscReal       *points;
1217     PetscInt         npoints = 0, Nc, Np, fpdim, i, k, p, d;
1218 
1219     /* Compose points from all dual basis functionals */
1220     ierr = PetscFEGetNumComponents(fe[fieldI], &Nc);CHKERRQ(ierr);
1221     ierr = PetscFEGetDualSpace(feRef[fieldI], &Qref);CHKERRQ(ierr);
1222     ierr = PetscDualSpaceGetDimension(Qref, &fpdim);CHKERRQ(ierr);
1223     for (i = 0; i < fpdim; ++i) {
1224       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1225       ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, NULL);CHKERRQ(ierr);
1226       npoints += Np;
1227     }
1228     ierr = PetscMalloc1(npoints*dim,&points);CHKERRQ(ierr);
1229     for (i = 0, k = 0; i < fpdim; ++i) {
1230       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1231       ierr = PetscQuadratureGetData(f, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr);
1232       for (p = 0; p < Np; ++p, ++k) for (d = 0; d < dim; ++d) points[k*dim+d] = qpoints[p*dim+d];
1233     }
1234 
1235     for (fieldJ = 0, offsetJ = 0; fieldJ < Nf; ++fieldJ) {
1236       PetscReal *B;
1237       PetscInt   NcJ, cpdim, j;
1238 
1239       /* Evaluate basis at points */
1240       ierr = PetscFEGetNumComponents(fe[fieldJ], &NcJ);CHKERRQ(ierr);
1241       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);
1242       ierr = PetscFEGetDimension(fe[fieldJ], &cpdim);CHKERRQ(ierr);
1243       /* For now, fields only interpolate themselves */
1244       if (fieldI == fieldJ) {
1245         ierr = PetscFEGetTabulation(fe[fieldJ], npoints, points, &B, NULL, NULL);CHKERRQ(ierr);
1246         for (i = 0, k = 0; i < fpdim; ++i) {
1247           ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1248           ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, &qweights);CHKERRQ(ierr);
1249           for (p = 0; p < Np; ++p, ++k) {
1250             for (j = 0; j < cpdim; ++j) {
1251               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];
1252             }
1253           }
1254         }
1255         ierr = PetscFERestoreTabulation(fe[fieldJ], npoints, points, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr);
1256       }
1257       offsetJ += cpdim*NcJ;
1258     }
1259     offsetI += fpdim*Nc;
1260     ierr = PetscFree(points);CHKERRQ(ierr);
1261   }
1262   if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(0, name, rCellDof, cCellDof, elemMat);CHKERRQ(ierr);}
1263   for (c = cStart; c < cEnd; ++c) {
1264     ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES);CHKERRQ(ierr);
1265   }
1266   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);}
1267   ierr = PetscFree(feRef);CHKERRQ(ierr);
1268   ierr = PetscFree(elemMat);CHKERRQ(ierr);
1269   ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1270   ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1271   if (mesh->printFEM) {
1272     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
1273     ierr = MatChop(In, 1.0e-10);CHKERRQ(ierr);
1274     ierr = MatView(In, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1275   }
1276 #if 0
1277   ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1278 #endif
1279   PetscFunctionReturn(0);
1280 }
1281 
1282 #undef __FUNCT__
1283 #define __FUNCT__ "DMPlexAddBoundary"
1284 /* The ids can be overridden by the command line option -bc_<boundary name> */
1285 PetscErrorCode DMPlexAddBoundary(DM dm, PetscBool isEssential, const char name[], PetscInt field, void (*bcFunc)(), PetscInt numids, const PetscInt *ids, void *ctx)
1286 {
1287   DM_Plex       *mesh = (DM_Plex *) dm->data;
1288   DMBoundary     b;
1289   PetscErrorCode ierr;
1290 
1291   PetscFunctionBegin;
1292   ierr = PetscNew(&b);CHKERRQ(ierr);
1293   ierr = PetscStrallocpy(name, (char **) &b->name);CHKERRQ(ierr);
1294   ierr = PetscMalloc1(numids, &b->ids);CHKERRQ(ierr);
1295   ierr = PetscMemcpy(b->ids, ids, numids*sizeof(PetscInt));CHKERRQ(ierr);
1296   b->essential   = isEssential;
1297   b->field       = field;
1298   b->func        = bcFunc;
1299   b->numids      = numids;
1300   b->ctx         = ctx;
1301   b->next        = mesh->boundary;
1302   mesh->boundary = b;
1303   PetscFunctionReturn(0);
1304 }
1305 
1306 #undef __FUNCT__
1307 #define __FUNCT__ "DMPlexGetNumBoundary"
1308 PetscErrorCode DMPlexGetNumBoundary(DM dm, PetscInt *numBd)
1309 {
1310   DM_Plex   *mesh = (DM_Plex *) dm->data;
1311   DMBoundary b    = mesh->boundary;
1312 
1313   PetscFunctionBegin;
1314   *numBd = 0;
1315   while (b) {++(*numBd); b = b->next;}
1316   PetscFunctionReturn(0);
1317 }
1318 
1319 #undef __FUNCT__
1320 #define __FUNCT__ "DMPlexGetBoundary"
1321 PetscErrorCode DMPlexGetBoundary(DM dm, PetscInt bd, PetscBool *isEssential, const char **name, PetscInt *field, void (**func)(), PetscInt *numids, const PetscInt **ids, void **ctx)
1322 {
1323   DM_Plex   *mesh = (DM_Plex *) dm->data;
1324   DMBoundary b    = mesh->boundary;
1325   PetscInt   n    = 0;
1326 
1327   PetscFunctionBegin;
1328   while (b) {
1329     if (n == bd) break;
1330     b = b->next;
1331     ++n;
1332   }
1333   if (n != bd) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Boundary %d is not in [0, %d)", bd, n);
1334   if (isEssential) {
1335     PetscValidPointer(isEssential, 3);
1336     *isEssential = b->essential;
1337   }
1338   if (name) {
1339     PetscValidPointer(name, 4);
1340     *name = b->name;
1341   }
1342   if (field) {
1343     PetscValidPointer(field, 5);
1344     *field = b->field;
1345   }
1346   if (func) {
1347     PetscValidPointer(func, 6);
1348     *func = b->func;
1349   }
1350   if (numids) {
1351     PetscValidPointer(numids, 7);
1352     *numids = b->numids;
1353   }
1354   if (ids) {
1355     PetscValidPointer(ids, 8);
1356     *ids = b->ids;
1357   }
1358   if (ctx) {
1359     PetscValidPointer(ctx, 9);
1360     *ctx = b->ctx;
1361   }
1362   PetscFunctionReturn(0);
1363 }
1364