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