1 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 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; 9 10 PetscFunctionBegin; 11 PetscCall(DMPlexGetDepth(dm, &depth)); 12 PetscCall(DMPlexGetDepthLabel(dm,&depthLabel)); 13 PetscCall(DMGetNumFields(dm, &Nf)); 14 PetscCall(PetscCalloc1(Nf, &isFE)); 15 for (f = 0; f < Nf; ++f) { 16 PetscObject obj; 17 PetscClassId id; 18 19 PetscCall(DMGetField(dm, f, NULL, &obj)); 20 PetscCall(PetscObjectGetClassId(obj, &id)); 21 if (id == PETSCFE_CLASSID) {isFE[f] = PETSC_TRUE;} 22 else if (id == PETSCFV_CLASSID) {isFE[f] = PETSC_FALSE;} 23 } 24 25 PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), section)); 26 if (Nf > 0) { 27 PetscCall(PetscSectionSetNumFields(*section, Nf)); 28 if (numComp) { 29 for (f = 0; f < Nf; ++f) { 30 PetscCall(PetscSectionSetFieldComponents(*section, f, numComp[f])); 31 if (isFE[f]) { 32 PetscFE fe; 33 PetscDualSpace dspace; 34 const PetscInt ***perms; 35 const PetscScalar ***flips; 36 const PetscInt *numDof; 37 38 PetscCall(DMGetField(dm,f,NULL,(PetscObject *) &fe)); 39 PetscCall(PetscFEGetDualSpace(fe,&dspace)); 40 PetscCall(PetscDualSpaceGetSymmetries(dspace,&perms,&flips)); 41 PetscCall(PetscDualSpaceGetNumDof(dspace,&numDof)); 42 if (perms || flips) { 43 DM K; 44 PetscInt sph, spdepth; 45 PetscSectionSym sym; 46 47 PetscCall(PetscDualSpaceGetDM(dspace,&K)); 48 PetscCall(DMPlexGetDepth(K, &spdepth)); 49 PetscCall(PetscSectionSymCreateLabel(PetscObjectComm((PetscObject)*section),depthLabel,&sym)); 50 for (sph = 0; sph <= spdepth; sph++) { 51 PetscDualSpace hspace; 52 PetscInt kStart, kEnd; 53 PetscInt kConeSize, h = sph + (depth - spdepth); 54 const PetscInt **perms0 = NULL; 55 const PetscScalar **flips0 = NULL; 56 57 PetscCall(PetscDualSpaceGetHeightSubspace(dspace, sph, &hspace)); 58 PetscCall(DMPlexGetHeightStratum(K, h, &kStart, &kEnd)); 59 if (!hspace) continue; 60 PetscCall(PetscDualSpaceGetSymmetries(hspace,&perms,&flips)); 61 if (perms) perms0 = perms[0]; 62 if (flips) flips0 = flips[0]; 63 if (!(perms0 || flips0)) continue; 64 { 65 DMPolytopeType ct; 66 /* The number of arrangements is no longer based on the number of faces */ 67 PetscCall(DMPlexGetCellType(K, kStart, &ct)); 68 kConeSize = DMPolytopeTypeGetNumArrangments(ct) / 2; 69 } 70 PetscCall(PetscSectionSymLabelSetStratum(sym,depth - h,numDof[depth - h],-kConeSize,kConeSize,PETSC_USE_POINTER,perms0 ? &perms0[-kConeSize] : NULL,flips0 ? &flips0[-kConeSize] : NULL)); 71 } 72 PetscCall(PetscSectionSetFieldSym(*section,f,sym)); 73 PetscCall(PetscSectionSymDestroy(&sym)); 74 } 75 } 76 } 77 } 78 } 79 PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 80 PetscCall(PetscSectionSetChart(*section, pStart, pEnd)); 81 PetscCall(PetscFree(isFE)); 82 PetscFunctionReturn(0); 83 } 84 85 /* Set the number of dof on each point and separate by fields */ 86 static PetscErrorCode DMPlexCreateSectionDof(DM dm, DMLabel label[],const PetscInt numDof[], PetscSection section) 87 { 88 DMLabel depthLabel; 89 DMPolytopeType ct; 90 PetscInt depth, cellHeight, pStart = 0, pEnd = 0; 91 PetscInt Nf, f, Nds, n, dim, d, dep, p; 92 PetscBool *isFE, hasCohesive = PETSC_FALSE; 93 94 PetscFunctionBegin; 95 PetscCall(DMGetDimension(dm, &dim)); 96 PetscCall(DMPlexGetDepth(dm, &depth)); 97 PetscCall(DMPlexGetDepthLabel(dm,&depthLabel)); 98 PetscCall(DMGetNumFields(dm, &Nf)); 99 PetscCall(DMGetNumDS(dm, &Nds)); 100 for (n = 0; n < Nds; ++n) { 101 PetscDS ds; 102 PetscBool isCohesive; 103 104 PetscCall(DMGetRegionNumDS(dm, n, NULL, NULL, &ds)); 105 PetscCall(PetscDSIsCohesive(ds, &isCohesive)); 106 if (isCohesive) {hasCohesive = PETSC_TRUE; break;} 107 } 108 PetscCall(PetscMalloc1(Nf, &isFE)); 109 for (f = 0; f < Nf; ++f) { 110 PetscObject obj; 111 PetscClassId id; 112 113 PetscCall(DMGetField(dm, f, NULL, &obj)); 114 PetscCall(PetscObjectGetClassId(obj, &id)); 115 /* User is allowed to put a "placeholder" field in (c.f. DMCreateDS) */ 116 isFE[f] = id == PETSCFE_CLASSID ? PETSC_TRUE : PETSC_FALSE; 117 } 118 119 PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight)); 120 for (f = 0; f < Nf; ++f) { 121 PetscBool avoidTensor; 122 123 PetscCall(DMGetFieldAvoidTensor(dm, f, &avoidTensor)); 124 avoidTensor = (avoidTensor || hasCohesive) ? PETSC_TRUE : PETSC_FALSE; 125 if (label && label[f]) { 126 IS pointIS; 127 const PetscInt *points; 128 PetscInt n; 129 130 PetscCall(DMLabelGetStratumIS(label[f], 1, &pointIS)); 131 if (!pointIS) continue; 132 PetscCall(ISGetLocalSize(pointIS, &n)); 133 PetscCall(ISGetIndices(pointIS, &points)); 134 for (p = 0; p < n; ++p) { 135 const PetscInt point = points[p]; 136 PetscInt dof, d; 137 138 PetscCall(DMPlexGetCellType(dm, point, &ct)); 139 PetscCall(DMLabelGetValue(depthLabel, point, &d)); 140 /* If this is a tensor prism point, use dof for one dimension lower */ 141 switch (ct) { 142 case DM_POLYTOPE_POINT_PRISM_TENSOR: 143 case DM_POLYTOPE_SEG_PRISM_TENSOR: 144 case DM_POLYTOPE_TRI_PRISM_TENSOR: 145 case DM_POLYTOPE_QUAD_PRISM_TENSOR: 146 if (hasCohesive) {--d;} break; 147 default: break; 148 } 149 dof = d < 0 ? 0 : numDof[f*(dim+1)+d]; 150 PetscCall(PetscSectionSetFieldDof(section, point, f, dof)); 151 PetscCall(PetscSectionAddDof(section, point, dof)); 152 } 153 PetscCall(ISRestoreIndices(pointIS, &points)); 154 PetscCall(ISDestroy(&pointIS)); 155 } else { 156 for (dep = 0; dep <= depth - cellHeight; ++dep) { 157 /* Cases: dim > depth (cell-vertex mesh), dim == depth (fully interpolated), dim < depth (interpolated submesh) */ 158 d = dim <= depth ? dep : (!dep ? 0 : dim); 159 PetscCall(DMPlexGetDepthStratum(dm, dep, &pStart, &pEnd)); 160 for (p = pStart; p < pEnd; ++p) { 161 const PetscInt dof = numDof[f*(dim+1)+d]; 162 163 PetscCall(DMPlexGetCellType(dm, p, &ct)); 164 switch (ct) { 165 case DM_POLYTOPE_POINT_PRISM_TENSOR: 166 case DM_POLYTOPE_SEG_PRISM_TENSOR: 167 case DM_POLYTOPE_TRI_PRISM_TENSOR: 168 case DM_POLYTOPE_QUAD_PRISM_TENSOR: 169 if (avoidTensor && isFE[f]) continue; 170 default: break; 171 } 172 PetscCall(PetscSectionSetFieldDof(section, p, f, dof)); 173 PetscCall(PetscSectionAddDof(section, p, dof)); 174 } 175 } 176 } 177 } 178 PetscCall(PetscFree(isFE)); 179 PetscFunctionReturn(0); 180 } 181 182 /* Set the number of dof on each point and separate by fields 183 If bcComps is NULL or the IS is NULL, constrain every dof on the point 184 */ 185 static PetscErrorCode DMPlexCreateSectionBCDof(DM dm, PetscInt numBC, const PetscInt bcField[], const IS bcComps[], const IS bcPoints[], PetscSection section) 186 { 187 PetscInt Nf; 188 PetscInt bc; 189 PetscSection aSec; 190 191 PetscFunctionBegin; 192 PetscCall(PetscSectionGetNumFields(section, &Nf)); 193 for (bc = 0; bc < numBC; ++bc) { 194 PetscInt field = 0; 195 const PetscInt *comp; 196 const PetscInt *idx; 197 PetscInt Nc = 0, cNc = -1, n, i; 198 199 if (Nf) { 200 field = bcField[bc]; 201 PetscCall(PetscSectionGetFieldComponents(section, field, &Nc)); 202 } 203 if (bcComps && bcComps[bc]) PetscCall(ISGetLocalSize(bcComps[bc], &cNc)); 204 if (bcComps && bcComps[bc]) PetscCall(ISGetIndices(bcComps[bc], &comp)); 205 PetscCall(ISGetLocalSize(bcPoints[bc], &n)); 206 PetscCall(ISGetIndices(bcPoints[bc], &idx)); 207 for (i = 0; i < n; ++i) { 208 const PetscInt p = idx[i]; 209 PetscInt numConst; 210 211 if (Nf) { 212 PetscCall(PetscSectionGetFieldDof(section, p, field, &numConst)); 213 } else { 214 PetscCall(PetscSectionGetDof(section, p, &numConst)); 215 } 216 /* If Nc <= 0, constrain every dof on the point */ 217 if (cNc > 0) { 218 /* We assume that a point may have multiple "nodes", which are collections of Nc dofs, 219 and that those dofs are numbered n*Nc+c */ 220 if (Nf) { 221 PetscCheckFalse(numConst % Nc,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Point %D has %D dof which is not divisible by %D field components", p, numConst, Nc); 222 numConst = (numConst/Nc) * cNc; 223 } else { 224 numConst = PetscMin(numConst, cNc); 225 } 226 } 227 if (Nf) PetscCall(PetscSectionAddFieldConstraintDof(section, p, field, numConst)); 228 PetscCall(PetscSectionAddConstraintDof(section, p, numConst)); 229 } 230 PetscCall(ISRestoreIndices(bcPoints[bc], &idx)); 231 if (bcComps && bcComps[bc]) PetscCall(ISRestoreIndices(bcComps[bc], &comp)); 232 } 233 PetscCall(DMPlexGetAnchors(dm, &aSec, NULL)); 234 if (aSec) { 235 PetscInt aStart, aEnd, a; 236 237 PetscCall(PetscSectionGetChart(aSec, &aStart, &aEnd)); 238 for (a = aStart; a < aEnd; a++) { 239 PetscInt dof, f; 240 241 PetscCall(PetscSectionGetDof(aSec, a, &dof)); 242 if (dof) { 243 /* if there are point-to-point constraints, then all dofs are constrained */ 244 PetscCall(PetscSectionGetDof(section, a, &dof)); 245 PetscCall(PetscSectionSetConstraintDof(section, a, dof)); 246 for (f = 0; f < Nf; f++) { 247 PetscCall(PetscSectionGetFieldDof(section, a, f, &dof)); 248 PetscCall(PetscSectionSetFieldConstraintDof(section, a, f, dof)); 249 } 250 } 251 } 252 } 253 PetscFunctionReturn(0); 254 } 255 256 /* Set the constrained field indices on each point 257 If bcComps is NULL or the IS is NULL, constrain every dof on the point 258 */ 259 static PetscErrorCode DMPlexCreateSectionBCIndicesField(DM dm, PetscInt numBC,const PetscInt bcField[], const IS bcComps[], const IS bcPoints[], PetscSection section) 260 { 261 PetscSection aSec; 262 PetscInt *indices; 263 PetscInt Nf, cdof, maxDof = 0, pStart, pEnd, p, bc, f, d; 264 265 PetscFunctionBegin; 266 PetscCall(PetscSectionGetNumFields(section, &Nf)); 267 if (!Nf) PetscFunctionReturn(0); 268 /* Initialize all field indices to -1 */ 269 PetscCall(PetscSectionGetChart(section, &pStart, &pEnd)); 270 for (p = pStart; p < pEnd; ++p) {PetscCall(PetscSectionGetConstraintDof(section, p, &cdof)); maxDof = PetscMax(maxDof, cdof);} 271 PetscCall(PetscMalloc1(maxDof, &indices)); 272 for (d = 0; d < maxDof; ++d) indices[d] = -1; 273 for (p = pStart; p < pEnd; ++p) for (f = 0; f < Nf; ++f) PetscCall(PetscSectionSetFieldConstraintIndices(section, p, f, indices)); 274 /* Handle BC constraints */ 275 for (bc = 0; bc < numBC; ++bc) { 276 const PetscInt field = bcField[bc]; 277 const PetscInt *comp, *idx; 278 PetscInt Nc, cNc = -1, n, i; 279 280 PetscCall(PetscSectionGetFieldComponents(section, field, &Nc)); 281 if (bcComps && bcComps[bc]) PetscCall(ISGetLocalSize(bcComps[bc], &cNc)); 282 if (bcComps && bcComps[bc]) PetscCall(ISGetIndices(bcComps[bc], &comp)); 283 PetscCall(ISGetLocalSize(bcPoints[bc], &n)); 284 PetscCall(ISGetIndices(bcPoints[bc], &idx)); 285 for (i = 0; i < n; ++i) { 286 const PetscInt p = idx[i]; 287 const PetscInt *find; 288 PetscInt fdof, fcdof, c, j; 289 290 PetscCall(PetscSectionGetFieldDof(section, p, field, &fdof)); 291 if (!fdof) continue; 292 if (cNc < 0) { 293 for (d = 0; d < fdof; ++d) indices[d] = d; 294 fcdof = fdof; 295 } else { 296 /* We assume that a point may have multiple "nodes", which are collections of Nc dofs, 297 and that those dofs are numbered n*Nc+c */ 298 PetscCall(PetscSectionGetFieldConstraintDof(section, p, field, &fcdof)); 299 PetscCall(PetscSectionGetFieldConstraintIndices(section, p, field, &find)); 300 /* Get indices constrained by previous bcs */ 301 for (d = 0; d < fcdof; ++d) {if (find[d] < 0) break; indices[d] = find[d];} 302 for (j = 0; j < fdof/Nc; ++j) for (c = 0; c < cNc; ++c) indices[d++] = j*Nc + comp[c]; 303 PetscCall(PetscSortRemoveDupsInt(&d, indices)); 304 for (c = d; c < fcdof; ++c) indices[c] = -1; 305 fcdof = d; 306 } 307 PetscCall(PetscSectionSetFieldConstraintDof(section, p, field, fcdof)); 308 PetscCall(PetscSectionSetFieldConstraintIndices(section, p, field, indices)); 309 } 310 if (bcComps && bcComps[bc]) PetscCall(ISRestoreIndices(bcComps[bc], &comp)); 311 PetscCall(ISRestoreIndices(bcPoints[bc], &idx)); 312 } 313 /* Handle anchors */ 314 PetscCall(DMPlexGetAnchors(dm, &aSec, NULL)); 315 if (aSec) { 316 PetscInt aStart, aEnd, a; 317 318 for (d = 0; d < maxDof; ++d) indices[d] = d; 319 PetscCall(PetscSectionGetChart(aSec, &aStart, &aEnd)); 320 for (a = aStart; a < aEnd; a++) { 321 PetscInt dof, f; 322 323 PetscCall(PetscSectionGetDof(aSec, a, &dof)); 324 if (dof) { 325 /* if there are point-to-point constraints, then all dofs are constrained */ 326 for (f = 0; f < Nf; f++) { 327 PetscCall(PetscSectionSetFieldConstraintIndices(section, a, f, indices)); 328 } 329 } 330 } 331 } 332 PetscCall(PetscFree(indices)); 333 PetscFunctionReturn(0); 334 } 335 336 /* Set the constrained indices on each point */ 337 static PetscErrorCode DMPlexCreateSectionBCIndices(DM dm, PetscSection section) 338 { 339 PetscInt *indices; 340 PetscInt Nf, maxDof, pStart, pEnd, p, f, d; 341 342 PetscFunctionBegin; 343 PetscCall(PetscSectionGetNumFields(section, &Nf)); 344 PetscCall(PetscSectionGetMaxDof(section, &maxDof)); 345 PetscCall(PetscSectionGetChart(section, &pStart, &pEnd)); 346 PetscCall(PetscMalloc1(maxDof, &indices)); 347 for (d = 0; d < maxDof; ++d) indices[d] = -1; 348 for (p = pStart; p < pEnd; ++p) { 349 PetscInt cdof, d; 350 351 PetscCall(PetscSectionGetConstraintDof(section, p, &cdof)); 352 if (cdof) { 353 if (Nf) { 354 PetscInt numConst = 0, foff = 0; 355 356 for (f = 0; f < Nf; ++f) { 357 const PetscInt *find; 358 PetscInt fcdof, fdof; 359 360 PetscCall(PetscSectionGetFieldDof(section, p, f, &fdof)); 361 PetscCall(PetscSectionGetFieldConstraintDof(section, p, f, &fcdof)); 362 /* Change constraint numbering from field component to local dof number */ 363 PetscCall(PetscSectionGetFieldConstraintIndices(section, p, f, &find)); 364 for (d = 0; d < fcdof; ++d) indices[numConst+d] = find[d] + foff; 365 numConst += fcdof; 366 foff += fdof; 367 } 368 if (cdof != numConst) PetscCall(PetscSectionSetConstraintDof(section, p, numConst)); 369 } else { 370 for (d = 0; d < cdof; ++d) indices[d] = d; 371 } 372 PetscCall(PetscSectionSetConstraintIndices(section, p, indices)); 373 } 374 } 375 PetscCall(PetscFree(indices)); 376 PetscFunctionReturn(0); 377 } 378 379 /*@C 380 DMPlexCreateSection - Create a PetscSection based upon the dof layout specification provided. 381 382 Not Collective 383 384 Input Parameters: 385 + dm - The DMPlex object 386 . label - The label indicating the mesh support of each field, or NULL for the whole mesh 387 . numComp - An array of size numFields that holds the number of components for each field 388 . numDof - An array of size numFields*(dim+1) which holds the number of dof for each field on a mesh piece of dimension d 389 . numBC - The number of boundary conditions 390 . bcField - An array of size numBC giving the field number for each boundry condition 391 . bcComps - [Optional] An array of size numBC giving an IS holding the field components to which each boundary condition applies 392 . bcPoints - An array of size numBC giving an IS holding the Plex points to which each boundary condition applies 393 - perm - Optional permutation of the chart, or NULL 394 395 Output Parameter: 396 . section - The PetscSection object 397 398 Notes: 399 numDof[f*(dim+1)+d] gives the number of dof for field f on points of dimension d. For instance, numDof[1] is the 400 number of dof for field 0 on each edge. 401 402 The chart permutation is the same one set using PetscSectionSetPermutation() 403 404 Level: developer 405 406 TODO: How is this related to DMCreateLocalSection() 407 408 .seealso: DMPlexCreate(), PetscSectionCreate(), PetscSectionSetPermutation() 409 @*/ 410 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) 411 { 412 PetscSection aSec; 413 414 PetscFunctionBegin; 415 PetscCall(DMPlexCreateSectionFields(dm, numComp, section)); 416 PetscCall(DMPlexCreateSectionDof(dm, label, numDof, *section)); 417 PetscCall(DMPlexCreateSectionBCDof(dm, numBC, bcField, bcComps, bcPoints, *section)); 418 if (perm) PetscCall(PetscSectionSetPermutation(*section, perm)); 419 PetscCall(PetscSectionSetFromOptions(*section)); 420 PetscCall(PetscSectionSetUp(*section)); 421 PetscCall(DMPlexGetAnchors(dm,&aSec,NULL)); 422 if (numBC || aSec) { 423 PetscCall(DMPlexCreateSectionBCIndicesField(dm, numBC, bcField, bcComps, bcPoints, *section)); 424 PetscCall(DMPlexCreateSectionBCIndices(dm, *section)); 425 } 426 PetscCall(PetscSectionViewFromOptions(*section,NULL,"-section_view")); 427 PetscFunctionReturn(0); 428 } 429 430 PetscErrorCode DMCreateLocalSection_Plex(DM dm) 431 { 432 PetscSection section; 433 DMLabel *labels; 434 IS *bcPoints, *bcComps; 435 PetscBool *isFE; 436 PetscInt *bcFields, *numComp, *numDof; 437 PetscInt depth, dim, numBC = 0, Nf, Nds, s, bc = 0, f; 438 PetscInt cStart, cEnd, cEndInterior; 439 440 PetscFunctionBegin; 441 PetscCall(DMGetNumFields(dm, &Nf)); 442 PetscCall(DMGetDimension(dm, &dim)); 443 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 444 /* FE and FV boundary conditions are handled slightly differently */ 445 PetscCall(PetscMalloc1(Nf, &isFE)); 446 for (f = 0; f < Nf; ++f) { 447 PetscObject obj; 448 PetscClassId id; 449 450 PetscCall(DMGetField(dm, f, NULL, &obj)); 451 PetscCall(PetscObjectGetClassId(obj, &id)); 452 if (id == PETSCFE_CLASSID) {isFE[f] = PETSC_TRUE;} 453 else if (id == PETSCFV_CLASSID) {isFE[f] = PETSC_FALSE;} 454 else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %D", f); 455 } 456 /* Allocate boundary point storage for FEM boundaries */ 457 PetscCall(DMGetNumDS(dm, &Nds)); 458 for (s = 0; s < Nds; ++s) { 459 PetscDS dsBC; 460 PetscInt numBd, bd; 461 462 PetscCall(DMGetRegionNumDS(dm, s, NULL, NULL, &dsBC)); 463 PetscCall(PetscDSGetNumBoundary(dsBC, &numBd)); 464 PetscCheckFalse(!Nf && numBd,PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "number of fields is zero and number of boundary conditions is nonzero (this should never happen)"); 465 for (bd = 0; bd < numBd; ++bd) { 466 PetscInt field; 467 DMBoundaryConditionType type; 468 DMLabel label; 469 470 PetscCall(PetscDSGetBoundary(dsBC, bd, NULL, &type, NULL, &label, NULL, NULL, &field, NULL, NULL, NULL, NULL, NULL)); 471 if (label && isFE[field] && (type & DM_BC_ESSENTIAL)) ++numBC; 472 } 473 } 474 /* Add ghost cell boundaries for FVM */ 475 PetscCall(DMPlexGetGhostCellStratum(dm, &cEndInterior, NULL)); 476 for (f = 0; f < Nf; ++f) if (!isFE[f] && cEndInterior >= 0) ++numBC; 477 PetscCall(PetscCalloc3(numBC, &bcFields, numBC, &bcPoints, numBC, &bcComps)); 478 /* Constrain ghost cells for FV */ 479 for (f = 0; f < Nf; ++f) { 480 PetscInt *newidx, c; 481 482 if (isFE[f] || cEndInterior < 0) continue; 483 PetscCall(PetscMalloc1(cEnd-cEndInterior,&newidx)); 484 for (c = cEndInterior; c < cEnd; ++c) newidx[c-cEndInterior] = c; 485 bcFields[bc] = f; 486 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, cEnd-cEndInterior, newidx, PETSC_OWN_POINTER, &bcPoints[bc++])); 487 } 488 /* Handle FEM Dirichlet boundaries */ 489 PetscCall(DMGetNumDS(dm, &Nds)); 490 for (s = 0; s < Nds; ++s) { 491 PetscDS dsBC; 492 PetscInt numBd, bd; 493 494 PetscCall(DMGetRegionNumDS(dm, s, NULL, NULL, &dsBC)); 495 PetscCall(PetscDSGetNumBoundary(dsBC, &numBd)); 496 for (bd = 0; bd < numBd; ++bd) { 497 DMLabel label; 498 const PetscInt *comps; 499 const PetscInt *values; 500 PetscInt bd2, field, numComps, numValues; 501 DMBoundaryConditionType type; 502 PetscBool duplicate = PETSC_FALSE; 503 504 PetscCall(PetscDSGetBoundary(dsBC, bd, NULL, &type, NULL, &label, &numValues, &values, &field, &numComps, &comps, NULL, NULL, NULL)); 505 if (!isFE[field] || !label) continue; 506 /* Only want to modify label once */ 507 for (bd2 = 0; bd2 < bd; ++bd2) { 508 DMLabel l; 509 510 PetscCall(PetscDSGetBoundary(dsBC, bd2, NULL, NULL, NULL, &l, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 511 duplicate = l == label ? PETSC_TRUE : PETSC_FALSE; 512 if (duplicate) break; 513 } 514 if (!duplicate && (isFE[field])) { 515 /* don't complete cells, which are just present to give orientation to the boundary */ 516 PetscCall(DMPlexLabelComplete(dm, label)); 517 } 518 /* Filter out cells, if you actually want to constrain cells you need to do things by hand right now */ 519 if (type & DM_BC_ESSENTIAL) { 520 PetscInt *newidx; 521 PetscInt n, newn = 0, p, v; 522 523 bcFields[bc] = field; 524 if (numComps) PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject) dm), numComps, comps, PETSC_COPY_VALUES, &bcComps[bc])); 525 for (v = 0; v < numValues; ++v) { 526 IS tmp; 527 const PetscInt *idx; 528 529 PetscCall(DMLabelGetStratumIS(label, values[v], &tmp)); 530 if (!tmp) continue; 531 PetscCall(ISGetLocalSize(tmp, &n)); 532 PetscCall(ISGetIndices(tmp, &idx)); 533 if (isFE[field]) { 534 for (p = 0; p < n; ++p) if ((idx[p] < cStart) || (idx[p] >= cEnd)) ++newn; 535 } else { 536 for (p = 0; p < n; ++p) if ((idx[p] >= cStart) || (idx[p] < cEnd)) ++newn; 537 } 538 PetscCall(ISRestoreIndices(tmp, &idx)); 539 PetscCall(ISDestroy(&tmp)); 540 } 541 PetscCall(PetscMalloc1(newn, &newidx)); 542 newn = 0; 543 for (v = 0; v < numValues; ++v) { 544 IS tmp; 545 const PetscInt *idx; 546 547 PetscCall(DMLabelGetStratumIS(label, values[v], &tmp)); 548 if (!tmp) continue; 549 PetscCall(ISGetLocalSize(tmp, &n)); 550 PetscCall(ISGetIndices(tmp, &idx)); 551 if (isFE[field]) { 552 for (p = 0; p < n; ++p) if ((idx[p] < cStart) || (idx[p] >= cEnd)) newidx[newn++] = idx[p]; 553 } else { 554 for (p = 0; p < n; ++p) if ((idx[p] >= cStart) || (idx[p] < cEnd)) newidx[newn++] = idx[p]; 555 } 556 PetscCall(ISRestoreIndices(tmp, &idx)); 557 PetscCall(ISDestroy(&tmp)); 558 } 559 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, newn, newidx, PETSC_OWN_POINTER, &bcPoints[bc++])); 560 } 561 } 562 } 563 /* Handle discretization */ 564 PetscCall(PetscCalloc3(Nf,&labels,Nf,&numComp,Nf*(dim+1),&numDof)); 565 for (f = 0; f < Nf; ++f) { 566 labels[f] = dm->fields[f].label; 567 if (isFE[f]) { 568 PetscFE fe = (PetscFE) dm->fields[f].disc; 569 const PetscInt *numFieldDof; 570 PetscInt fedim, d; 571 572 PetscCall(PetscFEGetNumComponents(fe, &numComp[f])); 573 PetscCall(PetscFEGetNumDof(fe, &numFieldDof)); 574 PetscCall(PetscFEGetSpatialDimension(fe, &fedim)); 575 for (d = 0; d < PetscMin(dim, fedim)+1; ++d) numDof[f*(dim+1)+d] = numFieldDof[d]; 576 } else { 577 PetscFV fv = (PetscFV) dm->fields[f].disc; 578 579 PetscCall(PetscFVGetNumComponents(fv, &numComp[f])); 580 numDof[f*(dim+1)+dim] = numComp[f]; 581 } 582 } 583 PetscCall(DMPlexGetDepth(dm, &depth)); 584 for (f = 0; f < Nf; ++f) { 585 PetscInt d; 586 for (d = 1; d < dim; ++d) { 587 PetscCheckFalse((numDof[f*(dim+1)+d] > 0) && (depth < dim),PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated when unknowns are specified on edges or faces."); 588 } 589 } 590 PetscCall(DMPlexCreateSection(dm, labels, numComp, numDof, numBC, bcFields, bcComps, bcPoints, NULL, §ion)); 591 for (f = 0; f < Nf; ++f) { 592 PetscFE fe; 593 const char *name; 594 595 PetscCall(DMGetField(dm, f, NULL, (PetscObject *) &fe)); 596 if (!((PetscObject) fe)->name) continue; 597 PetscCall(PetscObjectGetName((PetscObject) fe, &name)); 598 PetscCall(PetscSectionSetFieldName(section, f, name)); 599 } 600 PetscCall(DMSetLocalSection(dm, section)); 601 PetscCall(PetscSectionDestroy(§ion)); 602 for (bc = 0; bc < numBC; ++bc) {PetscCall(ISDestroy(&bcPoints[bc]));PetscCall(ISDestroy(&bcComps[bc]));} 603 PetscCall(PetscFree3(bcFields,bcPoints,bcComps)); 604 PetscCall(PetscFree3(labels,numComp,numDof)); 605 PetscCall(PetscFree(isFE)); 606 PetscFunctionReturn(0); 607 } 608