Actual source code: plexsection.c
petsc-3.11.4 2019-09-28
1: #include <petsc/private/dmpleximpl.h>
3: /* Set the number of dof on each point and separate by fields */
4: static PetscErrorCode DMPlexCreateSectionFields(DM dm, const PetscInt numComp[], PetscSection *section)
5: {
6: DMLabel depthLabel;
7: PetscInt depth, Nf, f, pStart, pEnd;
8: PetscBool *isFE;
12: DMPlexGetDepth(dm, &depth);
13: DMPlexGetDepthLabel(dm,&depthLabel);
14: DMGetNumFields(dm, &Nf);
15: PetscCalloc1(Nf, &isFE);
16: for (f = 0; f < Nf; ++f) {
17: PetscObject obj;
18: PetscClassId id;
20: DMGetField(dm, f, NULL, &obj);
21: PetscObjectGetClassId(obj, &id);
22: if (id == PETSCFE_CLASSID) {isFE[f] = PETSC_TRUE;}
23: else if (id == PETSCFV_CLASSID) {isFE[f] = PETSC_FALSE;}
24: }
26: PetscSectionCreate(PetscObjectComm((PetscObject)dm), section);
27: if (Nf > 0) {
28: PetscSectionSetNumFields(*section, Nf);
29: if (numComp) {
30: for (f = 0; f < Nf; ++f) {
31: PetscSectionSetFieldComponents(*section, f, numComp[f]);
32: if (isFE[f]) {
33: PetscFE fe;
34: PetscDualSpace dspace;
35: const PetscInt ***perms;
36: const PetscScalar ***flips;
37: const PetscInt *numDof;
39: DMGetField(dm,f,NULL,(PetscObject *) &fe);
40: PetscFEGetDualSpace(fe,&dspace);
41: PetscDualSpaceGetSymmetries(dspace,&perms,&flips);
42: PetscDualSpaceGetNumDof(dspace,&numDof);
43: if (perms || flips) {
44: DM K;
45: PetscInt h;
46: PetscSectionSym sym;
48: PetscDualSpaceGetDM(dspace,&K);
49: PetscSectionSymCreateLabel(PetscObjectComm((PetscObject)*section),depthLabel,&sym);
50: for (h = 0; h <= depth; h++) {
51: PetscDualSpace hspace;
52: PetscInt kStart, kEnd;
53: PetscInt kConeSize;
54: const PetscInt **perms0 = NULL;
55: const PetscScalar **flips0 = NULL;
57: PetscDualSpaceGetHeightSubspace(dspace,h,&hspace);
58: DMPlexGetHeightStratum(K,h,&kStart,&kEnd);
59: if (!hspace) continue;
60: PetscDualSpaceGetSymmetries(hspace,&perms,&flips);
61: if (perms) perms0 = perms[0];
62: if (flips) flips0 = flips[0];
63: if (!(perms0 || flips0)) continue;
64: DMPlexGetConeSize(K,kStart,&kConeSize);
65: PetscSectionSymLabelSetStratum(sym,depth - h,numDof[depth - h],-kConeSize,kConeSize,PETSC_USE_POINTER,perms0 ? &perms0[-kConeSize] : NULL,flips0 ? &flips0[-kConeSize] : NULL);
66: }
67: PetscSectionSetFieldSym(*section,f,sym);
68: PetscSectionSymDestroy(&sym);
69: }
70: }
71: }
72: }
73: }
74: DMPlexGetChart(dm, &pStart, &pEnd);
75: PetscSectionSetChart(*section, pStart, pEnd);
76: PetscFree(isFE);
77: return(0);
78: }
80: /* Set the number of dof on each point and separate by fields */
81: static PetscErrorCode DMPlexCreateSectionDof(DM dm, DMLabel label[],const PetscInt numDof[], PetscSection section)
82: {
83: DMLabel depthLabel;
84: PetscInt *pMax;
85: PetscInt depth, cellHeight, pStart = 0, pEnd = 0;
86: PetscInt Nf, f, dim, d, dep, p;
87: PetscBool *isFE;
91: DMGetDimension(dm, &dim);
92: DMPlexGetDepth(dm, &depth);
93: DMPlexGetDepthLabel(dm,&depthLabel);
94: DMGetNumFields(dm, &Nf);
95: PetscMalloc1(Nf, &isFE);
96: for (f = 0; f < Nf; ++f) {
97: PetscObject obj;
98: PetscClassId id;
100: DMGetField(dm, f, NULL, &obj);
101: PetscObjectGetClassId(obj, &id);
102: if (id == PETSCFE_CLASSID) {isFE[f] = PETSC_TRUE;}
103: else if (id == PETSCFV_CLASSID) {isFE[f] = PETSC_FALSE;}
104: }
106: PetscMalloc1(depth+1, &pMax);
107: DMPlexGetHybridBounds(dm, depth >= 0 ? &pMax[depth] : NULL, depth>1 ? &pMax[depth-1] : NULL, depth>2 ? &pMax[1] : NULL, &pMax[0]);
108: DMPlexGetVTKCellHeight(dm, &cellHeight);
109: for (f = 0; f < Nf; ++f) {
110: if (label && label[f]) {
111: IS pointIS;
112: const PetscInt *points;
113: PetscInt n;
115: DMLabelGetStratumIS(label[f], 1, &pointIS);
116: if (!pointIS) continue;
117: ISGetLocalSize(pointIS, &n);
118: ISGetIndices(pointIS, &points);
119: for (p = 0; p < n; ++p) {
120: const PetscInt point = points[p];
121: PetscInt dof, d;
123: DMLabelGetValue(depthLabel, point, &d);
124: DMPlexGetDepthStratum(dm, d, NULL, &pEnd);
125: /* If this is a hybrid point, use dof for one dimension lower */
126: if ((point >= pMax[d]) && (point < pEnd)) --d;
127: dof = d < 0 ? 0 : numDof[f*(dim+1)+d];
128: PetscSectionSetFieldDof(section, point, f, dof);
129: PetscSectionAddDof(section, point, dof);
130: }
131: ISRestoreIndices(pointIS, &points);
132: ISDestroy(&pointIS);
133: } else {
134: for (dep = 0; dep <= depth - cellHeight; ++dep) {
135: d = dim == depth ? dep : (!dep ? 0 : dim);
136: DMPlexGetDepthStratum(dm, dep, &pStart, &pEnd);
137: pMax[dep] = pMax[dep] < 0 ? pEnd : pMax[dep];
138: for (p = pStart; p < pEnd; ++p) {
139: const PetscInt dof = numDof[f*(dim+1)+d];
141: if (isFE[f] && p >= pMax[dep]) continue;
142: PetscSectionSetFieldDof(section, p, f, dof);
143: PetscSectionAddDof(section, p, dof);
144: }
145: }
146: }
147: }
148: PetscFree(pMax);
149: PetscFree(isFE);
150: return(0);
151: }
153: /* Set the number of dof on each point and separate by fields
154: If bcComps is NULL or the IS is NULL, constrain every dof on the point
155: */
156: static PetscErrorCode DMPlexCreateSectionBCDof(DM dm, PetscInt numBC, const PetscInt bcField[], const IS bcComps[], const IS bcPoints[], PetscSection section)
157: {
158: PetscInt Nf;
159: PetscInt bc;
160: PetscSection aSec;
164: PetscSectionGetNumFields(section, &Nf);
165: for (bc = 0; bc < numBC; ++bc) {
166: PetscInt field = 0;
167: const PetscInt *comp;
168: const PetscInt *idx;
169: PetscInt Nc = -1, n, i;
171: if (Nf) field = bcField[bc];
172: if (bcComps && bcComps[bc]) {ISGetLocalSize(bcComps[bc], &Nc);}
173: if (bcComps && bcComps[bc]) {ISGetIndices(bcComps[bc], &comp);}
174: ISGetLocalSize(bcPoints[bc], &n);
175: ISGetIndices(bcPoints[bc], &idx);
176: for (i = 0; i < n; ++i) {
177: const PetscInt p = idx[i];
178: PetscInt numConst;
180: if (Nf) {
181: PetscSectionGetFieldDof(section, p, field, &numConst);
182: } else {
183: PetscSectionGetDof(section, p, &numConst);
184: }
185: /* If Nc < 0, constrain every dof on the point */
186: /* TODO: Matt, this only works if there is one node on the point. We need to handle numDofs > NumComponents */
187: if (Nc > 0) numConst = PetscMin(numConst, Nc);
188: if (Nf) {PetscSectionAddFieldConstraintDof(section, p, field, numConst);}
189: PetscSectionAddConstraintDof(section, p, numConst);
190: }
191: ISRestoreIndices(bcPoints[bc], &idx);
192: if (bcComps && bcComps[bc]) {ISRestoreIndices(bcComps[bc], &comp);}
193: }
194: DMPlexGetAnchors(dm, &aSec, NULL);
195: if (aSec) {
196: PetscInt aStart, aEnd, a;
198: PetscSectionGetChart(aSec, &aStart, &aEnd);
199: for (a = aStart; a < aEnd; a++) {
200: PetscInt dof, f;
202: PetscSectionGetDof(aSec, a, &dof);
203: if (dof) {
204: /* if there are point-to-point constraints, then all dofs are constrained */
205: PetscSectionGetDof(section, a, &dof);
206: PetscSectionSetConstraintDof(section, a, dof);
207: for (f = 0; f < Nf; f++) {
208: PetscSectionGetFieldDof(section, a, f, &dof);
209: PetscSectionSetFieldConstraintDof(section, a, f, dof);
210: }
211: }
212: }
213: }
214: return(0);
215: }
217: /* Set the constrained field indices on each point
218: If bcComps is NULL or the IS is NULL, constrain every dof on the point
219: */
220: static PetscErrorCode DMPlexCreateSectionBCIndicesField(DM dm, PetscInt numBC,const PetscInt bcField[], const IS bcComps[], const IS bcPoints[], PetscSection section)
221: {
222: PetscSection aSec;
223: PetscInt *indices;
224: PetscInt Nf, cdof, maxDof = 0, pStart, pEnd, p, bc, f, d;
228: PetscSectionGetNumFields(section, &Nf);
229: if (!Nf) return(0);
230: /* Initialize all field indices to -1 */
231: PetscSectionGetChart(section, &pStart, &pEnd);
232: for (p = pStart; p < pEnd; ++p) {PetscSectionGetConstraintDof(section, p, &cdof); maxDof = PetscMax(maxDof, cdof);}
233: PetscMalloc1(maxDof, &indices);
234: for (d = 0; d < maxDof; ++d) indices[d] = -1;
235: for (p = pStart; p < pEnd; ++p) for (f = 0; f < Nf; ++f) {PetscSectionSetFieldConstraintIndices(section, p, f, indices);}
236: /* Handle BC constraints */
237: for (bc = 0; bc < numBC; ++bc) {
238: const PetscInt field = bcField[bc];
239: const PetscInt *comp, *idx;
240: PetscInt Nc = -1, n, i;
242: if (bcComps && bcComps[bc]) {ISGetLocalSize(bcComps[bc], &Nc);}
243: if (bcComps && bcComps[bc]) {ISGetIndices(bcComps[bc], &comp);}
244: ISGetLocalSize(bcPoints[bc], &n);
245: ISGetIndices(bcPoints[bc], &idx);
246: for (i = 0; i < n; ++i) {
247: const PetscInt p = idx[i];
248: const PetscInt *find;
249: PetscInt fdof, fcdof, c;
251: PetscSectionGetFieldDof(section, p, field, &fdof);
252: if (!fdof) continue;
253: if (Nc < 0) {
254: for (d = 0; d < fdof; ++d) indices[d] = d;
255: fcdof = fdof;
256: } else {
257: PetscSectionGetFieldConstraintDof(section, p, field, &fcdof);
258: PetscSectionGetFieldConstraintIndices(section, p, field, &find);
259: for (d = 0; d < fcdof; ++d) {if (find[d] < 0) break; indices[d] = find[d];}
260: for (c = 0; c < Nc; ++c) indices[d++] = comp[c];
261: PetscSortRemoveDupsInt(&d, indices);
262: for (c = d; c < fcdof; ++c) indices[c] = -1;
263: fcdof = d;
264: }
265: PetscSectionSetFieldConstraintDof(section, p, field, fcdof);
266: PetscSectionSetFieldConstraintIndices(section, p, field, indices);
267: }
268: if (bcComps && bcComps[bc]) {ISRestoreIndices(bcComps[bc], &comp);}
269: ISRestoreIndices(bcPoints[bc], &idx);
270: }
271: /* Handle anchors */
272: DMPlexGetAnchors(dm, &aSec, NULL);
273: if (aSec) {
274: PetscInt aStart, aEnd, a;
276: for (d = 0; d < maxDof; ++d) indices[d] = d;
277: PetscSectionGetChart(aSec, &aStart, &aEnd);
278: for (a = aStart; a < aEnd; a++) {
279: PetscInt dof, f;
281: PetscSectionGetDof(aSec, a, &dof);
282: if (dof) {
283: /* if there are point-to-point constraints, then all dofs are constrained */
284: for (f = 0; f < Nf; f++) {
285: PetscSectionSetFieldConstraintIndices(section, a, f, indices);
286: }
287: }
288: }
289: }
290: PetscFree(indices);
291: return(0);
292: }
294: /* Set the constrained indices on each point */
295: static PetscErrorCode DMPlexCreateSectionBCIndices(DM dm, PetscSection section)
296: {
297: PetscInt *indices;
298: PetscInt Nf, maxDof, pStart, pEnd, p, f, d;
302: PetscSectionGetNumFields(section, &Nf);
303: PetscSectionGetMaxDof(section, &maxDof);
304: PetscSectionGetChart(section, &pStart, &pEnd);
305: PetscMalloc1(maxDof, &indices);
306: for (d = 0; d < maxDof; ++d) indices[d] = -1;
307: for (p = pStart; p < pEnd; ++p) {
308: PetscInt cdof, d;
310: PetscSectionGetConstraintDof(section, p, &cdof);
311: if (cdof) {
312: if (Nf) {
313: PetscInt numConst = 0, foff = 0;
315: for (f = 0; f < Nf; ++f) {
316: const PetscInt *find;
317: PetscInt fcdof, fdof;
319: PetscSectionGetFieldDof(section, p, f, &fdof);
320: PetscSectionGetFieldConstraintDof(section, p, f, &fcdof);
321: /* Change constraint numbering from field component to local dof number */
322: PetscSectionGetFieldConstraintIndices(section, p, f, &find);
323: for (d = 0; d < fcdof; ++d) indices[numConst+d] = find[d] + foff;
324: numConst += fcdof;
325: foff += fdof;
326: }
327: if (cdof != numConst) {PetscSectionSetConstraintDof(section, p, numConst);}
328: } else {
329: for (d = 0; d < cdof; ++d) indices[d] = d;
330: }
331: PetscSectionSetConstraintIndices(section, p, indices);
332: }
333: }
334: PetscFree(indices);
335: return(0);
336: }
338: /*@C
339: DMPlexCreateSection - Create a PetscSection based upon the dof layout specification provided.
341: Not Collective
343: Input Parameters:
344: + dm - The DMPlex object
345: . numComp - An array of size numFields that holds the number of components for each field
346: . numDof - An array of size numFields*(dim+1) which holds the number of dof for each field on a mesh piece of dimension d
347: . label - The label indicating the mesh support of each field, or NULL for the whole mesh
348: . numBC - The number of boundary conditions
349: . bcField - An array of size numBC giving the field number for each boundry condition
350: . bcComps - [Optional] An array of size numBC giving an IS holding the field components to which each boundary condition applies
351: . bcPoints - An array of size numBC giving an IS holding the Plex points to which each boundary condition applies
352: - perm - Optional permutation of the chart, or NULL
354: Output Parameter:
355: . section - The PetscSection object
357: Notes:
358: numDof[f*(dim+1)+d] gives the number of dof for field f on points of dimension d. For instance, numDof[1] is the
359: number of dof for field 0 on each edge.
361: The chart permutation is the same one set using PetscSectionSetPermutation()
363: Level: developer
365: Fortran Notes:
366: A Fortran 90 version is available as DMPlexCreateSectionF90()
368: .keywords: mesh, elements
369: .seealso: DMPlexCreate(), PetscSectionCreate(), PetscSectionSetPermutation()
370: @*/
371: PetscErrorCode DMPlexCreateSection(DM dm, DMLabel label[], const PetscInt numComp[],const PetscInt numDof[], PetscInt numBC, const PetscInt bcField[], const IS bcComps[], const IS bcPoints[], IS perm, PetscSection *section)
372: {
373: PetscSection aSec;
377: DMPlexCreateSectionFields(dm, numComp, section);
378: DMPlexCreateSectionDof(dm, label, numDof, *section);
379: DMPlexCreateSectionBCDof(dm, numBC, bcField, bcComps, bcPoints, *section);
380: if (perm) {PetscSectionSetPermutation(*section, perm);}
381: PetscSectionSetFromOptions(*section);
382: PetscSectionSetUp(*section);
383: DMPlexGetAnchors(dm,&aSec,NULL);
384: if (numBC || aSec) {
385: DMPlexCreateSectionBCIndicesField(dm, numBC, bcField, bcComps, bcPoints, *section);
386: DMPlexCreateSectionBCIndices(dm, *section);
387: }
388: PetscSectionViewFromOptions(*section,NULL,"-section_view");
389: return(0);
390: }
392: PetscErrorCode DMCreateDefaultSection_Plex(DM dm)
393: {
394: PetscDS probBC;
395: PetscSection section;
396: DMLabel *labels;
397: IS *bcPoints, *bcComps;
398: PetscBool *isFE;
399: PetscInt *bcFields, *numComp, *numDof;
400: PetscInt depth, dim, numBd, numBC = 0, Nf, bd, bc = 0, f;
401: PetscInt cStart, cEnd, cEndInterior;
405: DMGetNumFields(dm, &Nf);
406: /* FE and FV boundary conditions are handled slightly differently */
407: PetscMalloc1(Nf, &isFE);
408: for (f = 0; f < Nf; ++f) {
409: PetscObject obj;
410: PetscClassId id;
412: DMGetField(dm, f, NULL, &obj);
413: PetscObjectGetClassId(obj, &id);
414: if (id == PETSCFE_CLASSID) {isFE[f] = PETSC_TRUE;}
415: else if (id == PETSCFV_CLASSID) {isFE[f] = PETSC_FALSE;}
416: else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %D", f);
417: }
418: /* Allocate boundary point storage for FEM boundaries */
419: DMPlexGetDepth(dm, &depth);
420: DMGetDimension(dm, &dim);
421: DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
422: DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);
423: DMGetDS(dm, &probBC);
424: PetscDSGetNumBoundary(probBC, &numBd);
425: if (!Nf && numBd) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "number of fields is zero and number of boundary conditions is nonzero (this should never happen)");
426: for (bd = 0; bd < numBd; ++bd) {
427: PetscInt field;
428: DMBoundaryConditionType type;
429: const char *labelName;
430: DMLabel label;
432: PetscDSGetBoundary(probBC, bd, &type, NULL, &labelName, &field, NULL, NULL, NULL, NULL, NULL, NULL);
433: DMGetLabel(dm,labelName,&label);
434: if (label && isFE[field] && (type & DM_BC_ESSENTIAL)) ++numBC;
435: }
436: /* Add ghost cell boundaries for FVM */
437: for (f = 0; f < Nf; ++f) if (!isFE[f] && cEndInterior >= 0) ++numBC;
438: PetscCalloc3(numBC,&bcFields,numBC,&bcPoints,numBC,&bcComps);
439: /* Constrain ghost cells for FV */
440: for (f = 0; f < Nf; ++f) {
441: PetscInt *newidx, c;
443: if (isFE[f] || cEndInterior < 0) continue;
444: PetscMalloc1(cEnd-cEndInterior,&newidx);
445: for (c = cEndInterior; c < cEnd; ++c) newidx[c-cEndInterior] = c;
446: bcFields[bc] = f;
447: ISCreateGeneral(PetscObjectComm((PetscObject) dm), cEnd-cEndInterior, newidx, PETSC_OWN_POINTER, &bcPoints[bc++]);
448: }
449: /* Handle FEM Dirichlet boundaries */
450: for (bd = 0; bd < numBd; ++bd) {
451: const char *bdLabel;
452: DMLabel label;
453: const PetscInt *comps;
454: const PetscInt *values;
455: PetscInt bd2, field, numComps, numValues;
456: DMBoundaryConditionType type;
457: PetscBool duplicate = PETSC_FALSE;
459: PetscDSGetBoundary(probBC, bd, &type, NULL, &bdLabel, &field, &numComps, &comps, NULL, &numValues, &values, NULL);
460: DMGetLabel(dm, bdLabel, &label);
461: if (!isFE[field] || !label) continue;
462: /* Only want to modify label once */
463: for (bd2 = 0; bd2 < bd; ++bd2) {
464: const char *bdname;
465: PetscDSGetBoundary(probBC, bd2, NULL, NULL, &bdname, NULL, NULL, NULL, NULL, NULL, NULL, NULL);
466: PetscStrcmp(bdname, bdLabel, &duplicate);
467: if (duplicate) break;
468: }
469: if (!duplicate && (isFE[field])) {
470: /* don't complete cells, which are just present to give orientation to the boundary */
471: DMPlexLabelComplete(dm, label);
472: }
473: /* Filter out cells, if you actually want to constrain cells you need to do things by hand right now */
474: if (type & DM_BC_ESSENTIAL) {
475: PetscInt *newidx;
476: PetscInt n, newn = 0, p, v;
478: bcFields[bc] = field;
479: if (numComps) {ISCreateGeneral(PetscObjectComm((PetscObject) dm), numComps, comps, PETSC_COPY_VALUES, &bcComps[bc]);}
480: for (v = 0; v < numValues; ++v) {
481: IS tmp;
482: const PetscInt *idx;
484: DMGetStratumIS(dm, bdLabel, values[v], &tmp);
485: if (!tmp) continue;
486: ISGetLocalSize(tmp, &n);
487: ISGetIndices(tmp, &idx);
488: if (isFE[field]) {
489: for (p = 0; p < n; ++p) if ((idx[p] < cStart) || (idx[p] >= cEnd)) ++newn;
490: } else {
491: for (p = 0; p < n; ++p) if ((idx[p] >= cStart) || (idx[p] < cEnd)) ++newn;
492: }
493: ISRestoreIndices(tmp, &idx);
494: ISDestroy(&tmp);
495: }
496: PetscMalloc1(newn,&newidx);
497: newn = 0;
498: for (v = 0; v < numValues; ++v) {
499: IS tmp;
500: const PetscInt *idx;
502: DMGetStratumIS(dm, bdLabel, values[v], &tmp);
503: if (!tmp) continue;
504: ISGetLocalSize(tmp, &n);
505: ISGetIndices(tmp, &idx);
506: if (isFE[field]) {
507: for (p = 0; p < n; ++p) if ((idx[p] < cStart) || (idx[p] >= cEnd)) newidx[newn++] = idx[p];
508: } else {
509: for (p = 0; p < n; ++p) if ((idx[p] >= cStart) || (idx[p] < cEnd)) newidx[newn++] = idx[p];
510: }
511: ISRestoreIndices(tmp, &idx);
512: ISDestroy(&tmp);
513: }
514: ISCreateGeneral(PetscObjectComm((PetscObject) dm), newn, newidx, PETSC_OWN_POINTER, &bcPoints[bc++]);
515: }
516: }
517: /* Handle discretization */
518: PetscCalloc3(Nf,&labels,Nf,&numComp,Nf*(dim+1),&numDof);
519: for (f = 0; f < Nf; ++f) {
520: labels[f] = dm->fields[f].label;
521: if (isFE[f]) {
522: PetscFE fe = (PetscFE) dm->fields[f].disc;
523: const PetscInt *numFieldDof;
524: PetscInt fedim, d;
526: PetscFEGetNumComponents(fe, &numComp[f]);
527: PetscFEGetNumDof(fe, &numFieldDof);
528: PetscFEGetSpatialDimension(fe, &fedim);
529: for (d = 0; d < PetscMin(dim, fedim)+1; ++d) numDof[f*(dim+1)+d] = numFieldDof[d];
530: } else {
531: PetscFV fv = (PetscFV) dm->fields[f].disc;
533: PetscFVGetNumComponents(fv, &numComp[f]);
534: numDof[f*(dim+1)+dim] = numComp[f];
535: }
536: }
537: for (f = 0; f < Nf; ++f) {
538: PetscInt d;
539: for (d = 1; d < dim; ++d) {
540: if ((numDof[f*(dim+1)+d] > 0) && (depth < dim)) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated when unknowns are specified on edges or faces.");
541: }
542: }
543: DMPlexCreateSection(dm, labels, numComp, numDof, numBC, bcFields, bcComps, bcPoints, NULL, §ion);
544: for (f = 0; f < Nf; ++f) {
545: PetscFE fe;
546: const char *name;
548: if (isFE[f]) {
549: DMGetField(dm, f, NULL, (PetscObject *) &fe);
550: PetscObjectGetName((PetscObject) fe, &name);
551: PetscSectionSetFieldName(section, f, name);
552: }
553: }
554: DMSetSection(dm, section);
555: PetscSectionDestroy(§ion);
556: for (bc = 0; bc < numBC; ++bc) {ISDestroy(&bcPoints[bc]);ISDestroy(&bcComps[bc]);}
557: PetscFree3(bcFields,bcPoints,bcComps);
558: PetscFree3(labels,numComp,numDof);
559: PetscFree(isFE);
560: return(0);
561: }