xref: /petsc/src/dm/impls/plex/plexfem.c (revision f9426fe092dba0ba2fdf65dfec8d938c4b10a31c)
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 = DMPlexGetCoordinateSection(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 *), 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;
192   PetscErrorCode  ierr;
193 
194   PetscFunctionBegin;
195   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
196   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
197   ierr = PetscMalloc(numFields * sizeof(PetscDualSpace), &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,PetscReal,&v0,dim*dim,PetscReal,&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       ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr);
219       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
220       for (d = 0; d < spDim; ++d) {
221         ierr = PetscDualSpaceApply(sp[f], d, geom, numComp, funcs[f], &values[v]);CHKERRQ(ierr);
222         v += numComp;
223       }
224     }
225     ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr);
226   }
227   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
228   ierr = PetscFree(sp);CHKERRQ(ierr);
229   PetscFunctionReturn(0);
230 }
231 
232 #undef __FUNCT__
233 #define __FUNCT__ "DMPlexProjectFunction"
234 /*@C
235   DMPlexProjectFunction - This projects the given function into the function space provided.
236 
237   Input Parameters:
238 + dm      - The DM
239 . fe      - The PetscFE associated with the field
240 . funcs   - The coordinate functions to evaluate, one per field
241 - mode    - The insertion mode for values
242 
243   Output Parameter:
244 . X - vector
245 
246   Level: developer
247 
248 .seealso: DMPlexComputeL2Diff()
249 @*/
250 PetscErrorCode DMPlexProjectFunction(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *), InsertMode mode, Vec X)
251 {
252   Vec            localX;
253   PetscErrorCode ierr;
254 
255   PetscFunctionBegin;
256   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
257   ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, mode, localX);CHKERRQ(ierr);
258   ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr);
259   ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr);
260   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
261   PetscFunctionReturn(0);
262 }
263 
264 #undef __FUNCT__
265 #define __FUNCT__ "DMPlexComputeL2Diff"
266 /*@C
267   DMPlexComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.
268 
269   Input Parameters:
270 + dm    - The DM
271 . fe    - The PetscFE object for each field
272 . funcs - The functions to evaluate for each field component
273 - X     - The coefficient vector u_h
274 
275   Output Parameter:
276 . diff - The diff ||u - u_h||_2
277 
278   Level: developer
279 
280 .seealso: DMPlexProjectFunction()
281 @*/
282 PetscErrorCode DMPlexComputeL2Diff(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *), Vec X, PetscReal *diff)
283 {
284   const PetscInt  debug = 0;
285   PetscSection    section;
286   PetscQuadrature quad;
287   Vec             localX;
288   PetscScalar    *funcVal;
289   PetscReal      *coords, *v0, *J, *invJ, detJ;
290   PetscReal       localDiff = 0.0;
291   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp;
292   PetscErrorCode  ierr;
293 
294   PetscFunctionBegin;
295   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
296   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
297   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
298   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
299   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
300   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
301   for (field = 0; field < numFields; ++field) {
302     PetscInt Nc;
303 
304     ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr);
305     numComponents += Nc;
306   }
307   ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
308   ierr = PetscMalloc5(numComponents,PetscScalar,&funcVal,dim,PetscReal,&coords,dim,PetscReal,&v0,dim*dim,PetscReal,&J,dim*dim,PetscReal,&invJ);CHKERRQ(ierr);
309   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
310   ierr = PetscFEGetQuadrature(fe[0], &quad);CHKERRQ(ierr);
311   for (c = cStart; c < cEnd; ++c) {
312     PetscScalar *x = NULL;
313     PetscReal    elemDiff = 0.0;
314 
315     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
316     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
317     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
318 
319     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
320       const PetscInt   numQuadPoints = quad.numQuadPoints;
321       const PetscReal *quadPoints    = quad.quadPoints;
322       const PetscReal *quadWeights   = quad.quadWeights;
323       PetscReal       *basis;
324       PetscInt         numBasisFuncs, numBasisComps, q, d, e, fc, f;
325 
326       ierr = PetscFEGetDimension(fe[field], &numBasisFuncs);CHKERRQ(ierr);
327       ierr = PetscFEGetNumComponents(fe[field], &numBasisComps);CHKERRQ(ierr);
328       ierr = PetscFEGetDefaultTabulation(fe[field], &basis, NULL, NULL);CHKERRQ(ierr);
329       if (debug) {
330         char title[1024];
331         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
332         ierr = DMPrintCellVector(c, title, numBasisFuncs*numBasisComps, &x[fieldOffset]);CHKERRQ(ierr);
333       }
334       for (q = 0; q < numQuadPoints; ++q) {
335         for (d = 0; d < dim; d++) {
336           coords[d] = v0[d];
337           for (e = 0; e < dim; e++) {
338             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
339           }
340         }
341         (*funcs[field])(coords, funcVal);
342         for (fc = 0; fc < numBasisComps; ++fc) {
343           PetscScalar interpolant = 0.0;
344 
345           for (f = 0; f < numBasisFuncs; ++f) {
346             const PetscInt fidx = f*numBasisComps+fc;
347             interpolant += x[fieldOffset+fidx]*basis[q*numBasisFuncs*numBasisComps+fidx];
348           }
349           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);}
350           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ;
351         }
352       }
353       comp        += numBasisComps;
354       fieldOffset += numBasisFuncs*numBasisComps;
355     }
356     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
357     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
358     localDiff += elemDiff;
359   }
360   ierr  = PetscFree5(funcVal,coords,v0,J,invJ);CHKERRQ(ierr);
361   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
362   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
363   *diff = PetscSqrtReal(*diff);
364   PetscFunctionReturn(0);
365 }
366 
367 #if 0
368 
369 #undef __FUNCT__
370 #define __FUNCT__ "DMPlexComputeResidualFEM"
371 PetscErrorCode DMPlexComputeResidualFEM(DM dm, Vec X, Vec F, void *user)
372 {
373   DM_Plex         *mesh   = (DM_Plex *) dm->data;
374   PetscFEM        *fem    = (PetscFEM *) user;
375   PetscQuadrature *quad   = fem->quad;
376   PetscQuadrature *quadBd = fem->quadBd;
377   PetscSection     section;
378   PetscReal       *v0, *n, *J, *invJ, *detJ;
379   PetscScalar     *elemVec, *u;
380   PetscInt         dim, numFields, field, numBatchesTmp = 1, numCells, cStart, cEnd, c;
381   PetscInt         cellDof, numComponents;
382   PetscBool        has;
383   PetscErrorCode   ierr;
384 
385   PetscFunctionBegin;
386   if (has && quadBd) {
387     DMLabel         label;
388     IS              pointIS;
389     const PetscInt *points;
390     PetscInt        numPoints, p;
391 
392     ierr = DMPlexGetLabel(dm, "boundary", &label);CHKERRQ(ierr);
393     ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr);
394     ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr);
395     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
396     for (field = 0, cellDof = 0, numComponents = 0; field < numFields; ++field) {
397       cellDof       += quadBd[field].numBasisFuncs*quadBd[field].numComponents;
398       numComponents += quadBd[field].numComponents;
399     }
400     ierr = PetscMalloc7(numPoints*cellDof,PetscScalar,&u,numPoints*dim,PetscReal,&v0,numPoints*dim,PetscReal,&n,numPoints*dim*dim,PetscReal,&J,numPoints*dim*dim,PetscReal,&invJ,numPoints,PetscReal,&detJ,numPoints*cellDof,PetscScalar,&elemVec);CHKERRQ(ierr);
401     for (p = 0; p < numPoints; ++p) {
402       const PetscInt point = points[p];
403       PetscScalar   *x;
404       PetscInt       i;
405 
406       /* TODO: Add normal determination here */
407       ierr = DMPlexComputeCellGeometry(dm, point, &v0[p*dim], &J[p*dim*dim], &invJ[p*dim*dim], &detJ[p]);CHKERRQ(ierr);
408       if (detJ[p] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[p], point);
409       ierr = DMPlexVecGetClosure(dm, NULL, X, point, NULL, &x);CHKERRQ(ierr);
410 
411       for (i = 0; i < cellDof; ++i) u[p*cellDof+i] = x[i];
412       ierr = DMPlexVecRestoreClosure(dm, NULL, X, point, NULL, &x);CHKERRQ(ierr);
413     }
414     for (field = 0; field < numFields; ++field) {
415       const PetscInt numQuadPoints = quadBd[field].numQuadPoints;
416       const PetscInt numBasisFuncs = quadBd[field].numBasisFuncs;
417       void           (*f0)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f0BdFuncs[field];
418       void           (*f1)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f1BdFuncs[field];
419       /* Conforming batches */
420       PetscInt blockSize  = numBasisFuncs*numQuadPoints;
421       PetscInt numBlocks  = 1;
422       PetscInt batchSize  = numBlocks * blockSize;
423       PetscInt numBatches = numBatchesTmp;
424       PetscInt numChunks  = numPoints / (numBatches*batchSize);
425       /* Remainder */
426       PetscInt numRemainder = numPoints % (numBatches * batchSize);
427       PetscInt offset       = numPoints - numRemainder;
428 
429       ierr = (*mesh->integrateBdResidualFEM)(numChunks*numBatches*batchSize, numFields, field, quadBd, u, v0, n, J, invJ, detJ, f0, f1, elemVec);CHKERRQ(ierr);
430       ierr = (*mesh->integrateBdResidualFEM)(numRemainder, numFields, field, quadBd, &u[offset*cellDof], &v0[offset*dim], &n[offset*dim], &J[offset*dim*dim], &invJ[offset*dim*dim], &detJ[offset],
431                                              f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr);
432     }
433     for (p = 0; p < numPoints; ++p) {
434       const PetscInt point = points[p];
435 
436       if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "Residual", cellDof, &elemVec[p*cellDof]);CHKERRQ(ierr);}
437       ierr = DMPlexVecSetClosure(dm, NULL, F, point, &elemVec[p*cellDof], ADD_VALUES);CHKERRQ(ierr);
438     }
439     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
440     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
441     ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemVec);CHKERRQ(ierr);
442   }
443   PetscFunctionReturn(0);
444 }
445 
446 #else
447 
448 #undef __FUNCT__
449 #define __FUNCT__ "DMPlexComputeResidualFEM"
450 /*@
451   DMPlexComputeResidualFEM - Form the local residual F from the local input X using pointwise functions specified by the user
452 
453   Input Parameters:
454 + dm - The mesh
455 . X  - Local input vector
456 - user - The user context
457 
458   Output Parameter:
459 . F  - Local output vector
460 
461   Note:
462   The second member of the user context must be an FEMContext.
463 
464   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
465   like a GPU, or vectorize on a multicore machine.
466 
467   Level: developer
468 
469 .seealso: DMPlexComputeJacobianActionFEM()
470 @*/
471 PetscErrorCode DMPlexComputeResidualFEM(DM dm, Vec X, Vec F, void *user)
472 {
473   DM_Plex          *mesh  = (DM_Plex *) dm->data;
474   PetscFEM         *fem   = (PetscFEM *) user;
475   PetscFE          *fe    = fem->fe;
476   PetscFE          *feAux = fem->feAux;
477   PetscFE          *feBd  = fem->feBd;
478   const char       *name  = "Residual";
479   DM                dmAux;
480   Vec               A;
481   PetscQuadrature   q;
482   PetscCellGeometry geom;
483   PetscSection      section, sectionAux;
484   PetscReal        *v0, *J, *invJ, *detJ;
485   PetscScalar      *elemVec, *u, *a;
486   PetscInt          dim, Nf, NfAux = 0, f, numCells, cStart, cEnd, c;
487   PetscInt          cellDof = 0, numComponents = 0;
488   PetscInt          cellDofAux = 0, numComponentsAux = 0;
489   PetscErrorCode    ierr;
490 
491   PetscFunctionBegin;
492   ierr = PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
493   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
494   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
495   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
496   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
497   numCells = cEnd - cStart;
498   for (f = 0; f < Nf; ++f) {
499     PetscInt Nb, Nc;
500 
501     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
502     ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr);
503     cellDof       += Nb*Nc;
504     numComponents += Nc;
505   }
506   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
507   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
508   if (dmAux) {
509     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
510     ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr);
511   }
512   for (f = 0; f < NfAux; ++f) {
513     PetscInt Nb, Nc;
514 
515     ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr);
516     ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr);
517     cellDofAux       += Nb*Nc;
518     numComponentsAux += Nc;
519   }
520   ierr = DMPlexProjectFunctionLocal(dm, fe, fem->bcFuncs, INSERT_BC_VALUES, X);CHKERRQ(ierr);
521   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
522   ierr = PetscMalloc6(numCells*cellDof,PetscScalar,&u,numCells*dim,PetscReal,&v0,numCells*dim*dim,PetscReal,&J,numCells*dim*dim,PetscReal,&invJ,numCells,PetscReal,&detJ,numCells*cellDof,PetscScalar,&elemVec);CHKERRQ(ierr);
523   if (dmAux) {ierr = PetscMalloc(numCells*cellDofAux * sizeof(PetscScalar), &a);CHKERRQ(ierr);}
524   for (c = cStart; c < cEnd; ++c) {
525     PetscScalar *x = NULL;
526     PetscInt     i;
527 
528     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
529     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
530     ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
531     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
532     ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
533     if (dmAux) {
534       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
535       for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x[i];
536       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
537     }
538   }
539   for (f = 0; f < Nf; ++f) {
540     void   (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f0Funcs[f];
541     void   (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f1Funcs[f];
542     PetscInt Nb;
543     /* Conforming batches */
544     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
545     /* Remainder */
546     PetscInt Nr, offset;
547 
548     ierr = PetscFEGetQuadrature(fe[f], &q);CHKERRQ(ierr);
549     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
550     ierr = PetscFEGetTileSizes(fe[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
551     blockSize = Nb*q.numQuadPoints;
552     batchSize = numBlocks * blockSize;
553     ierr =  PetscFESetTileSizes(fe[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
554     numChunks = numCells / (numBatches*batchSize);
555     Ne        = numChunks*numBatches*batchSize;
556     Nr        = numCells % (numBatches*batchSize);
557     offset    = numCells - Nr;
558     geom.v0   = v0;
559     geom.J    = J;
560     geom.invJ = invJ;
561     geom.detJ = detJ;
562     ierr = PetscFEIntegrateResidual(fe[f], Ne, Nf, fe, f, geom, u, NfAux, feAux, a, f0, f1, elemVec);CHKERRQ(ierr);
563     geom.v0   = &v0[offset*dim];
564     geom.J    = &J[offset*dim*dim];
565     geom.invJ = &invJ[offset*dim*dim];
566     geom.detJ = &detJ[offset];
567     ierr = PetscFEIntegrateResidual(fe[f], Nr, Nf, fe, f, geom, &u[offset*cellDof], NfAux, feAux, &a[offset*cellDofAux], f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr);
568   }
569   for (c = cStart; c < cEnd; ++c) {
570     if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, name, cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);}
571     ierr = DMPlexVecSetClosure(dm, section, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr);
572   }
573   ierr = PetscFree6(u,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr);
574   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
575   if (feBd) {
576     DMLabel         label;
577     IS              pointIS;
578     const PetscInt *points;
579     PetscInt        numPoints, p;
580     PetscReal      *n;
581 
582     ierr = DMPlexGetLabel(dm, "boundary", &label);CHKERRQ(ierr);
583     ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr);
584     ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr);
585     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
586     for (f = 0, cellDof = 0, numComponents = 0; f < Nf; ++f) {
587       PetscInt Nb, Nc;
588 
589       ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr);
590       ierr = PetscFEGetNumComponents(feBd[f], &Nc);CHKERRQ(ierr);
591       cellDof       += Nb*Nc;
592       numComponents += Nc;
593     }
594     ierr = PetscMalloc7(numPoints*cellDof,PetscScalar,&u,numPoints*dim,PetscReal,&v0,numPoints*dim,PetscReal,&n,numPoints*dim*dim,PetscReal,&J,numPoints*dim*dim,PetscReal,&invJ,numPoints,PetscReal,&detJ,numPoints*cellDof,PetscScalar,&elemVec);CHKERRQ(ierr);
595     for (p = 0; p < numPoints; ++p) {
596       const PetscInt point = points[p];
597       PetscScalar   *x     = NULL;
598       PetscInt       i;
599 
600       /* TODO: Add normal determination here */
601       ierr = DMPlexComputeCellGeometry(dm, point, &v0[p*dim], &J[p*dim*dim], &invJ[p*dim*dim], &detJ[p]);CHKERRQ(ierr);
602       if (detJ[p] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[p], point);
603       ierr = DMPlexVecGetClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
604       for (i = 0; i < cellDof; ++i) u[p*cellDof+i] = x[i];
605       ierr = DMPlexVecRestoreClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr);
606     }
607     for (f = 0; f < Nf; ++f) {
608       void   (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f0BdFuncs[f];
609       void   (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f1BdFuncs[f];
610       PetscInt Nb;
611       /* Conforming batches */
612       PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
613       /* Remainder */
614       PetscInt Nr, offset;
615 
616       ierr = PetscFEGetQuadrature(feBd[f], &q);CHKERRQ(ierr);
617       ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr);
618       ierr = PetscFEGetTileSizes(feBd[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
619       blockSize = Nb*q.numQuadPoints;
620       batchSize = numBlocks * blockSize;
621       ierr =  PetscFESetTileSizes(feBd[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
622       numChunks = numPoints / (numBatches*batchSize);
623       Ne        = numChunks*numBatches*batchSize;
624       Nr        = numPoints % (numBatches*batchSize);
625       offset    = numPoints - Nr;
626       geom.v0   = v0;
627       geom.n    = n;
628       geom.J    = J;
629       geom.invJ = invJ;
630       geom.detJ = detJ;
631       ierr = PetscFEIntegrateBdResidual(feBd[f], Ne, Nf, feBd, f, geom, u, 0, NULL, NULL, f0, f1, elemVec);CHKERRQ(ierr);
632       geom.v0   = &v0[offset*dim];
633       geom.n    = &n[offset*dim];
634       geom.J    = &J[offset*dim*dim];
635       geom.invJ = &invJ[offset*dim*dim];
636       geom.detJ = &detJ[offset];
637       ierr = PetscFEIntegrateBdResidual(feBd[f], Nr, Nf, feBd, f, geom, &u[offset*cellDof], 0, NULL, NULL, f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr);
638     }
639     for (p = 0; p < numPoints; ++p) {
640       const PetscInt point = points[p];
641 
642       if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "BdResidual", cellDof, &elemVec[p*cellDof]);CHKERRQ(ierr);}
643       ierr = DMPlexVecSetClosure(dm, NULL, F, point, &elemVec[p*cellDof], ADD_VALUES);CHKERRQ(ierr);
644     }
645     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
646     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
647     ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemVec);CHKERRQ(ierr);
648   }
649   if (mesh->printFEM) {ierr = DMPrintLocalVec(dm, name, mesh->printTol, F);CHKERRQ(ierr);}
650   ierr = PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr);
651   PetscFunctionReturn(0);
652 }
653 
654 #endif
655 
656 #undef __FUNCT__
657 #define __FUNCT__ "DMPlexComputeJacobianActionFEM"
658 /*@C
659   DMPlexComputeJacobianActionFEM - Form the local action of Jacobian J(u) on the local input X using pointwise functions specified by the user
660 
661   Input Parameters:
662 + dm - The mesh
663 . J  - The Jacobian shell matrix
664 . X  - Local input vector
665 - user - The user context
666 
667   Output Parameter:
668 . F  - Local output vector
669 
670   Note:
671   The second member of the user context must be an FEMContext.
672 
673   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
674   like a GPU, or vectorize on a multicore machine.
675 
676   Level: developer
677 
678 .seealso: DMPlexComputeResidualFEM()
679 @*/
680 PetscErrorCode DMPlexComputeJacobianActionFEM(DM dm, Mat Jac, Vec X, Vec F, void *user)
681 {
682   DM_Plex          *mesh = (DM_Plex *) dm->data;
683   PetscFEM         *fem  = (PetscFEM *) user;
684   PetscFE          *fe   = fem->fe;
685   PetscQuadrature   quad;
686   PetscCellGeometry geom;
687   PetscSection      section;
688   JacActionCtx     *jctx;
689   PetscReal        *v0, *J, *invJ, *detJ;
690   PetscScalar      *elemVec, *u, *a;
691   PetscInt          dim, numFields, field, numCells, cStart, cEnd, c;
692   PetscInt          cellDof = 0;
693   PetscErrorCode    ierr;
694 
695   PetscFunctionBegin;
696   /* ierr = PetscLogEventBegin(DMPLEX_JacobianActionFEM,dm,0,0,0);CHKERRQ(ierr); */
697   ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr);
698   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
699   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
700   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
701   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
702   numCells = cEnd - cStart;
703   for (field = 0; field < numFields; ++field) {
704     PetscInt Nb, Nc;
705 
706     ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr);
707     ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr);
708     cellDof += Nb*Nc;
709   }
710   ierr = VecSet(F, 0.0);CHKERRQ(ierr);
711   ierr = PetscMalloc7(numCells*cellDof,PetscScalar,&u,numCells*cellDof,PetscScalar,&a,numCells*dim,PetscReal,&v0,numCells*dim*dim,PetscReal,&J,numCells*dim*dim,PetscReal,&invJ,numCells,PetscReal,&detJ,numCells*cellDof,PetscScalar,&elemVec);CHKERRQ(ierr);
712   for (c = cStart; c < cEnd; ++c) {
713     PetscScalar *x = NULL;
714     PetscInt     i;
715 
716     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
717     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
718     ierr = DMPlexVecGetClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr);
719     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
720     ierr = DMPlexVecRestoreClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr);
721     ierr = DMPlexVecGetClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr);
722     for (i = 0; i < cellDof; ++i) a[c*cellDof+i] = x[i];
723     ierr = DMPlexVecRestoreClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr);
724   }
725   for (field = 0; field < numFields; ++field) {
726     PetscInt Nb;
727     /* Conforming batches */
728     PetscInt numBlocks  = 1;
729     PetscInt numBatches = 1;
730     PetscInt numChunks, Ne, blockSize, batchSize;
731     /* Remainder */
732     PetscInt Nr, offset;
733 
734     ierr = PetscFEGetQuadrature(fe[field], &quad);CHKERRQ(ierr);
735     ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr);
736     blockSize = Nb*quad.numQuadPoints;
737     batchSize = numBlocks * blockSize;
738     numChunks = numCells / (numBatches*batchSize);
739     Ne        = numChunks*numBatches*batchSize;
740     Nr        = numCells % (numBatches*batchSize);
741     offset    = numCells - Nr;
742     geom.v0   = v0;
743     geom.J    = J;
744     geom.invJ = invJ;
745     geom.detJ = detJ;
746     ierr = PetscFEIntegrateJacobianAction(fe[field], Ne, numFields, fe, field, geom, u, a, fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, elemVec);CHKERRQ(ierr);
747     geom.v0   = &v0[offset*dim];
748     geom.J    = &J[offset*dim*dim];
749     geom.invJ = &invJ[offset*dim*dim];
750     geom.detJ = &detJ[offset];
751     ierr = PetscFEIntegrateJacobianAction(fe[field], Nr, numFields, fe, field, geom, &u[offset*cellDof], &a[offset*cellDof],
752                                           fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, &elemVec[offset*cellDof]);CHKERRQ(ierr);
753   }
754   for (c = cStart; c < cEnd; ++c) {
755     if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, "Jacobian Action", cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);}
756     ierr = DMPlexVecSetClosure(dm, NULL, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr);
757   }
758   ierr = PetscFree7(u,a,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr);
759   if (mesh->printFEM) {
760     PetscMPIInt rank, numProcs;
761     PetscInt    p;
762 
763     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
764     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm), &numProcs);CHKERRQ(ierr);
765     ierr = PetscPrintf(PetscObjectComm((PetscObject)dm), "Jacobian Action:\n");CHKERRQ(ierr);
766     for (p = 0; p < numProcs; ++p) {
767       if (p == rank) {ierr = VecView(F, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);}
768       ierr = PetscBarrier((PetscObject) dm);CHKERRQ(ierr);
769     }
770   }
771   /* ierr = PetscLogEventEnd(DMPLEX_JacobianActionFEM,dm,0,0,0);CHKERRQ(ierr); */
772   PetscFunctionReturn(0);
773 }
774 
775 #undef __FUNCT__
776 #define __FUNCT__ "DMPlexComputeJacobianFEM"
777 /*@
778   DMPlexComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user.
779 
780   Input Parameters:
781 + dm - The mesh
782 . X  - Local input vector
783 - user - The user context
784 
785   Output Parameter:
786 . Jac  - Jacobian matrix
787 
788   Note:
789   The second member of the user context must be an FEMContext.
790 
791   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
792   like a GPU, or vectorize on a multicore machine.
793 
794   Level: developer
795 
796 .seealso: FormFunctionLocal()
797 @*/
798 PetscErrorCode DMPlexComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP, MatStructure *str,void *user)
799 {
800   DM_Plex          *mesh  = (DM_Plex *) dm->data;
801   PetscFEM         *fem   = (PetscFEM *) user;
802   PetscFE          *fe    = fem->fe;
803   PetscFE          *feAux = fem->feAux;
804   const char       *name  = "Jacobian";
805   DM                dmAux;
806   Vec               A;
807   PetscQuadrature   quad;
808   PetscCellGeometry geom;
809   PetscSection      section, globalSection, sectionAux;
810   PetscReal        *v0, *J, *invJ, *detJ;
811   PetscScalar      *elemMat, *u, *a;
812   PetscInt          dim, Nf, NfAux = 0, f, fieldI, fieldJ, numCells, cStart, cEnd, c;
813   PetscInt          cellDof = 0, numComponents = 0;
814   PetscInt          cellDofAux = 0, numComponentsAux = 0;
815   PetscBool         isShell;
816   PetscErrorCode    ierr;
817 
818   PetscFunctionBegin;
819   ierr = PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
820   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
821   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
822   ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);
823   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
824   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
825   numCells = cEnd - cStart;
826   for (f = 0; f < Nf; ++f) {
827     PetscInt Nb, Nc;
828 
829     ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr);
830     ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr);
831     cellDof       += Nb*Nc;
832     numComponents += Nc;
833   }
834   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
835   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
836   if (dmAux) {
837     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
838     ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr);
839   }
840   for (f = 0; f < NfAux; ++f) {
841     PetscInt Nb, Nc;
842 
843     ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr);
844     ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr);
845     cellDofAux       += Nb*Nc;
846     numComponentsAux += Nc;
847   }
848   ierr = DMPlexProjectFunctionLocal(dm, fe, fem->bcFuncs, INSERT_BC_VALUES, X);CHKERRQ(ierr);
849   ierr = MatZeroEntries(JacP);CHKERRQ(ierr);
850   ierr = PetscMalloc6(numCells*cellDof,PetscScalar,&u,numCells*dim,PetscReal,&v0,numCells*dim*dim,PetscReal,&J,numCells*dim*dim,PetscReal,&invJ,numCells,PetscReal,&detJ,numCells*cellDof*cellDof,PetscScalar,&elemMat);CHKERRQ(ierr);
851   if (dmAux) {ierr = PetscMalloc(numCells*cellDofAux * sizeof(PetscScalar), &a);CHKERRQ(ierr);}
852   for (c = cStart; c < cEnd; ++c) {
853     PetscScalar *x = NULL;
854     PetscInt     i;
855 
856     ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr);
857     if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c);
858     ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
859     for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i];
860     ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr);
861     if (dmAux) {
862       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
863       for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x[i];
864       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
865     }
866   }
867   ierr = PetscMemzero(elemMat, numCells*cellDof*cellDof * sizeof(PetscScalar));CHKERRQ(ierr);
868   for (fieldI = 0; fieldI < Nf; ++fieldI) {
869     PetscInt Nb;
870     /* Conforming batches */
871     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
872     /* Remainder */
873     PetscInt Nr, offset;
874 
875     ierr = PetscFEGetQuadrature(fe[fieldI], &quad);CHKERRQ(ierr);
876     ierr = PetscFEGetDimension(fe[fieldI], &Nb);CHKERRQ(ierr);
877     ierr = PetscFEGetTileSizes(fe[fieldI], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
878     blockSize = Nb*quad.numQuadPoints;
879     batchSize = numBlocks * blockSize;
880     ierr = PetscFESetTileSizes(fe[fieldI], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
881     numChunks = numCells / (numBatches*batchSize);
882     Ne        = numChunks*numBatches*batchSize;
883     Nr        = numCells % (numBatches*batchSize);
884     offset    = numCells - Nr;
885     for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
886       void   (*g0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g0Funcs[fieldI*Nf+fieldJ];
887       void   (*g1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g1Funcs[fieldI*Nf+fieldJ];
888       void   (*g2)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g2Funcs[fieldI*Nf+fieldJ];
889       void   (*g3)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g3Funcs[fieldI*Nf+fieldJ];
890 
891       geom.v0   = v0;
892       geom.J    = J;
893       geom.invJ = invJ;
894       geom.detJ = detJ;
895       ierr = PetscFEIntegrateJacobian(fe[fieldI], Ne, Nf, fe, fieldI, fieldJ, geom, u, NfAux, feAux, a, g0, g1, g2, g3, elemMat);CHKERRQ(ierr);
896       geom.v0   = &v0[offset*dim];
897       geom.J    = &J[offset*dim*dim];
898       geom.invJ = &invJ[offset*dim*dim];
899       geom.detJ = &detJ[offset];
900       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);
901     }
902   }
903   for (c = cStart; c < cEnd; ++c) {
904     if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, name, cellDof, cellDof, &elemMat[c*cellDof*cellDof]);CHKERRQ(ierr);}
905     ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, c, &elemMat[c*cellDof*cellDof], ADD_VALUES);CHKERRQ(ierr);
906   }
907   ierr = PetscFree6(u,v0,J,invJ,detJ,elemMat);CHKERRQ(ierr);
908   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
909   ierr = MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
910   ierr = MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
911   if (mesh->printFEM) {
912     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
913     ierr = MatChop(JacP, 1.0e-10);CHKERRQ(ierr);
914     ierr = MatView(JacP, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
915   }
916   ierr = PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
917   ierr = PetscObjectTypeCompare((PetscObject) Jac, MATSHELL, &isShell);CHKERRQ(ierr);
918   if (isShell) {
919     JacActionCtx *jctx;
920 
921     ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr);
922     ierr = VecCopy(X, jctx->u);CHKERRQ(ierr);
923   }
924   *str = SAME_NONZERO_PATTERN;
925   PetscFunctionReturn(0);
926 }
927