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