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