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