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