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, PetscBool useComp, PetscInt numFields, const PetscInt fields[], IS *is) 89 { 90 PetscInt *subIndices; 91 PetscInt bs, bsLocal[2], bsMinMax[2]; 92 PetscInt pStart, pEnd, subSize = 0, subOff = 0; 93 94 PetscFunctionBegin; 95 PetscCall(PetscSectionGetChart(gs, &pStart, &pEnd)); 96 if (useComp) { 97 PetscInt Nc; 98 99 PetscCall(PetscSectionGetFieldComponents(s, 0, &Nc)); 100 for (PetscInt f = 0; f < numFields; ++f) PetscCheck(fields[f] < Nc, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "FV component[%" PetscInt_FMT "]: %" PetscInt_FMT " not in [0, %" PetscInt_FMT ")", f, fields[f], Nc); 101 bs = numFields; 102 for (PetscInt p = pStart; p < pEnd; ++p) { 103 PetscInt gdof, fcdof; 104 105 PetscCall(PetscSectionGetDof(gs, p, &gdof)); 106 if (gdof > 0) { 107 PetscCheck(gdof == Nc, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "This is collocated FV, so there should be %" PetscInt_FMT " dofs in each cell, not %" PetscInt_FMT, Nc, gdof); 108 PetscCall(PetscSectionGetFieldConstraintDof(s, p, 0, &fcdof)); 109 PetscCheck(!fcdof, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "This is collocated FV, so there should be 0 constraints in each cell, not %" PetscInt_FMT, fcdof); 110 subSize += numFields; 111 } 112 } 113 PetscCall(PetscMalloc1(subSize, &subIndices)); 114 for (PetscInt p = pStart; p < pEnd; ++p) { 115 PetscInt gdof, goff; 116 117 PetscCall(PetscSectionGetDof(gs, p, &gdof)); 118 if (gdof > 0) { 119 PetscCall(PetscSectionGetOffset(gs, p, &goff)); 120 for (PetscInt f = 0; f < numFields; ++f, ++subOff) subIndices[subOff] = goff + fields[f]; 121 } 122 } 123 PetscCheck(subSize == subOff, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "The offset array size %" PetscInt_FMT " != %" PetscInt_FMT " the number of indices", subSize, subOff); 124 PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)gs), subSize, subIndices, PETSC_OWN_POINTER, is)); 125 PetscCall(ISSetBlockSize(*is, bs)); 126 PetscFunctionReturn(PETSC_SUCCESS); 127 } 128 bs = 0; 129 for (PetscInt f = 0; f < numFields; ++f) { 130 PetscInt Nc; 131 132 PetscCall(PetscSectionGetFieldComponents(s, fields[f], &Nc)); 133 bs += Nc; 134 } 135 for (PetscInt p = pStart; p < pEnd; ++p) { 136 PetscInt gdof, pSubSize = 0; 137 138 PetscCall(PetscSectionGetDof(gs, p, &gdof)); 139 if (gdof > 0) { 140 for (PetscInt f = 0; f < numFields; ++f) { 141 PetscInt fdof, fcdof; 142 143 PetscCall(PetscSectionGetFieldDof(s, p, fields[f], &fdof)); 144 PetscCall(PetscSectionGetFieldConstraintDof(s, p, fields[f], &fcdof)); 145 pSubSize += fdof - fcdof; 146 } 147 subSize += pSubSize; 148 if (pSubSize && bs != pSubSize) { 149 // Layout does not admit a pointwise block size 150 bs = 1; 151 } 152 } 153 } 154 // Must have same blocksize on all procs (some might have no points) 155 bsLocal[0] = bs < 0 ? PETSC_MAX_INT : bs; 156 bsLocal[1] = bs; 157 PetscCall(PetscGlobalMinMaxInt(PetscObjectComm((PetscObject)gs), bsLocal, bsMinMax)); 158 if (bsMinMax[0] != bsMinMax[1]) { 159 bs = 1; 160 } else { 161 bs = bsMinMax[0]; 162 } 163 PetscCall(PetscMalloc1(subSize, &subIndices)); 164 for (PetscInt p = pStart; p < pEnd; ++p) { 165 PetscInt gdof, goff; 166 167 PetscCall(PetscSectionGetDof(gs, p, &gdof)); 168 if (gdof > 0) { 169 PetscCall(PetscSectionGetOffset(gs, p, &goff)); 170 for (PetscInt f = 0; f < numFields; ++f) { 171 PetscInt fdof, fcdof, poff = 0; 172 173 /* Can get rid of this loop by storing field information in the global section */ 174 for (PetscInt f2 = 0; f2 < fields[f]; ++f2) { 175 PetscCall(PetscSectionGetFieldDof(s, p, f2, &fdof)); 176 PetscCall(PetscSectionGetFieldConstraintDof(s, p, f2, &fcdof)); 177 poff += fdof - fcdof; 178 } 179 PetscCall(PetscSectionGetFieldDof(s, p, fields[f], &fdof)); 180 PetscCall(PetscSectionGetFieldConstraintDof(s, p, fields[f], &fcdof)); 181 for (PetscInt fc = 0; fc < fdof - fcdof; ++fc, ++subOff) subIndices[subOff] = goff + poff + fc; 182 } 183 } 184 } 185 PetscCheck(subSize == subOff, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "The offset array size %" PetscInt_FMT " != %" PetscInt_FMT " the number of indices", subSize, subOff); 186 PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)gs), subSize, subIndices, PETSC_OWN_POINTER, is)); 187 if (bs > 1) { 188 // We need to check that the block size does not come from non-contiguous fields 189 PetscInt set = 1, rset = 1; 190 for (PetscInt i = 0; i < subSize; i += bs) { 191 for (PetscInt j = 0; j < bs; ++j) { 192 if (subIndices[i + j] != subIndices[i] + j) { 193 set = 0; 194 break; 195 } 196 } 197 } 198 PetscCall(MPIU_Allreduce(&set, &rset, 1, MPIU_INT, MPI_PROD, PetscObjectComm((PetscObject)gs))); 199 if (rset) PetscCall(ISSetBlockSize(*is, bs)); 200 } 201 PetscFunctionReturn(PETSC_SUCCESS); 202 } 203 204 static PetscErrorCode DMSelectFields_Private(DM dm, PetscSection section, PetscBool useComp, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm) 205 { 206 PetscSection subsection; 207 PetscBool haveNull = PETSC_FALSE; 208 PetscInt nf = 0, of = 0; 209 210 PetscFunctionBegin; 211 if (useComp) { 212 PetscCall(PetscSectionCreateComponentSubsection(section, numFields, fields, &subsection)); 213 PetscCall(DMSetLocalSection(*subdm, subsection)); 214 PetscCall(PetscSectionDestroy(&subsection)); 215 (*subdm)->nullspaceConstructors[0] = dm->nullspaceConstructors[0]; 216 if (dm->probs) { 217 PetscFV fv, fvNew; 218 PetscInt fnum[1] = {0}; 219 220 PetscCall(DMSetNumFields(*subdm, 1)); 221 PetscCall(DMGetField(dm, 0, NULL, (PetscObject *)&fv)); 222 PetscCall(PetscFVClone(fv, &fvNew)); 223 PetscCall(PetscFVSetNumComponents(fvNew, numFields)); 224 PetscCall(DMSetField(*subdm, 0, NULL, (PetscObject)fvNew)); 225 PetscCall(PetscFVDestroy(&fvNew)); 226 PetscCall(DMCreateDS(*subdm)); 227 if (numFields == 1 && is) { 228 PetscObject disc, space, pmat; 229 230 PetscCall(DMGetField(*subdm, 0, NULL, &disc)); 231 PetscCall(PetscObjectQuery(disc, "nullspace", &space)); 232 if (space) PetscCall(PetscObjectCompose((PetscObject)*is, "nullspace", space)); 233 PetscCall(PetscObjectQuery(disc, "nearnullspace", &space)); 234 if (space) PetscCall(PetscObjectCompose((PetscObject)*is, "nearnullspace", space)); 235 PetscCall(PetscObjectQuery(disc, "pmat", &pmat)); 236 if (pmat) PetscCall(PetscObjectCompose((PetscObject)*is, "pmat", pmat)); 237 } 238 PetscCall(PetscDSCopyConstants(dm->probs[0].ds, (*subdm)->probs[0].ds)); 239 PetscCall(PetscDSCopyBoundary(dm->probs[0].ds, 1, fnum, (*subdm)->probs[0].ds)); 240 PetscCall(PetscDSSelectEquations(dm->probs[0].ds, 1, fnum, (*subdm)->probs[0].ds)); 241 } 242 if ((*subdm)->nullspaceConstructors[0] && is) { 243 MatNullSpace nullSpace; 244 245 PetscCall((*(*subdm)->nullspaceConstructors[0])(*subdm, 0, 0, &nullSpace)); 246 PetscCall(PetscObjectCompose((PetscObject)*is, "nullspace", (PetscObject)nullSpace)); 247 PetscCall(MatNullSpaceDestroy(&nullSpace)); 248 } 249 if (dm->coarseMesh) PetscCall(DMCreateSubDM(dm->coarseMesh, numFields, fields, NULL, &(*subdm)->coarseMesh)); 250 PetscFunctionReturn(PETSC_SUCCESS); 251 } 252 PetscCall(PetscSectionCreateSubsection(section, numFields, fields, &subsection)); 253 PetscCall(DMSetLocalSection(*subdm, subsection)); 254 PetscCall(PetscSectionDestroy(&subsection)); 255 for (PetscInt f = 0; f < numFields; ++f) { 256 (*subdm)->nullspaceConstructors[f] = dm->nullspaceConstructors[fields[f]]; 257 if ((*subdm)->nullspaceConstructors[f]) { 258 haveNull = PETSC_TRUE; 259 nf = f; 260 of = fields[f]; 261 } 262 } 263 if (dm->probs) { 264 PetscCall(DMSetNumFields(*subdm, numFields)); 265 for (PetscInt f = 0; f < numFields; ++f) { 266 PetscObject disc; 267 268 PetscCall(DMGetField(dm, fields[f], NULL, &disc)); 269 PetscCall(DMSetField(*subdm, f, NULL, disc)); 270 } 271 // TODO: if only FV, then cut down the components 272 PetscCall(DMCreateDS(*subdm)); 273 if (numFields == 1 && is) { 274 PetscObject disc, space, pmat; 275 276 PetscCall(DMGetField(*subdm, 0, NULL, &disc)); 277 PetscCall(PetscObjectQuery(disc, "nullspace", &space)); 278 if (space) PetscCall(PetscObjectCompose((PetscObject)*is, "nullspace", space)); 279 PetscCall(PetscObjectQuery(disc, "nearnullspace", &space)); 280 if (space) PetscCall(PetscObjectCompose((PetscObject)*is, "nearnullspace", space)); 281 PetscCall(PetscObjectQuery(disc, "pmat", &pmat)); 282 if (pmat) PetscCall(PetscObjectCompose((PetscObject)*is, "pmat", pmat)); 283 } 284 // Check if DSes record their DM fields 285 if (dm->probs[0].fields) { 286 PetscInt d, e; 287 288 for (d = 0, e = 0; d < dm->Nds && e < (*subdm)->Nds; ++d) { 289 const PetscInt Nf = dm->probs[d].ds->Nf; 290 const PetscInt *fld; 291 PetscInt f, g; 292 293 PetscCall(ISGetIndices(dm->probs[d].fields, &fld)); 294 for (f = 0; f < Nf; ++f) { 295 for (g = 0; g < numFields; ++g) 296 if (fld[f] == fields[g]) break; 297 if (g < numFields) break; 298 } 299 PetscCall(ISRestoreIndices(dm->probs[d].fields, &fld)); 300 if (f == Nf) continue; 301 PetscCall(PetscDSCopyConstants(dm->probs[d].ds, (*subdm)->probs[e].ds)); 302 PetscCall(PetscDSCopyBoundary(dm->probs[d].ds, numFields, fields, (*subdm)->probs[e].ds)); 303 // Translate DM fields to DS fields 304 { 305 IS infields, dsfields; 306 const PetscInt *fld, *ofld; 307 PetscInt *fidx; 308 PetscInt onf, nf; 309 310 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, numFields, fields, PETSC_USE_POINTER, &infields)); 311 PetscCall(ISIntersect(infields, dm->probs[d].fields, &dsfields)); 312 PetscCall(ISDestroy(&infields)); 313 PetscCall(ISGetLocalSize(dsfields, &nf)); 314 PetscCheck(nf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "DS cannot be supported on 0 fields"); 315 PetscCall(ISGetIndices(dsfields, &fld)); 316 PetscCall(ISGetLocalSize(dm->probs[d].fields, &onf)); 317 PetscCall(ISGetIndices(dm->probs[d].fields, &ofld)); 318 PetscCall(PetscMalloc1(nf, &fidx)); 319 for (PetscInt f = 0, g = 0; f < onf && g < nf; ++f) { 320 if (ofld[f] == fld[g]) fidx[g++] = f; 321 } 322 PetscCall(ISRestoreIndices(dm->probs[d].fields, &ofld)); 323 PetscCall(ISRestoreIndices(dsfields, &fld)); 324 PetscCall(ISDestroy(&dsfields)); 325 PetscCall(PetscDSSelectDiscretizations(dm->probs[0].ds, nf, fidx, (*subdm)->probs[0].ds)); 326 PetscCall(PetscDSSelectEquations(dm->probs[0].ds, nf, fidx, (*subdm)->probs[0].ds)); 327 PetscCall(PetscFree(fidx)); 328 } 329 ++e; 330 } 331 } else { 332 PetscCall(PetscDSCopyConstants(dm->probs[0].ds, (*subdm)->probs[0].ds)); 333 PetscCall(PetscDSCopyBoundary(dm->probs[0].ds, PETSC_DETERMINE, NULL, (*subdm)->probs[0].ds)); 334 PetscCall(PetscDSSelectDiscretizations(dm->probs[0].ds, numFields, fields, (*subdm)->probs[0].ds)); 335 PetscCall(PetscDSSelectEquations(dm->probs[0].ds, numFields, fields, (*subdm)->probs[0].ds)); 336 } 337 } 338 if (haveNull && is) { 339 MatNullSpace nullSpace; 340 341 PetscCall((*(*subdm)->nullspaceConstructors[nf])(*subdm, of, nf, &nullSpace)); 342 PetscCall(PetscObjectCompose((PetscObject)*is, "nullspace", (PetscObject)nullSpace)); 343 PetscCall(MatNullSpaceDestroy(&nullSpace)); 344 } 345 if (dm->coarseMesh) PetscCall(DMCreateSubDM(dm->coarseMesh, numFields, fields, NULL, &(*subdm)->coarseMesh)); 346 PetscFunctionReturn(PETSC_SUCCESS); 347 } 348 349 /*@C 350 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`. 351 352 Not Collective 353 354 Input Parameters: 355 + dm - The `DM` object 356 . numFields - The number of fields to incorporate into `subdm` 357 - fields - The field numbers of the selected fields 358 359 Output Parameters: 360 + is - The global indices for the subproblem or `NULL` 361 - subdm - The `DM` for the subproblem, which must already have be cloned from `dm` or `NULL` 362 363 Level: intermediate 364 365 Notes: 366 If `is` and `subdm` are both `NULL` this does nothing 367 368 If there is only one field, and it is a `PetscFV`, then the selected fields will be interpreted as selected components. 369 370 .seealso: `DMCreateSubDM()`, `DMGetLocalSection()`, `DMPlexSetMigrationSF()`, `DMView()` 371 @*/ 372 PetscErrorCode DMCreateSectionSubDM(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm) 373 { 374 PetscSection section, sectionGlobal; 375 PetscInt Nf; 376 PetscBool useComp = PETSC_FALSE; 377 378 PetscFunctionBegin; 379 if (!numFields) PetscFunctionReturn(PETSC_SUCCESS); 380 PetscCall(DMGetLocalSection(dm, §ion)); 381 PetscCall(DMGetGlobalSection(dm, §ionGlobal)); 382 PetscCheck(section, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields"); 383 PetscCheck(sectionGlobal, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields"); 384 PetscCall(PetscSectionGetNumFields(section, &Nf)); 385 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); 386 387 // If there is only one field and it is a PetscFV, then interpret selection fields as components 388 { 389 PetscInt Nf; 390 391 PetscCall(DMGetNumFields(dm, &Nf)); 392 if (Nf == 1) { 393 PetscObject obj; 394 PetscClassId id; 395 396 PetscCall(DMGetField(dm, 0, NULL, &obj)); 397 PetscCall(PetscObjectGetClassId(obj, &id)); 398 if (id == PETSCFV_CLASSID) useComp = PETSC_TRUE; 399 } 400 } 401 if (is) PetscCall(PetscSectionSelectFields_Private(section, sectionGlobal, useComp, numFields, fields, is)); 402 if (subdm) PetscCall(DMSelectFields_Private(dm, section, useComp, numFields, fields, is, subdm)); 403 PetscFunctionReturn(PETSC_SUCCESS); 404 } 405 406 /*@C 407 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` 408 409 Not Collective 410 411 Input Parameters: 412 + dms - The `DM` objects, the must all have the same topology; for example obtained with `DMClone()` 413 - len - The number of `DM` in `dms` 414 415 Output Parameters: 416 + is - The global indices for the subproblem, or `NULL` 417 - superdm - The `DM` for the superproblem, which must already have be cloned and contain the same topology as the `dms` 418 419 Level: intermediate 420 421 .seealso: `DMCreateSuperDM()`, `DMGetLocalSection()`, `DMPlexSetMigrationSF()`, `DMView()` 422 @*/ 423 PetscErrorCode DMCreateSectionSuperDM(DM dms[], PetscInt len, IS **is, DM *superdm) 424 { 425 MPI_Comm comm; 426 PetscSection supersection, *sections, *sectionGlobals; 427 PetscInt *Nfs, Nf = 0, f, supf, oldf = -1, nullf = -1, i; 428 PetscBool haveNull = PETSC_FALSE; 429 430 PetscFunctionBegin; 431 PetscCall(PetscObjectGetComm((PetscObject)dms[0], &comm)); 432 /* Pull out local and global sections */ 433 PetscCall(PetscMalloc3(len, &Nfs, len, §ions, len, §ionGlobals)); 434 for (i = 0; i < len; ++i) { 435 PetscCall(DMGetLocalSection(dms[i], §ions[i])); 436 PetscCall(DMGetGlobalSection(dms[i], §ionGlobals[i])); 437 PetscCheck(sections[i], comm, PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields"); 438 PetscCheck(sectionGlobals[i], comm, PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields"); 439 PetscCall(PetscSectionGetNumFields(sections[i], &Nfs[i])); 440 Nf += Nfs[i]; 441 } 442 /* Create the supersection */ 443 PetscCall(PetscSectionCreateSupersection(sections, len, &supersection)); 444 PetscCall(DMSetLocalSection(*superdm, supersection)); 445 /* Create ISes */ 446 if (is) { 447 PetscSection supersectionGlobal; 448 PetscInt bs = -1, startf = 0; 449 450 PetscCall(PetscMalloc1(len, is)); 451 PetscCall(DMGetGlobalSection(*superdm, &supersectionGlobal)); 452 for (i = 0; i < len; startf += Nfs[i], ++i) { 453 PetscInt *subIndices; 454 PetscInt subSize, subOff, pStart, pEnd, p, start, end, dummy; 455 456 PetscCall(PetscSectionGetChart(sectionGlobals[i], &pStart, &pEnd)); 457 PetscCall(PetscSectionGetConstrainedStorageSize(sectionGlobals[i], &subSize)); 458 PetscCall(PetscMalloc1(subSize, &subIndices)); 459 for (p = pStart, subOff = 0; p < pEnd; ++p) { 460 PetscInt gdof, gcdof, gtdof, d; 461 462 PetscCall(PetscSectionGetDof(sectionGlobals[i], p, &gdof)); 463 PetscCall(PetscSectionGetConstraintDof(sections[i], p, &gcdof)); 464 gtdof = gdof - gcdof; 465 if (gdof > 0 && gtdof) { 466 if (bs < 0) { 467 bs = gtdof; 468 } else if (bs != gtdof) { 469 bs = 1; 470 } 471 PetscCall(DMGetGlobalFieldOffset_Private(*superdm, p, startf, &start, &dummy)); 472 PetscCall(DMGetGlobalFieldOffset_Private(*superdm, p, startf + Nfs[i] - 1, &dummy, &end)); 473 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); 474 for (d = start; d < end; ++d, ++subOff) subIndices[subOff] = d; 475 } 476 } 477 PetscCall(ISCreateGeneral(comm, subSize, subIndices, PETSC_OWN_POINTER, &(*is)[i])); 478 /* Must have same blocksize on all procs (some might have no points) */ 479 { 480 PetscInt bs = -1, bsLocal[2], bsMinMax[2]; 481 482 bsLocal[0] = bs < 0 ? PETSC_MAX_INT : bs; 483 bsLocal[1] = bs; 484 PetscCall(PetscGlobalMinMaxInt(comm, bsLocal, bsMinMax)); 485 if (bsMinMax[0] != bsMinMax[1]) { 486 bs = 1; 487 } else { 488 bs = bsMinMax[0]; 489 } 490 PetscCall(ISSetBlockSize((*is)[i], bs)); 491 } 492 } 493 } 494 /* Preserve discretizations */ 495 if (len && dms[0]->probs) { 496 PetscCall(DMSetNumFields(*superdm, Nf)); 497 for (i = 0, supf = 0; i < len; ++i) { 498 for (f = 0; f < Nfs[i]; ++f, ++supf) { 499 PetscObject disc; 500 501 PetscCall(DMGetField(dms[i], f, NULL, &disc)); 502 PetscCall(DMSetField(*superdm, supf, NULL, disc)); 503 } 504 } 505 PetscCall(DMCreateDS(*superdm)); 506 } 507 /* Preserve nullspaces */ 508 for (i = 0, supf = 0; i < len; ++i) { 509 for (f = 0; f < Nfs[i]; ++f, ++supf) { 510 (*superdm)->nullspaceConstructors[supf] = dms[i]->nullspaceConstructors[f]; 511 if ((*superdm)->nullspaceConstructors[supf]) { 512 haveNull = PETSC_TRUE; 513 nullf = supf; 514 oldf = f; 515 } 516 } 517 } 518 /* Attach nullspace to IS */ 519 if (haveNull && is) { 520 MatNullSpace nullSpace; 521 522 PetscCall((*(*superdm)->nullspaceConstructors[nullf])(*superdm, oldf, nullf, &nullSpace)); 523 PetscCall(PetscObjectCompose((PetscObject)(*is)[nullf], "nullspace", (PetscObject)nullSpace)); 524 PetscCall(MatNullSpaceDestroy(&nullSpace)); 525 } 526 PetscCall(PetscSectionDestroy(&supersection)); 527 PetscCall(PetscFree3(Nfs, sections, sectionGlobals)); 528 PetscFunctionReturn(PETSC_SUCCESS); 529 } 530