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