1 #include <petsc/private/dmimpl.h> /*I "petscdm.h" I*/ 2 #include <petscds.h> 3 4 // Greatest common divisor of two nonnegative integers 5 PetscInt PetscGCD(PetscInt a, PetscInt b) 6 { 7 while (b != 0) { 8 PetscInt tmp = b; 9 b = a % b; 10 a = tmp; 11 } 12 return a; 13 } 14 15 PetscErrorCode DMCreateGlobalVector_Section_Private(DM dm, Vec *vec) 16 { 17 PetscSection gSection; 18 PetscInt localSize, bs, blockSize = -1, pStart, pEnd, p; 19 PetscInt in[2], out[2]; 20 21 PetscFunctionBegin; 22 PetscCall(DMGetGlobalSection(dm, &gSection)); 23 PetscCall(PetscSectionGetChart(gSection, &pStart, &pEnd)); 24 for (p = pStart; p < pEnd; ++p) { 25 PetscInt dof, cdof; 26 27 PetscCall(PetscSectionGetDof(gSection, p, &dof)); 28 PetscCall(PetscSectionGetConstraintDof(gSection, p, &cdof)); 29 30 if (dof - cdof > 0) { 31 if (blockSize < 0) { 32 /* set blockSize */ 33 blockSize = dof - cdof; 34 } else { 35 blockSize = PetscGCD(dof - cdof, blockSize); 36 } 37 } 38 } 39 40 in[0] = blockSize < 0 ? PETSC_MIN_INT : -blockSize; 41 in[1] = blockSize; 42 PetscCall(MPIU_Allreduce(in, out, 2, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm))); 43 /* -out[0] = min(blockSize), out[1] = max(blockSize) */ 44 if (-out[0] == out[1]) { 45 bs = out[1]; 46 } else bs = 1; 47 48 if (bs < 0) { /* Everyone was empty */ 49 blockSize = 1; 50 bs = 1; 51 } 52 53 PetscCall(PetscSectionGetConstrainedStorageSize(gSection, &localSize)); 54 PetscCheck(localSize % blockSize == 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mismatch between blocksize %" PetscInt_FMT " and local storage size %" PetscInt_FMT, blockSize, localSize); 55 PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), vec)); 56 PetscCall(VecSetSizes(*vec, localSize, PETSC_DETERMINE)); 57 PetscCall(VecSetBlockSize(*vec, bs)); 58 PetscCall(VecSetType(*vec, dm->vectype)); 59 PetscCall(VecSetDM(*vec, dm)); 60 /* PetscCall(VecSetLocalToGlobalMapping(*vec, dm->ltogmap)); */ 61 PetscFunctionReturn(PETSC_SUCCESS); 62 } 63 64 PetscErrorCode DMCreateLocalVector_Section_Private(DM dm, Vec *vec) 65 { 66 PetscSection section; 67 PetscInt localSize, blockSize = -1, pStart, pEnd, p; 68 69 PetscFunctionBegin; 70 PetscCall(DMGetLocalSection(dm, §ion)); 71 PetscCall(PetscSectionGetChart(section, &pStart, &pEnd)); 72 for (p = pStart; p < pEnd; ++p) { 73 PetscInt dof; 74 75 PetscCall(PetscSectionGetDof(section, p, &dof)); 76 if ((blockSize < 0) && (dof > 0)) blockSize = dof; 77 if (dof > 0) blockSize = PetscGCD(dof, blockSize); 78 } 79 PetscCall(PetscSectionGetStorageSize(section, &localSize)); 80 PetscCall(VecCreate(PETSC_COMM_SELF, vec)); 81 PetscCall(VecSetSizes(*vec, localSize, localSize)); 82 PetscCall(VecSetBlockSize(*vec, blockSize)); 83 PetscCall(VecSetType(*vec, dm->vectype)); 84 PetscCall(VecSetDM(*vec, dm)); 85 PetscFunctionReturn(PETSC_SUCCESS); 86 } 87 88 static PetscErrorCode PetscSectionSelectFields_Private(PetscSection s, PetscSection gs, PetscInt numFields, const PetscInt fields[], const PetscInt numComps[], const PetscInt comps[], IS *is) 89 { 90 PetscInt *subIndices; 91 PetscInt bs = 0, bsLocal[2], bsMinMax[2]; 92 PetscInt pStart, pEnd, Nc, subSize = 0, subOff = 0; 93 94 PetscFunctionBegin; 95 PetscCall(PetscSectionGetChart(gs, &pStart, &pEnd)); 96 if (numComps) { 97 for (PetscInt f = 0, off = 0; f < numFields; ++f) { 98 PetscInt Nc; 99 100 if (numComps[f] < 0) continue; 101 PetscCall(PetscSectionGetFieldComponents(s, f, &Nc)); 102 PetscCheck(numComps[f] <= Nc, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field %" PetscInt_FMT ": Number of selected components %" PetscInt_FMT " > %" PetscInt_FMT " number of field components", f, numComps[f], Nc); 103 for (PetscInt c = 0; c < numComps[f]; ++c, ++off) PetscCheck(comps[off] < Nc, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field %" PetscInt_FMT ": component %" PetscInt_FMT " not in [0, %" PetscInt_FMT ")", f, comps[off], Nc); 104 bs += numComps[f]; 105 } 106 } else { 107 for (PetscInt f = 0; f < numFields; ++f) { 108 PetscInt Nc; 109 110 PetscCall(PetscSectionGetFieldComponents(s, fields[f], &Nc)); 111 bs += Nc; 112 } 113 } 114 for (PetscInt p = pStart; p < pEnd; ++p) { 115 PetscInt gdof, pSubSize = 0; 116 117 PetscCall(PetscSectionGetDof(gs, p, &gdof)); 118 if (gdof > 0) { 119 PetscInt off = 0; 120 121 for (PetscInt f = 0; f < numFields; ++f) { 122 PetscInt fdof, fcdof, sfdof, sfcdof = 0; 123 124 PetscCall(PetscSectionGetFieldComponents(s, f, &Nc)); 125 PetscCall(PetscSectionGetFieldDof(s, p, fields[f], &fdof)); 126 PetscCall(PetscSectionGetFieldConstraintDof(s, p, fields[f], &fcdof)); 127 if (numComps && numComps[f] >= 0) { 128 const PetscInt *ind; 129 130 // Assume sets of dofs on points are of size Nc 131 PetscCheck(!(fdof % Nc), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components %" PetscInt_FMT " should evenly divide the dofs %" PetscInt_FMT " on point %" PetscInt_FMT, Nc, fdof, p); 132 sfdof = (fdof / Nc) * numComps[f]; 133 PetscCall(PetscSectionGetFieldConstraintIndices(s, p, fields[f], &ind)); 134 for (PetscInt i = 0; i < (fdof / Nc); ++i) { 135 for (PetscInt c = 0, fcc = 0; c < Nc; ++c) { 136 if (c == comps[off + fcc]) { 137 ++fcc; 138 ++sfcdof; 139 } 140 } 141 } 142 pSubSize += sfdof - sfcdof; 143 off += numComps[f]; 144 } else { 145 pSubSize += fdof - fcdof; 146 } 147 } 148 subSize += pSubSize; 149 if (pSubSize && bs != pSubSize) { 150 // Layout does not admit a pointwise block size 151 bs = 1; 152 } 153 } 154 } 155 // Must have same blocksize on all procs (some might have no points) 156 bsLocal[0] = bs < 0 ? PETSC_MAX_INT : bs; 157 bsLocal[1] = bs; 158 PetscCall(PetscGlobalMinMaxInt(PetscObjectComm((PetscObject)gs), bsLocal, bsMinMax)); 159 if (bsMinMax[0] != bsMinMax[1]) { 160 bs = 1; 161 } else { 162 bs = bsMinMax[0]; 163 } 164 PetscCall(PetscMalloc1(subSize, &subIndices)); 165 for (PetscInt p = pStart; p < pEnd; ++p) { 166 PetscInt gdof, goff; 167 168 PetscCall(PetscSectionGetDof(gs, p, &gdof)); 169 if (gdof > 0) { 170 PetscInt off = 0; 171 172 PetscCall(PetscSectionGetOffset(gs, p, &goff)); 173 for (PetscInt f = 0; f < numFields; ++f) { 174 PetscInt fdof, fcdof, poff = 0; 175 176 /* Can get rid of this loop by storing field information in the global section */ 177 for (PetscInt f2 = 0; f2 < fields[f]; ++f2) { 178 PetscCall(PetscSectionGetFieldDof(s, p, f2, &fdof)); 179 PetscCall(PetscSectionGetFieldConstraintDof(s, p, f2, &fcdof)); 180 poff += fdof - fcdof; 181 } 182 PetscCall(PetscSectionGetFieldDof(s, p, fields[f], &fdof)); 183 PetscCall(PetscSectionGetFieldConstraintDof(s, p, fields[f], &fcdof)); 184 185 if (numComps && numComps[f] >= 0) { 186 const PetscInt *ind; 187 188 // Assume sets of dofs on points are of size Nc 189 PetscCall(PetscSectionGetFieldConstraintIndices(s, p, fields[f], &ind)); 190 for (PetscInt i = 0, fcoff = 0, pfoff = 0; i < (fdof / Nc); ++i) { 191 for (PetscInt c = 0, fcc = 0; c < Nc; ++c) { 192 const PetscInt k = i * Nc + c; 193 194 if (ind[fcoff] == k) { 195 ++fcoff; 196 continue; 197 } 198 if (c == comps[off + fcc]) { 199 ++fcc; 200 subIndices[subOff++] = goff + poff + pfoff; 201 } 202 ++pfoff; 203 } 204 } 205 off += numComps[f]; 206 } else { 207 for (PetscInt fc = 0; fc < fdof - fcdof; ++fc, ++subOff) subIndices[subOff] = goff + poff + fc; 208 } 209 } 210 } 211 } 212 PetscCheck(subSize == subOff, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "The offset array size %" PetscInt_FMT " != %" PetscInt_FMT " the number of indices", subSize, subOff); 213 PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)gs), subSize, subIndices, PETSC_OWN_POINTER, is)); 214 if (bs > 1) { 215 // We need to check that the block size does not come from non-contiguous fields 216 PetscInt set = 1, rset = 1; 217 for (PetscInt i = 0; i < subSize; i += bs) { 218 for (PetscInt j = 0; j < bs; ++j) { 219 if (subIndices[i + j] != subIndices[i] + j) { 220 set = 0; 221 break; 222 } 223 } 224 } 225 PetscCall(MPIU_Allreduce(&set, &rset, 1, MPIU_INT, MPI_PROD, PetscObjectComm((PetscObject)gs))); 226 if (rset) PetscCall(ISSetBlockSize(*is, bs)); 227 } 228 PetscFunctionReturn(PETSC_SUCCESS); 229 } 230 231 static PetscErrorCode DMSelectFields_Private(DM dm, PetscSection section, PetscInt numFields, const PetscInt fields[], const PetscInt numComps[], const PetscInt comps[], IS *is, DM *subdm) 232 { 233 PetscSection subsection; 234 PetscBool haveNull = PETSC_FALSE; 235 PetscInt nf = 0, of = 0; 236 237 PetscFunctionBegin; 238 if (numComps) { 239 const PetscInt field = fields[0]; 240 241 PetscCheck(numFields == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "We only support a single field for component selection right now"); 242 PetscCall(PetscSectionCreateComponentSubsection(section, numComps[field], comps, &subsection)); 243 PetscCall(DMSetLocalSection(*subdm, subsection)); 244 PetscCall(PetscSectionDestroy(&subsection)); 245 (*subdm)->nullspaceConstructors[field] = dm->nullspaceConstructors[field]; 246 if (dm->probs) { 247 PetscFV fv, fvNew; 248 PetscInt fnum[1] = {field}; 249 250 PetscCall(DMSetNumFields(*subdm, 1)); 251 PetscCall(DMGetField(dm, field, NULL, (PetscObject *)&fv)); 252 PetscCall(PetscFVClone(fv, &fvNew)); 253 PetscCall(PetscFVSetNumComponents(fvNew, numComps[0])); 254 PetscCall(DMSetField(*subdm, 0, NULL, (PetscObject)fvNew)); 255 PetscCall(PetscFVDestroy(&fvNew)); 256 PetscCall(DMCreateDS(*subdm)); 257 if (numComps[0] == 1 && is) { 258 PetscObject disc, space, pmat; 259 260 PetscCall(DMGetField(*subdm, field, NULL, &disc)); 261 PetscCall(PetscObjectQuery(disc, "nullspace", &space)); 262 if (space) PetscCall(PetscObjectCompose((PetscObject)*is, "nullspace", space)); 263 PetscCall(PetscObjectQuery(disc, "nearnullspace", &space)); 264 if (space) PetscCall(PetscObjectCompose((PetscObject)*is, "nearnullspace", space)); 265 PetscCall(PetscObjectQuery(disc, "pmat", &pmat)); 266 if (pmat) PetscCall(PetscObjectCompose((PetscObject)*is, "pmat", pmat)); 267 } 268 PetscCall(PetscDSCopyConstants(dm->probs[field].ds, (*subdm)->probs[0].ds)); 269 PetscCall(PetscDSCopyBoundary(dm->probs[field].ds, 1, fnum, (*subdm)->probs[0].ds)); 270 PetscCall(PetscDSSelectEquations(dm->probs[field].ds, 1, fnum, (*subdm)->probs[0].ds)); 271 } 272 if ((*subdm)->nullspaceConstructors[0] && is) { 273 MatNullSpace nullSpace; 274 275 PetscCall((*(*subdm)->nullspaceConstructors[0])(*subdm, 0, 0, &nullSpace)); 276 PetscCall(PetscObjectCompose((PetscObject)*is, "nullspace", (PetscObject)nullSpace)); 277 PetscCall(MatNullSpaceDestroy(&nullSpace)); 278 } 279 if (dm->coarseMesh) PetscCall(DMCreateSubDM(dm->coarseMesh, numFields, fields, NULL, &(*subdm)->coarseMesh)); 280 PetscFunctionReturn(PETSC_SUCCESS); 281 } 282 PetscCall(PetscSectionCreateSubsection(section, numFields, fields, &subsection)); 283 PetscCall(DMSetLocalSection(*subdm, subsection)); 284 PetscCall(PetscSectionDestroy(&subsection)); 285 for (PetscInt f = 0; f < numFields; ++f) { 286 (*subdm)->nullspaceConstructors[f] = dm->nullspaceConstructors[fields[f]]; 287 if ((*subdm)->nullspaceConstructors[f]) { 288 haveNull = PETSC_TRUE; 289 nf = f; 290 of = fields[f]; 291 } 292 } 293 if (dm->probs) { 294 PetscCall(DMSetNumFields(*subdm, numFields)); 295 for (PetscInt f = 0; f < numFields; ++f) { 296 PetscObject disc; 297 298 PetscCall(DMGetField(dm, fields[f], NULL, &disc)); 299 PetscCall(DMSetField(*subdm, f, NULL, disc)); 300 } 301 // TODO: if only FV, then cut down the components 302 PetscCall(DMCreateDS(*subdm)); 303 if (numFields == 1 && is) { 304 PetscObject disc, space, pmat; 305 306 PetscCall(DMGetField(*subdm, 0, NULL, &disc)); 307 PetscCall(PetscObjectQuery(disc, "nullspace", &space)); 308 if (space) PetscCall(PetscObjectCompose((PetscObject)*is, "nullspace", space)); 309 PetscCall(PetscObjectQuery(disc, "nearnullspace", &space)); 310 if (space) PetscCall(PetscObjectCompose((PetscObject)*is, "nearnullspace", space)); 311 PetscCall(PetscObjectQuery(disc, "pmat", &pmat)); 312 if (pmat) PetscCall(PetscObjectCompose((PetscObject)*is, "pmat", pmat)); 313 } 314 // Check if DSes record their DM fields 315 if (dm->probs[0].fields) { 316 PetscInt d, e; 317 318 for (d = 0, e = 0; d < dm->Nds && e < (*subdm)->Nds; ++d) { 319 const PetscInt Nf = dm->probs[d].ds->Nf; 320 const PetscInt *fld; 321 PetscInt f, g; 322 323 PetscCall(ISGetIndices(dm->probs[d].fields, &fld)); 324 for (f = 0; f < Nf; ++f) { 325 for (g = 0; g < numFields; ++g) 326 if (fld[f] == fields[g]) break; 327 if (g < numFields) break; 328 } 329 PetscCall(ISRestoreIndices(dm->probs[d].fields, &fld)); 330 if (f == Nf) continue; 331 PetscCall(PetscDSCopyConstants(dm->probs[d].ds, (*subdm)->probs[e].ds)); 332 PetscCall(PetscDSCopyBoundary(dm->probs[d].ds, numFields, fields, (*subdm)->probs[e].ds)); 333 // Translate DM fields to DS fields 334 { 335 IS infields, dsfields; 336 const PetscInt *fld, *ofld; 337 PetscInt *fidx; 338 PetscInt onf, nf; 339 340 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, numFields, fields, PETSC_USE_POINTER, &infields)); 341 PetscCall(ISIntersect(infields, dm->probs[d].fields, &dsfields)); 342 PetscCall(ISDestroy(&infields)); 343 PetscCall(ISGetLocalSize(dsfields, &nf)); 344 PetscCheck(nf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "DS cannot be supported on 0 fields"); 345 PetscCall(ISGetIndices(dsfields, &fld)); 346 PetscCall(ISGetLocalSize(dm->probs[d].fields, &onf)); 347 PetscCall(ISGetIndices(dm->probs[d].fields, &ofld)); 348 PetscCall(PetscMalloc1(nf, &fidx)); 349 for (PetscInt f = 0, g = 0; f < onf && g < nf; ++f) { 350 if (ofld[f] == fld[g]) fidx[g++] = f; 351 } 352 PetscCall(ISRestoreIndices(dm->probs[d].fields, &ofld)); 353 PetscCall(ISRestoreIndices(dsfields, &fld)); 354 PetscCall(ISDestroy(&dsfields)); 355 PetscCall(PetscDSSelectDiscretizations(dm->probs[0].ds, nf, fidx, (*subdm)->probs[0].ds)); 356 PetscCall(PetscDSSelectEquations(dm->probs[0].ds, nf, fidx, (*subdm)->probs[0].ds)); 357 PetscCall(PetscFree(fidx)); 358 } 359 ++e; 360 } 361 } else { 362 PetscCall(PetscDSCopyConstants(dm->probs[0].ds, (*subdm)->probs[0].ds)); 363 PetscCall(PetscDSCopyBoundary(dm->probs[0].ds, PETSC_DETERMINE, NULL, (*subdm)->probs[0].ds)); 364 PetscCall(PetscDSSelectDiscretizations(dm->probs[0].ds, numFields, fields, (*subdm)->probs[0].ds)); 365 PetscCall(PetscDSSelectEquations(dm->probs[0].ds, numFields, fields, (*subdm)->probs[0].ds)); 366 } 367 } 368 if (haveNull && is) { 369 MatNullSpace nullSpace; 370 371 PetscCall((*(*subdm)->nullspaceConstructors[nf])(*subdm, of, nf, &nullSpace)); 372 PetscCall(PetscObjectCompose((PetscObject)*is, "nullspace", (PetscObject)nullSpace)); 373 PetscCall(MatNullSpaceDestroy(&nullSpace)); 374 } 375 if (dm->coarseMesh) PetscCall(DMCreateSubDM(dm->coarseMesh, numFields, fields, NULL, &(*subdm)->coarseMesh)); 376 PetscFunctionReturn(PETSC_SUCCESS); 377 } 378 379 /*@C 380 DMCreateSectionSubDM - Returns an `IS` and `subDM` containing a `PetscSection` that encapsulates a subproblem defined by a subset of the fields in a `PetscSection` in the `DM`. 381 382 Not Collective 383 384 Input Parameters: 385 + dm - The `DM` object 386 . numFields - The number of fields to incorporate into `subdm` 387 . fields - The field numbers of the selected fields 388 . numComps - The number of components from each field to incorporate into `subdm`, or PETSC_DECIDE for all components 389 - comps - The component numbers of the selected fields (omitted for PTESC_DECIDE fields) 390 391 Output Parameters: 392 + is - The global indices for the subproblem or `NULL` 393 - subdm - The `DM` for the subproblem, which must already have be cloned from `dm` or `NULL` 394 395 Level: intermediate 396 397 Notes: 398 If `is` and `subdm` are both `NULL` this does nothing 399 400 .seealso: `DMCreateSubDM()`, `DMGetLocalSection()`, `DMPlexSetMigrationSF()`, `DMView()` 401 @*/ 402 PetscErrorCode DMCreateSectionSubDM(DM dm, PetscInt numFields, const PetscInt fields[], const PetscInt numComps[], const PetscInt comps[], IS *is, DM *subdm) 403 { 404 PetscSection section, sectionGlobal; 405 PetscInt Nf; 406 407 PetscFunctionBegin; 408 if (!numFields) PetscFunctionReturn(PETSC_SUCCESS); 409 PetscCall(DMGetLocalSection(dm, §ion)); 410 PetscCall(DMGetGlobalSection(dm, §ionGlobal)); 411 PetscCheck(section, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields"); 412 PetscCheck(sectionGlobal, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields"); 413 PetscCall(PetscSectionGetNumFields(section, &Nf)); 414 PetscCheck(numFields <= Nf, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Number of requested fields %" PetscInt_FMT " greater than number of DM fields %" PetscInt_FMT, numFields, Nf); 415 416 if (is) PetscCall(PetscSectionSelectFields_Private(section, sectionGlobal, numFields, fields, numComps, comps, is)); 417 if (subdm) PetscCall(DMSelectFields_Private(dm, section, numFields, fields, numComps, comps, is, subdm)); 418 PetscFunctionReturn(PETSC_SUCCESS); 419 } 420 421 /*@C 422 DMCreateSectionSuperDM - Returns an arrays of `IS` and a `DM` containing a `PetscSection` that encapsulates a superproblem defined by the array of `DM` and their `PetscSection` 423 424 Not Collective 425 426 Input Parameters: 427 + dms - The `DM` objects, the must all have the same topology; for example obtained with `DMClone()` 428 - len - The number of `DM` in `dms` 429 430 Output Parameters: 431 + is - The global indices for the subproblem, or `NULL` 432 - superdm - The `DM` for the superproblem, which must already have be cloned and contain the same topology as the `dms` 433 434 Level: intermediate 435 436 .seealso: `DMCreateSuperDM()`, `DMGetLocalSection()`, `DMPlexSetMigrationSF()`, `DMView()` 437 @*/ 438 PetscErrorCode DMCreateSectionSuperDM(DM dms[], PetscInt len, IS **is, DM *superdm) 439 { 440 MPI_Comm comm; 441 PetscSection supersection, *sections, *sectionGlobals; 442 PetscInt *Nfs, Nf = 0, f, supf, oldf = -1, nullf = -1, i; 443 PetscBool haveNull = PETSC_FALSE; 444 445 PetscFunctionBegin; 446 PetscCall(PetscObjectGetComm((PetscObject)dms[0], &comm)); 447 /* Pull out local and global sections */ 448 PetscCall(PetscMalloc3(len, &Nfs, len, §ions, len, §ionGlobals)); 449 for (i = 0; i < len; ++i) { 450 PetscCall(DMGetLocalSection(dms[i], §ions[i])); 451 PetscCall(DMGetGlobalSection(dms[i], §ionGlobals[i])); 452 PetscCheck(sections[i], comm, PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields"); 453 PetscCheck(sectionGlobals[i], comm, PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields"); 454 PetscCall(PetscSectionGetNumFields(sections[i], &Nfs[i])); 455 Nf += Nfs[i]; 456 } 457 /* Create the supersection */ 458 PetscCall(PetscSectionCreateSupersection(sections, len, &supersection)); 459 PetscCall(DMSetLocalSection(*superdm, supersection)); 460 /* Create ISes */ 461 if (is) { 462 PetscSection supersectionGlobal; 463 PetscInt bs = -1, startf = 0; 464 465 PetscCall(PetscMalloc1(len, is)); 466 PetscCall(DMGetGlobalSection(*superdm, &supersectionGlobal)); 467 for (i = 0; i < len; startf += Nfs[i], ++i) { 468 PetscInt *subIndices; 469 PetscInt subSize, subOff, pStart, pEnd, p, start, end, dummy; 470 471 PetscCall(PetscSectionGetChart(sectionGlobals[i], &pStart, &pEnd)); 472 PetscCall(PetscSectionGetConstrainedStorageSize(sectionGlobals[i], &subSize)); 473 PetscCall(PetscMalloc1(subSize, &subIndices)); 474 for (p = pStart, subOff = 0; p < pEnd; ++p) { 475 PetscInt gdof, gcdof, gtdof, d; 476 477 PetscCall(PetscSectionGetDof(sectionGlobals[i], p, &gdof)); 478 PetscCall(PetscSectionGetConstraintDof(sections[i], p, &gcdof)); 479 gtdof = gdof - gcdof; 480 if (gdof > 0 && gtdof) { 481 if (bs < 0) { 482 bs = gtdof; 483 } else if (bs != gtdof) { 484 bs = 1; 485 } 486 PetscCall(DMGetGlobalFieldOffset_Private(*superdm, p, startf, &start, &dummy)); 487 PetscCall(DMGetGlobalFieldOffset_Private(*superdm, p, startf + Nfs[i] - 1, &dummy, &end)); 488 PetscCheck(end - start == gtdof, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid number of global dofs %" PetscInt_FMT " != %" PetscInt_FMT " for dm %" PetscInt_FMT " on point %" PetscInt_FMT, end - start, gtdof, i, p); 489 for (d = start; d < end; ++d, ++subOff) subIndices[subOff] = d; 490 } 491 } 492 PetscCall(ISCreateGeneral(comm, subSize, subIndices, PETSC_OWN_POINTER, &(*is)[i])); 493 /* Must have same blocksize on all procs (some might have no points) */ 494 { 495 PetscInt bs = -1, bsLocal[2], bsMinMax[2]; 496 497 bsLocal[0] = bs < 0 ? PETSC_MAX_INT : bs; 498 bsLocal[1] = bs; 499 PetscCall(PetscGlobalMinMaxInt(comm, bsLocal, bsMinMax)); 500 if (bsMinMax[0] != bsMinMax[1]) { 501 bs = 1; 502 } else { 503 bs = bsMinMax[0]; 504 } 505 PetscCall(ISSetBlockSize((*is)[i], bs)); 506 } 507 } 508 } 509 /* Preserve discretizations */ 510 if (len && dms[0]->probs) { 511 PetscCall(DMSetNumFields(*superdm, Nf)); 512 for (i = 0, supf = 0; i < len; ++i) { 513 for (f = 0; f < Nfs[i]; ++f, ++supf) { 514 PetscObject disc; 515 516 PetscCall(DMGetField(dms[i], f, NULL, &disc)); 517 PetscCall(DMSetField(*superdm, supf, NULL, disc)); 518 } 519 } 520 PetscCall(DMCreateDS(*superdm)); 521 } 522 /* Preserve nullspaces */ 523 for (i = 0, supf = 0; i < len; ++i) { 524 for (f = 0; f < Nfs[i]; ++f, ++supf) { 525 (*superdm)->nullspaceConstructors[supf] = dms[i]->nullspaceConstructors[f]; 526 if ((*superdm)->nullspaceConstructors[supf]) { 527 haveNull = PETSC_TRUE; 528 nullf = supf; 529 oldf = f; 530 } 531 } 532 } 533 /* Attach nullspace to IS */ 534 if (haveNull && is) { 535 MatNullSpace nullSpace; 536 537 PetscCall((*(*superdm)->nullspaceConstructors[nullf])(*superdm, oldf, nullf, &nullSpace)); 538 PetscCall(PetscObjectCompose((PetscObject)(*is)[nullf], "nullspace", (PetscObject)nullSpace)); 539 PetscCall(MatNullSpaceDestroy(&nullSpace)); 540 } 541 PetscCall(PetscSectionDestroy(&supersection)); 542 PetscCall(PetscFree3(Nfs, sections, sectionGlobals)); 543 PetscFunctionReturn(PETSC_SUCCESS); 544 } 545