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 if (bcPoints[bc]) { 206 PetscCall(ISGetLocalSize(bcPoints[bc], &n)); 207 PetscCall(ISGetIndices(bcPoints[bc], &idx)); 208 for (i = 0; i < n; ++i) { 209 const PetscInt p = idx[i]; 210 PetscInt numConst; 211 212 if (Nf) { 213 PetscCall(PetscSectionGetFieldDof(section, p, field, &numConst)); 214 } else { 215 PetscCall(PetscSectionGetDof(section, p, &numConst)); 216 } 217 /* If Nc <= 0, constrain every dof on the point */ 218 if (cNc > 0) { 219 /* We assume that a point may have multiple "nodes", which are collections of Nc dofs, 220 and that those dofs are numbered n*Nc+c */ 221 if (Nf) { 222 PetscCheck(numConst % Nc == 0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Point %" PetscInt_FMT " has %" PetscInt_FMT " dof which is not divisible by %" PetscInt_FMT " field components", p, numConst, Nc); 223 numConst = (numConst/Nc) * cNc; 224 } else { 225 numConst = PetscMin(numConst, cNc); 226 } 227 } 228 if (Nf) PetscCall(PetscSectionAddFieldConstraintDof(section, p, field, numConst)); 229 PetscCall(PetscSectionAddConstraintDof(section, p, numConst)); 230 } 231 PetscCall(ISRestoreIndices(bcPoints[bc], &idx)); 232 } 233 if (bcComps && bcComps[bc]) PetscCall(ISRestoreIndices(bcComps[bc], &comp)); 234 } 235 PetscCall(DMPlexGetAnchors(dm, &aSec, NULL)); 236 if (aSec) { 237 PetscInt aStart, aEnd, a; 238 239 PetscCall(PetscSectionGetChart(aSec, &aStart, &aEnd)); 240 for (a = aStart; a < aEnd; a++) { 241 PetscInt dof, f; 242 243 PetscCall(PetscSectionGetDof(aSec, a, &dof)); 244 if (dof) { 245 /* if there are point-to-point constraints, then all dofs are constrained */ 246 PetscCall(PetscSectionGetDof(section, a, &dof)); 247 PetscCall(PetscSectionSetConstraintDof(section, a, dof)); 248 for (f = 0; f < Nf; f++) { 249 PetscCall(PetscSectionGetFieldDof(section, a, f, &dof)); 250 PetscCall(PetscSectionSetFieldConstraintDof(section, a, f, dof)); 251 } 252 } 253 } 254 } 255 PetscFunctionReturn(0); 256 } 257 258 /* Set the constrained field indices on each point 259 If bcComps is NULL or the IS is NULL, constrain every dof on the point 260 */ 261 static PetscErrorCode DMPlexCreateSectionBCIndicesField(DM dm, PetscInt numBC,const PetscInt bcField[], const IS bcComps[], const IS bcPoints[], PetscSection section) 262 { 263 PetscSection aSec; 264 PetscInt *indices; 265 PetscInt Nf, cdof, maxDof = 0, pStart, pEnd, p, bc, f, d; 266 267 PetscFunctionBegin; 268 PetscCall(PetscSectionGetNumFields(section, &Nf)); 269 if (!Nf) PetscFunctionReturn(0); 270 /* Initialize all field indices to -1 */ 271 PetscCall(PetscSectionGetChart(section, &pStart, &pEnd)); 272 for (p = pStart; p < pEnd; ++p) {PetscCall(PetscSectionGetConstraintDof(section, p, &cdof)); maxDof = PetscMax(maxDof, cdof);} 273 PetscCall(PetscMalloc1(maxDof, &indices)); 274 for (d = 0; d < maxDof; ++d) indices[d] = -1; 275 for (p = pStart; p < pEnd; ++p) for (f = 0; f < Nf; ++f) PetscCall(PetscSectionSetFieldConstraintIndices(section, p, f, indices)); 276 /* Handle BC constraints */ 277 for (bc = 0; bc < numBC; ++bc) { 278 const PetscInt field = bcField[bc]; 279 const PetscInt *comp, *idx; 280 PetscInt Nc, cNc = -1, n, i; 281 282 PetscCall(PetscSectionGetFieldComponents(section, field, &Nc)); 283 if (bcComps && bcComps[bc]) PetscCall(ISGetLocalSize(bcComps[bc], &cNc)); 284 if (bcComps && bcComps[bc]) PetscCall(ISGetIndices(bcComps[bc], &comp)); 285 if (bcPoints[bc]) { 286 PetscCall(ISGetLocalSize(bcPoints[bc], &n)); 287 PetscCall(ISGetIndices(bcPoints[bc], &idx)); 288 for (i = 0; i < n; ++i) { 289 const PetscInt p = idx[i]; 290 const PetscInt *find; 291 PetscInt fdof, fcdof, c, j; 292 293 PetscCall(PetscSectionGetFieldDof(section, p, field, &fdof)); 294 if (!fdof) continue; 295 if (cNc < 0) { 296 for (d = 0; d < fdof; ++d) indices[d] = d; 297 fcdof = fdof; 298 } else { 299 /* We assume that a point may have multiple "nodes", which are collections of Nc dofs, 300 and that those dofs are numbered n*Nc+c */ 301 PetscCall(PetscSectionGetFieldConstraintDof(section, p, field, &fcdof)); 302 PetscCall(PetscSectionGetFieldConstraintIndices(section, p, field, &find)); 303 /* Get indices constrained by previous bcs */ 304 for (d = 0; d < fcdof; ++d) {if (find[d] < 0) break; indices[d] = find[d];} 305 for (j = 0; j < fdof/Nc; ++j) for (c = 0; c < cNc; ++c) indices[d++] = j*Nc + comp[c]; 306 PetscCall(PetscSortRemoveDupsInt(&d, indices)); 307 for (c = d; c < fcdof; ++c) indices[c] = -1; 308 fcdof = d; 309 } 310 PetscCall(PetscSectionSetFieldConstraintDof(section, p, field, fcdof)); 311 PetscCall(PetscSectionSetFieldConstraintIndices(section, p, field, indices)); 312 } 313 PetscCall(ISRestoreIndices(bcPoints[bc], &idx)); 314 } 315 if (bcComps && bcComps[bc]) PetscCall(ISRestoreIndices(bcComps[bc], &comp)); 316 } 317 /* Handle anchors */ 318 PetscCall(DMPlexGetAnchors(dm, &aSec, NULL)); 319 if (aSec) { 320 PetscInt aStart, aEnd, a; 321 322 for (d = 0; d < maxDof; ++d) indices[d] = d; 323 PetscCall(PetscSectionGetChart(aSec, &aStart, &aEnd)); 324 for (a = aStart; a < aEnd; a++) { 325 PetscInt dof, f; 326 327 PetscCall(PetscSectionGetDof(aSec, a, &dof)); 328 if (dof) { 329 /* if there are point-to-point constraints, then all dofs are constrained */ 330 for (f = 0; f < Nf; f++) { 331 PetscCall(PetscSectionSetFieldConstraintIndices(section, a, f, indices)); 332 } 333 } 334 } 335 } 336 PetscCall(PetscFree(indices)); 337 PetscFunctionReturn(0); 338 } 339 340 /* Set the constrained indices on each point */ 341 static PetscErrorCode DMPlexCreateSectionBCIndices(DM dm, PetscSection section) 342 { 343 PetscInt *indices; 344 PetscInt Nf, maxDof, pStart, pEnd, p, f, d; 345 346 PetscFunctionBegin; 347 PetscCall(PetscSectionGetNumFields(section, &Nf)); 348 PetscCall(PetscSectionGetMaxDof(section, &maxDof)); 349 PetscCall(PetscSectionGetChart(section, &pStart, &pEnd)); 350 PetscCall(PetscMalloc1(maxDof, &indices)); 351 for (d = 0; d < maxDof; ++d) indices[d] = -1; 352 for (p = pStart; p < pEnd; ++p) { 353 PetscInt cdof, d; 354 355 PetscCall(PetscSectionGetConstraintDof(section, p, &cdof)); 356 if (cdof) { 357 if (Nf) { 358 PetscInt numConst = 0, foff = 0; 359 360 for (f = 0; f < Nf; ++f) { 361 const PetscInt *find; 362 PetscInt fcdof, fdof; 363 364 PetscCall(PetscSectionGetFieldDof(section, p, f, &fdof)); 365 PetscCall(PetscSectionGetFieldConstraintDof(section, p, f, &fcdof)); 366 /* Change constraint numbering from field component to local dof number */ 367 PetscCall(PetscSectionGetFieldConstraintIndices(section, p, f, &find)); 368 for (d = 0; d < fcdof; ++d) indices[numConst+d] = find[d] + foff; 369 numConst += fcdof; 370 foff += fdof; 371 } 372 if (cdof != numConst) PetscCall(PetscSectionSetConstraintDof(section, p, numConst)); 373 } else { 374 for (d = 0; d < cdof; ++d) indices[d] = d; 375 } 376 PetscCall(PetscSectionSetConstraintIndices(section, p, indices)); 377 } 378 } 379 PetscCall(PetscFree(indices)); 380 PetscFunctionReturn(0); 381 } 382 383 /*@C 384 DMPlexCreateSection - Create a PetscSection based upon the dof layout specification provided. 385 386 Not Collective 387 388 Input Parameters: 389 + dm - The DMPlex object 390 . label - The label indicating the mesh support of each field, or NULL for the whole mesh 391 . numComp - An array of size numFields that holds the number of components for each field 392 . numDof - An array of size numFields*(dim+1) which holds the number of dof for each field on a mesh piece of dimension d 393 . numBC - The number of boundary conditions 394 . bcField - An array of size numBC giving the field number for each boundry condition 395 . bcComps - [Optional] An array of size numBC giving an IS holding the field components to which each boundary condition applies 396 . bcPoints - An array of size numBC giving an IS holding the Plex points to which each boundary condition applies 397 - perm - Optional permutation of the chart, or NULL 398 399 Output Parameter: 400 . section - The PetscSection object 401 402 Notes: 403 numDof[f*(dim+1)+d] gives the number of dof for field f on points of dimension d. For instance, numDof[1] is the 404 number of dof for field 0 on each edge. 405 406 The chart permutation is the same one set using PetscSectionSetPermutation() 407 408 Level: developer 409 410 TODO: How is this related to DMCreateLocalSection() 411 412 .seealso: `DMPlexCreate()`, `PetscSectionCreate()`, `PetscSectionSetPermutation()` 413 @*/ 414 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) 415 { 416 PetscSection aSec; 417 418 PetscFunctionBegin; 419 PetscCall(DMPlexCreateSectionFields(dm, numComp, section)); 420 PetscCall(DMPlexCreateSectionDof(dm, label, numDof, *section)); 421 PetscCall(DMPlexCreateSectionBCDof(dm, numBC, bcField, bcComps, bcPoints, *section)); 422 if (perm) PetscCall(PetscSectionSetPermutation(*section, perm)); 423 PetscCall(PetscSectionSetFromOptions(*section)); 424 PetscCall(PetscSectionSetUp(*section)); 425 PetscCall(DMPlexGetAnchors(dm,&aSec,NULL)); 426 if (numBC || aSec) { 427 PetscCall(DMPlexCreateSectionBCIndicesField(dm, numBC, bcField, bcComps, bcPoints, *section)); 428 PetscCall(DMPlexCreateSectionBCIndices(dm, *section)); 429 } 430 PetscCall(PetscSectionViewFromOptions(*section,NULL,"-section_view")); 431 PetscFunctionReturn(0); 432 } 433 434 static PetscErrorCode DMCompleteBCLabels_Private(DM dm, const PetscBool isFE[]) 435 { 436 DMLabel *labels, *glabels; 437 const char **names; 438 char *sendNames, *recvNames; 439 PetscInt Nds, s, maxLabels = 0, maxLen = 0, gmaxLen, Nl = 0, gNl, l, gl, m; 440 size_t len; 441 MPI_Comm comm; 442 PetscMPIInt rank, size, p, *counts, *displs; 443 444 PetscFunctionBegin; 445 PetscCall(PetscObjectGetComm((PetscObject) dm, &comm)); 446 PetscCallMPI(MPI_Comm_size(comm, &size)); 447 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 448 PetscCall(DMGetNumDS(dm, &Nds)); 449 for (s = 0; s < Nds; ++s) { 450 PetscDS dsBC; 451 PetscInt numBd; 452 453 PetscCall(DMGetRegionNumDS(dm, s, NULL, NULL, &dsBC)); 454 PetscCall(PetscDSGetNumBoundary(dsBC, &numBd)); 455 maxLabels += numBd; 456 } 457 PetscCall(PetscCalloc1(maxLabels, &labels)); 458 /* Get list of labels to be completed */ 459 for (s = 0; s < Nds; ++s) { 460 PetscDS dsBC; 461 PetscInt numBd, bd; 462 463 PetscCall(DMGetRegionNumDS(dm, s, NULL, NULL, &dsBC)); 464 PetscCall(PetscDSGetNumBoundary(dsBC, &numBd)); 465 for (bd = 0; bd < numBd; ++bd) { 466 DMLabel label; 467 PetscInt field; 468 469 PetscCall(PetscDSGetBoundary(dsBC, bd, NULL, NULL, NULL, &label, NULL, NULL, &field, NULL, NULL, NULL, NULL, NULL)); 470 if (!isFE[field] || !label) continue; 471 for (l = 0; l < Nl; ++l) if (labels[l] == label) break; 472 if (l == Nl) labels[Nl++] = label; 473 } 474 } 475 /* Get label names */ 476 PetscCall(PetscMalloc1(Nl, &names)); 477 for (l = 0; l < Nl; ++l) PetscCall(PetscObjectGetName((PetscObject) labels[l], &names[l])); 478 for (l = 0; l < Nl; ++l) {PetscCall(PetscStrlen(names[l], &len)); maxLen = PetscMax(maxLen, len+2);} 479 PetscCall(PetscFree(labels)); 480 PetscCallMPI(MPI_Allreduce(&maxLen, &gmaxLen, 1, MPIU_INT, MPI_MAX, comm)); 481 PetscCall(PetscMalloc1(Nl * gmaxLen, &sendNames)); 482 for (l = 0; l < Nl; ++l) PetscCall(PetscStrcpy(&sendNames[gmaxLen*l], names[l])); 483 PetscCall(PetscFree(names)); 484 /* Put all names on all processes */ 485 PetscCall(PetscCalloc2(size, &counts, size+1, &displs)); 486 PetscCallMPI(MPI_Allgather(&Nl, 1, MPIU_INT, counts, 1, MPIU_INT, comm)); 487 for (p = 0; p < size; ++p) displs[p+1] = displs[p] + counts[p]; 488 gNl = displs[size]; 489 for (p = 0; p < size; ++p) {counts[p] *= gmaxLen; displs[p] *= gmaxLen;} 490 PetscCall(PetscMalloc2(gNl * gmaxLen, &recvNames, gNl, &glabels)); 491 PetscCallMPI(MPI_Allgatherv(sendNames, counts[rank], MPI_CHAR, recvNames, counts, displs, MPI_CHAR, comm)); 492 PetscCall(PetscFree2(counts, displs)); 493 PetscCall(PetscFree(sendNames)); 494 for (l = 0, gl = 0; l < gNl; ++l) { 495 PetscCall(DMGetLabel(dm, &recvNames[l*gmaxLen], &glabels[gl])); 496 PetscCheck(glabels[gl], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Label %s missing on rank %d", &recvNames[l*gmaxLen], rank); 497 for (m = 0; m < gl; ++m) if (glabels[m] == glabels[gl]) continue; 498 PetscCall(DMPlexLabelComplete(dm, glabels[gl])); 499 ++gl; 500 } 501 PetscCall(PetscFree2(recvNames, glabels)); 502 PetscFunctionReturn(0); 503 } 504 505 PetscErrorCode DMCreateLocalSection_Plex(DM dm) 506 { 507 PetscSection section; 508 DMLabel *labels; 509 IS *bcPoints, *bcComps; 510 PetscBool *isFE; 511 PetscInt *bcFields, *numComp, *numDof; 512 PetscInt depth, dim, numBC = 0, Nf, Nds, s, bc = 0, f; 513 PetscInt cStart, cEnd, cEndInterior; 514 515 PetscFunctionBegin; 516 PetscCall(DMGetNumFields(dm, &Nf)); 517 PetscCall(DMGetDimension(dm, &dim)); 518 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 519 /* FE and FV boundary conditions are handled slightly differently */ 520 PetscCall(PetscMalloc1(Nf, &isFE)); 521 for (f = 0; f < Nf; ++f) { 522 PetscObject obj; 523 PetscClassId id; 524 525 PetscCall(DMGetField(dm, f, NULL, &obj)); 526 PetscCall(PetscObjectGetClassId(obj, &id)); 527 if (id == PETSCFE_CLASSID) {isFE[f] = PETSC_TRUE;} 528 else if (id == PETSCFV_CLASSID) {isFE[f] = PETSC_FALSE;} 529 else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, f); 530 } 531 /* Allocate boundary point storage for FEM boundaries */ 532 PetscCall(DMGetNumDS(dm, &Nds)); 533 for (s = 0; s < Nds; ++s) { 534 PetscDS dsBC; 535 PetscInt numBd, bd; 536 537 PetscCall(DMGetRegionNumDS(dm, s, NULL, NULL, &dsBC)); 538 PetscCall(PetscDSGetNumBoundary(dsBC, &numBd)); 539 PetscCheck(Nf || !numBd,PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "number of fields is zero and number of boundary conditions is nonzero (this should never happen)"); 540 for (bd = 0; bd < numBd; ++bd) { 541 PetscInt field; 542 DMBoundaryConditionType type; 543 DMLabel label; 544 545 PetscCall(PetscDSGetBoundary(dsBC, bd, NULL, &type, NULL, &label, NULL, NULL, &field, NULL, NULL, NULL, NULL, NULL)); 546 if (label && isFE[field] && (type & DM_BC_ESSENTIAL)) ++numBC; 547 } 548 } 549 /* Add ghost cell boundaries for FVM */ 550 PetscCall(DMPlexGetGhostCellStratum(dm, &cEndInterior, NULL)); 551 for (f = 0; f < Nf; ++f) if (!isFE[f] && cEndInterior >= 0) ++numBC; 552 PetscCall(PetscCalloc3(numBC, &bcFields, numBC, &bcPoints, numBC, &bcComps)); 553 /* Constrain ghost cells for FV */ 554 for (f = 0; f < Nf; ++f) { 555 PetscInt *newidx, c; 556 557 if (isFE[f] || cEndInterior < 0) continue; 558 PetscCall(PetscMalloc1(cEnd-cEndInterior,&newidx)); 559 for (c = cEndInterior; c < cEnd; ++c) newidx[c-cEndInterior] = c; 560 bcFields[bc] = f; 561 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, cEnd-cEndInterior, newidx, PETSC_OWN_POINTER, &bcPoints[bc++])); 562 } 563 /* Complete labels for boundary conditions */ 564 PetscCall(DMCompleteBCLabels_Private(dm, isFE)); 565 /* Handle FEM Dirichlet boundaries */ 566 PetscCall(DMGetNumDS(dm, &Nds)); 567 for (s = 0; s < Nds; ++s) { 568 PetscDS dsBC; 569 PetscInt numBd, bd; 570 571 PetscCall(DMGetRegionNumDS(dm, s, NULL, NULL, &dsBC)); 572 PetscCall(PetscDSGetNumBoundary(dsBC, &numBd)); 573 for (bd = 0; bd < numBd; ++bd) { 574 DMLabel label; 575 const PetscInt *comps; 576 const PetscInt *values; 577 PetscInt bd2, field, numComps, numValues; 578 DMBoundaryConditionType type; 579 PetscBool duplicate = PETSC_FALSE; 580 581 PetscCall(PetscDSGetBoundary(dsBC, bd, NULL, &type, NULL, &label, &numValues, &values, &field, &numComps, &comps, NULL, NULL, NULL)); 582 if (!isFE[field] || !label) continue; 583 /* Only want to modify label once */ 584 for (bd2 = 0; bd2 < bd; ++bd2) { 585 DMLabel l; 586 587 PetscCall(PetscDSGetBoundary(dsBC, bd2, NULL, NULL, NULL, &l, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 588 duplicate = l == label ? PETSC_TRUE : PETSC_FALSE; 589 if (duplicate) break; 590 } 591 /* Filter out cells, if you actually want to constrain cells you need to do things by hand right now */ 592 if (type & DM_BC_ESSENTIAL) { 593 PetscInt *newidx; 594 PetscInt n, newn = 0, p, v; 595 596 bcFields[bc] = field; 597 if (numComps) PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject) dm), numComps, comps, PETSC_COPY_VALUES, &bcComps[bc])); 598 for (v = 0; v < numValues; ++v) { 599 IS tmp; 600 const PetscInt *idx; 601 602 PetscCall(DMLabelGetStratumIS(label, values[v], &tmp)); 603 if (!tmp) continue; 604 PetscCall(ISGetLocalSize(tmp, &n)); 605 PetscCall(ISGetIndices(tmp, &idx)); 606 if (isFE[field]) { 607 for (p = 0; p < n; ++p) if ((idx[p] < cStart) || (idx[p] >= cEnd)) ++newn; 608 } else { 609 for (p = 0; p < n; ++p) if ((idx[p] >= cStart) || (idx[p] < cEnd)) ++newn; 610 } 611 PetscCall(ISRestoreIndices(tmp, &idx)); 612 PetscCall(ISDestroy(&tmp)); 613 } 614 PetscCall(PetscMalloc1(newn, &newidx)); 615 newn = 0; 616 for (v = 0; v < numValues; ++v) { 617 IS tmp; 618 const PetscInt *idx; 619 620 PetscCall(DMLabelGetStratumIS(label, values[v], &tmp)); 621 if (!tmp) continue; 622 PetscCall(ISGetLocalSize(tmp, &n)); 623 PetscCall(ISGetIndices(tmp, &idx)); 624 if (isFE[field]) { 625 for (p = 0; p < n; ++p) if ((idx[p] < cStart) || (idx[p] >= cEnd)) newidx[newn++] = idx[p]; 626 } else { 627 for (p = 0; p < n; ++p) if ((idx[p] >= cStart) || (idx[p] < cEnd)) newidx[newn++] = idx[p]; 628 } 629 PetscCall(ISRestoreIndices(tmp, &idx)); 630 PetscCall(ISDestroy(&tmp)); 631 } 632 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, newn, newidx, PETSC_OWN_POINTER, &bcPoints[bc++])); 633 } 634 } 635 } 636 /* Handle discretization */ 637 PetscCall(PetscCalloc3(Nf,&labels,Nf,&numComp,Nf*(dim+1),&numDof)); 638 for (f = 0; f < Nf; ++f) { 639 labels[f] = dm->fields[f].label; 640 if (isFE[f]) { 641 PetscFE fe = (PetscFE) dm->fields[f].disc; 642 const PetscInt *numFieldDof; 643 PetscInt fedim, d; 644 645 PetscCall(PetscFEGetNumComponents(fe, &numComp[f])); 646 PetscCall(PetscFEGetNumDof(fe, &numFieldDof)); 647 PetscCall(PetscFEGetSpatialDimension(fe, &fedim)); 648 for (d = 0; d < PetscMin(dim, fedim)+1; ++d) numDof[f*(dim+1)+d] = numFieldDof[d]; 649 } else { 650 PetscFV fv = (PetscFV) dm->fields[f].disc; 651 652 PetscCall(PetscFVGetNumComponents(fv, &numComp[f])); 653 numDof[f*(dim+1)+dim] = numComp[f]; 654 } 655 } 656 PetscCall(DMPlexGetDepth(dm, &depth)); 657 for (f = 0; f < Nf; ++f) { 658 PetscInt d; 659 for (d = 1; d < dim; ++d) { 660 PetscCheck(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."); 661 } 662 } 663 PetscCall(DMPlexCreateSection(dm, labels, numComp, numDof, numBC, bcFields, bcComps, bcPoints, NULL, §ion)); 664 for (f = 0; f < Nf; ++f) { 665 PetscFE fe; 666 const char *name; 667 668 PetscCall(DMGetField(dm, f, NULL, (PetscObject *) &fe)); 669 if (!((PetscObject) fe)->name) continue; 670 PetscCall(PetscObjectGetName((PetscObject) fe, &name)); 671 PetscCall(PetscSectionSetFieldName(section, f, name)); 672 } 673 PetscCall(DMSetLocalSection(dm, section)); 674 PetscCall(PetscSectionDestroy(§ion)); 675 for (bc = 0; bc < numBC; ++bc) {PetscCall(ISDestroy(&bcPoints[bc]));PetscCall(ISDestroy(&bcComps[bc]));} 676 PetscCall(PetscFree3(bcFields,bcPoints,bcComps)); 677 PetscCall(PetscFree3(labels,numComp,numDof)); 678 PetscCall(PetscFree(isFE)); 679 PetscFunctionReturn(0); 680 } 681