Actual source code: plexsection.c

petsc-master 2019-08-22
Report Typos and Errors
  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:         /* Cases: dim > depth (cell-vertex mesh), dim == depth (fully interpolated), dim < depth (interpolated submesh) */
136:         d    = dim <= depth ? dep : (!dep ? 0 : dim);
137:         DMPlexGetDepthStratum(dm, dep, &pStart, &pEnd);
138:         pMax[dep] = pMax[dep] < 0 ? pEnd : pMax[dep];
139:         for (p = pStart; p < pEnd; ++p) {
140:           const PetscInt dof = numDof[f*(dim+1)+d];

142:           if (isFE[f] && p >= pMax[dep]) continue;
143:           PetscSectionSetFieldDof(section, p, f, dof);
144:           PetscSectionAddDof(section, p, dof);
145:         }
146:       }
147:     }
148:   }
149:   PetscFree(pMax);
150:   PetscFree(isFE);
151:   return(0);
152: }

154: /* Set the number of dof on each point and separate by fields
155:    If bcComps is NULL or the IS is NULL, constrain every dof on the point
156: */
157: static PetscErrorCode DMPlexCreateSectionBCDof(DM dm, PetscInt numBC, const PetscInt bcField[], const IS bcComps[], const IS bcPoints[], PetscSection section)
158: {
159:   PetscInt       Nf;
160:   PetscInt       bc;
161:   PetscSection   aSec;

165:   PetscSectionGetNumFields(section, &Nf);
166:   for (bc = 0; bc < numBC; ++bc) {
167:     PetscInt        field = 0;
168:     const PetscInt *comp;
169:     const PetscInt *idx;
170:     PetscInt        Nc = -1, n, i;

172:     if (Nf) field = bcField[bc];
173:     if (bcComps && bcComps[bc]) {ISGetLocalSize(bcComps[bc], &Nc);}
174:     if (bcComps && bcComps[bc]) {ISGetIndices(bcComps[bc], &comp);}
175:     ISGetLocalSize(bcPoints[bc], &n);
176:     ISGetIndices(bcPoints[bc], &idx);
177:     for (i = 0; i < n; ++i) {
178:       const PetscInt p = idx[i];
179:       PetscInt       numConst;

181:       if (Nf) {
182:         PetscSectionGetFieldDof(section, p, field, &numConst);
183:       } else {
184:         PetscSectionGetDof(section, p, &numConst);
185:       }
186:       /* If Nc < 0, constrain every dof on the point */
187:       /* TODO: Matt, this only works if there is one node on the point.  We need to handle numDofs > NumComponents */
188:       if (Nc > 0) numConst = PetscMin(numConst, Nc);
189:       if (Nf) {PetscSectionAddFieldConstraintDof(section, p, field, numConst);}
190:       PetscSectionAddConstraintDof(section, p, numConst);
191:     }
192:     ISRestoreIndices(bcPoints[bc], &idx);
193:     if (bcComps && bcComps[bc]) {ISRestoreIndices(bcComps[bc], &comp);}
194:   }
195:   DMPlexGetAnchors(dm, &aSec, NULL);
196:   if (aSec) {
197:     PetscInt aStart, aEnd, a;

199:     PetscSectionGetChart(aSec, &aStart, &aEnd);
200:     for (a = aStart; a < aEnd; a++) {
201:       PetscInt dof, f;

203:       PetscSectionGetDof(aSec, a, &dof);
204:       if (dof) {
205:         /* if there are point-to-point constraints, then all dofs are constrained */
206:         PetscSectionGetDof(section, a, &dof);
207:         PetscSectionSetConstraintDof(section, a, dof);
208:         for (f = 0; f < Nf; f++) {
209:           PetscSectionGetFieldDof(section, a, f, &dof);
210:           PetscSectionSetFieldConstraintDof(section, a, f, dof);
211:         }
212:       }
213:     }
214:   }
215:   return(0);
216: }

218: /* Set the constrained field indices on each point
219:    If bcComps is NULL or the IS is NULL, constrain every dof on the point
220: */
221: static PetscErrorCode DMPlexCreateSectionBCIndicesField(DM dm, PetscInt numBC,const PetscInt bcField[], const IS bcComps[], const IS bcPoints[], PetscSection section)
222: {
223:   PetscSection   aSec;
224:   PetscInt      *indices;
225:   PetscInt       Nf, cdof, maxDof = 0, pStart, pEnd, p, bc, f, d;

229:   PetscSectionGetNumFields(section, &Nf);
230:   if (!Nf) return(0);
231:   /* Initialize all field indices to -1 */
232:   PetscSectionGetChart(section, &pStart, &pEnd);
233:   for (p = pStart; p < pEnd; ++p) {PetscSectionGetConstraintDof(section, p, &cdof); maxDof = PetscMax(maxDof, cdof);}
234:   PetscMalloc1(maxDof, &indices);
235:   for (d = 0; d < maxDof; ++d) indices[d] = -1;
236:   for (p = pStart; p < pEnd; ++p) for (f = 0; f < Nf; ++f) {PetscSectionSetFieldConstraintIndices(section, p, f, indices);}
237:   /* Handle BC constraints */
238:   for (bc = 0; bc < numBC; ++bc) {
239:     const PetscInt  field = bcField[bc];
240:     const PetscInt *comp, *idx;
241:     PetscInt        Nc = -1, n, i;

243:     if (bcComps && bcComps[bc]) {ISGetLocalSize(bcComps[bc], &Nc);}
244:     if (bcComps && bcComps[bc]) {ISGetIndices(bcComps[bc], &comp);}
245:     ISGetLocalSize(bcPoints[bc], &n);
246:     ISGetIndices(bcPoints[bc], &idx);
247:     for (i = 0; i < n; ++i) {
248:       const PetscInt  p = idx[i];
249:       const PetscInt *find;
250:       PetscInt        fdof, fcdof, c;

252:       PetscSectionGetFieldDof(section, p, field, &fdof);
253:       if (!fdof) continue;
254:       if (Nc < 0) {
255:         for (d = 0; d < fdof; ++d) indices[d] = d;
256:         fcdof = fdof;
257:       } else {
258:         PetscSectionGetFieldConstraintDof(section, p, field, &fcdof);
259:         PetscSectionGetFieldConstraintIndices(section, p, field, &find);
260:         for (d = 0; d < fcdof; ++d) {if (find[d] < 0) break; indices[d] = find[d];}
261:         for (c = 0; c < Nc; ++c) indices[d++] = comp[c];
262:         PetscSortRemoveDupsInt(&d, indices);
263:         for (c = d; c < fcdof; ++c) indices[c] = -1;
264:         fcdof = d;
265:       }
266:       PetscSectionSetFieldConstraintDof(section, p, field, fcdof);
267:       PetscSectionSetFieldConstraintIndices(section, p, field, indices);
268:     }
269:     if (bcComps && bcComps[bc]) {ISRestoreIndices(bcComps[bc], &comp);}
270:     ISRestoreIndices(bcPoints[bc], &idx);
271:   }
272:   /* Handle anchors */
273:   DMPlexGetAnchors(dm, &aSec, NULL);
274:   if (aSec) {
275:     PetscInt aStart, aEnd, a;

277:     for (d = 0; d < maxDof; ++d) indices[d] = d;
278:     PetscSectionGetChart(aSec, &aStart, &aEnd);
279:     for (a = aStart; a < aEnd; a++) {
280:       PetscInt dof, f;

282:       PetscSectionGetDof(aSec, a, &dof);
283:       if (dof) {
284:         /* if there are point-to-point constraints, then all dofs are constrained */
285:         for (f = 0; f < Nf; f++) {
286:           PetscSectionSetFieldConstraintIndices(section, a, f, indices);
287:         }
288:       }
289:     }
290:   }
291:   PetscFree(indices);
292:   return(0);
293: }

295: /* Set the constrained indices on each point */
296: static PetscErrorCode DMPlexCreateSectionBCIndices(DM dm, PetscSection section)
297: {
298:   PetscInt      *indices;
299:   PetscInt       Nf, maxDof, pStart, pEnd, p, f, d;

303:   PetscSectionGetNumFields(section, &Nf);
304:   PetscSectionGetMaxDof(section, &maxDof);
305:   PetscSectionGetChart(section, &pStart, &pEnd);
306:   PetscMalloc1(maxDof, &indices);
307:   for (d = 0; d < maxDof; ++d) indices[d] = -1;
308:   for (p = pStart; p < pEnd; ++p) {
309:     PetscInt cdof, d;

311:     PetscSectionGetConstraintDof(section, p, &cdof);
312:     if (cdof) {
313:       if (Nf) {
314:         PetscInt numConst = 0, foff = 0;

316:         for (f = 0; f < Nf; ++f) {
317:           const PetscInt *find;
318:           PetscInt        fcdof, fdof;

320:           PetscSectionGetFieldDof(section, p, f, &fdof);
321:           PetscSectionGetFieldConstraintDof(section, p, f, &fcdof);
322:           /* Change constraint numbering from field component to local dof number */
323:           PetscSectionGetFieldConstraintIndices(section, p, f, &find);
324:           for (d = 0; d < fcdof; ++d) indices[numConst+d] = find[d] + foff;
325:           numConst += fcdof;
326:           foff     += fdof;
327:         }
328:         if (cdof != numConst) {PetscSectionSetConstraintDof(section, p, numConst);}
329:       } else {
330:         for (d = 0; d < cdof; ++d) indices[d] = d;
331:       }
332:       PetscSectionSetConstraintIndices(section, p, indices);
333:     }
334:   }
335:   PetscFree(indices);
336:   return(0);
337: }

339: /*@C
340:   DMPlexCreateSection - Create a PetscSection based upon the dof layout specification provided.

342:   Not Collective

344:   Input Parameters:
345: + dm        - The DMPlex object
346: . numComp   - An array of size numFields that holds the number of components for each field
347: . numDof    - An array of size numFields*(dim+1) which holds the number of dof for each field on a mesh piece of dimension d
348: . label     - The label indicating the mesh support of each field, or NULL for the whole mesh
349: . numBC     - The number of boundary conditions
350: . bcField   - An array of size numBC giving the field number for each boundry condition
351: . bcComps   - [Optional] An array of size numBC giving an IS holding the field components to which each boundary condition applies
352: . bcPoints  - An array of size numBC giving an IS holding the Plex points to which each boundary condition applies
353: - perm      - Optional permutation of the chart, or NULL

355:   Output Parameter:
356: . section - The PetscSection object

358:   Notes:
359:     numDof[f*(dim+1)+d] gives the number of dof for field f on points of dimension d. For instance, numDof[1] is the
360:   number of dof for field 0 on each edge.

362:   The chart permutation is the same one set using PetscSectionSetPermutation()

364:   Level: developer

366:   Fortran Notes:
367:   A Fortran 90 version is available as DMPlexCreateSectionF90()

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, &section);
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:   DMSetLocalSection(dm, section);
555:   PetscSectionDestroy(&section);
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: }