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