1af0996ceSBarry Smith #include <petsc/private/dmimpl.h> /*I "petscdm.h" I*/ 22764a2aaSMatthew G. Knepley #include <petscds.h> 311689aebSJed Brown 45b8243e1SJed Brown // Greatest common divisor of two nonnegative integers 5d71ae5a4SJacob Faibussowitsch PetscInt PetscGCD(PetscInt a, PetscInt b) 6d71ae5a4SJacob Faibussowitsch { 75b8243e1SJed Brown while (b != 0) { 85b8243e1SJed Brown PetscInt tmp = b; 95b8243e1SJed Brown b = a % b; 105b8243e1SJed Brown a = tmp; 115b8243e1SJed Brown } 125b8243e1SJed Brown return a; 135b8243e1SJed Brown } 145b8243e1SJed Brown 15d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateGlobalVector_Section_Private(DM dm, Vec *vec) 16d71ae5a4SJacob Faibussowitsch { 1711689aebSJed Brown PetscSection gSection; 18cc30b04dSMatthew G Knepley PetscInt localSize, bs, blockSize = -1, pStart, pEnd, p; 195d1c0279SHong Zhang PetscInt in[2], out[2]; 2011689aebSJed Brown 2111689aebSJed Brown PetscFunctionBegin; 229566063dSJacob Faibussowitsch PetscCall(DMGetGlobalSection(dm, &gSection)); 239566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(gSection, &pStart, &pEnd)); 2411689aebSJed Brown for (p = pStart; p < pEnd; ++p) { 2511689aebSJed Brown PetscInt dof, cdof; 2611689aebSJed Brown 279566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(gSection, p, &dof)); 289566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(gSection, p, &cdof)); 290680ce0aSHong Zhang 305b8243e1SJed Brown if (dof - cdof > 0) { 315b8243e1SJed Brown if (blockSize < 0) { 320680ce0aSHong Zhang /* set blockSize */ 330680ce0aSHong Zhang blockSize = dof - cdof; 345b8243e1SJed Brown } else { 355b8243e1SJed Brown blockSize = PetscGCD(dof - cdof, blockSize); 3611689aebSJed Brown } 3711689aebSJed Brown } 380680ce0aSHong Zhang } 395d1c0279SHong Zhang 402c8be454SMatthew G. Knepley in[0] = blockSize < 0 ? PETSC_MIN_INT : -blockSize; 415d1c0279SHong Zhang in[1] = blockSize; 421c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(in, out, 2, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm))); 435d1c0279SHong Zhang /* -out[0] = min(blockSize), out[1] = max(blockSize) */ 445d1c0279SHong Zhang if (-out[0] == out[1]) { 455d1c0279SHong Zhang bs = out[1]; 465d1c0279SHong Zhang } else bs = 1; 475d1c0279SHong Zhang 485d1c0279SHong Zhang if (bs < 0) { /* Everyone was empty */ 495d1c0279SHong Zhang blockSize = 1; 505d1c0279SHong Zhang bs = 1; 515d1c0279SHong Zhang } 525d1c0279SHong Zhang 539566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstrainedStorageSize(gSection, &localSize)); 547a8be351SBarry Smith PetscCheck(localSize % blockSize == 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mismatch between blocksize %" PetscInt_FMT " and local storage size %" PetscInt_FMT, blockSize, localSize); 559566063dSJacob Faibussowitsch PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), vec)); 569566063dSJacob Faibussowitsch PetscCall(VecSetSizes(*vec, localSize, PETSC_DETERMINE)); 579566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(*vec, bs)); 589566063dSJacob Faibussowitsch PetscCall(VecSetType(*vec, dm->vectype)); 599566063dSJacob Faibussowitsch PetscCall(VecSetDM(*vec, dm)); 609566063dSJacob Faibussowitsch /* PetscCall(VecSetLocalToGlobalMapping(*vec, dm->ltogmap)); */ 613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6211689aebSJed Brown } 6311689aebSJed Brown 64d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateLocalVector_Section_Private(DM dm, Vec *vec) 65d71ae5a4SJacob Faibussowitsch { 66fad22124SMatthew G Knepley PetscSection section; 6711689aebSJed Brown PetscInt localSize, blockSize = -1, pStart, pEnd, p; 6811689aebSJed Brown 6911689aebSJed Brown PetscFunctionBegin; 709566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, §ion)); 719566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(section, &pStart, &pEnd)); 7211689aebSJed Brown for (p = pStart; p < pEnd; ++p) { 7311689aebSJed Brown PetscInt dof; 7411689aebSJed Brown 759566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, p, &dof)); 7611689aebSJed Brown if ((blockSize < 0) && (dof > 0)) blockSize = dof; 775b8243e1SJed Brown if (dof > 0) blockSize = PetscGCD(dof, blockSize); 7811689aebSJed Brown } 799566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(section, &localSize)); 809566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, vec)); 819566063dSJacob Faibussowitsch PetscCall(VecSetSizes(*vec, localSize, localSize)); 829566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(*vec, blockSize)); 839566063dSJacob Faibussowitsch PetscCall(VecSetType(*vec, dm->vectype)); 849566063dSJacob Faibussowitsch PetscCall(VecSetDM(*vec, dm)); 853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8611689aebSJed Brown } 874d9407bcSMatthew G. Knepley 88792b654fSMatthew G. Knepley /*@C 89*0b3275a6SBarry Smith 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`. 90792b654fSMatthew G. Knepley 9120f4b53cSBarry Smith Not Collective 92792b654fSMatthew G. Knepley 93792b654fSMatthew G. Knepley Input Parameters: 9420f4b53cSBarry Smith + dm - The `DM` object 95*0b3275a6SBarry Smith . numFields - The number of fields to incorporate into `subdm` 96792b654fSMatthew G. Knepley - fields - The field numbers of the selected fields 97792b654fSMatthew G. Knepley 98792b654fSMatthew G. Knepley Output Parameters: 99*0b3275a6SBarry Smith + is - The global indices for the subproblem or `NULL` 100*0b3275a6SBarry Smith - subdm - The `DM` for the subproblem, which must already have be cloned from `dm` or `NULL` 101792b654fSMatthew G. Knepley 102792b654fSMatthew G. Knepley Level: intermediate 103792b654fSMatthew G. Knepley 10420f4b53cSBarry Smith Note: 105*0b3275a6SBarry Smith If `is` and `subdm` are both `NULL` this does nothing 10620f4b53cSBarry Smith 107a4e35b19SJacob Faibussowitsch .seealso: `DMCreateSubDM()`, `DMGetLocalSection()`, `DMPlexSetMigrationSF()`, `DMView()` 108792b654fSMatthew G. Knepley @*/ 109d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateSectionSubDM(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm) 110d71ae5a4SJacob Faibussowitsch { 1114d9407bcSMatthew G. Knepley PetscSection section, sectionGlobal; 1124d9407bcSMatthew G. Knepley PetscInt *subIndices; 113696188eaSMatthew G. Knepley PetscInt subSize = 0, subOff = 0, Nf, f, pStart, pEnd, p; 1144d9407bcSMatthew G. Knepley 1154d9407bcSMatthew G. Knepley PetscFunctionBegin; 1163ba16761SJacob Faibussowitsch if (!numFields) PetscFunctionReturn(PETSC_SUCCESS); 1179566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, §ion)); 1189566063dSJacob Faibussowitsch PetscCall(DMGetGlobalSection(dm, §ionGlobal)); 1197a8be351SBarry Smith PetscCheck(section, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields"); 1207a8be351SBarry Smith PetscCheck(sectionGlobal, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields"); 1219566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section, &Nf)); 1227a8be351SBarry Smith 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); 1234d9407bcSMatthew G. Knepley if (is) { 1243a544194SStefano Zampini PetscInt bs, bsLocal[2], bsMinMax[2]; 1253a544194SStefano Zampini 1263a544194SStefano Zampini for (f = 0, bs = 0; f < numFields; ++f) { 1273a544194SStefano Zampini PetscInt Nc; 1283a544194SStefano Zampini 1299566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldComponents(section, fields[f], &Nc)); 1303a544194SStefano Zampini bs += Nc; 1313a544194SStefano Zampini } 1329566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(sectionGlobal, &pStart, &pEnd)); 1334d9407bcSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 134b9eca265Ssarens PetscInt gdof, pSubSize = 0; 1354d9407bcSMatthew G. Knepley 1369566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(sectionGlobal, p, &gdof)); 1374d9407bcSMatthew G. Knepley if (gdof > 0) { 1384d9407bcSMatthew G. Knepley for (f = 0; f < numFields; ++f) { 1394d9407bcSMatthew G. Knepley PetscInt fdof, fcdof; 1404d9407bcSMatthew G. Knepley 1419566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(section, p, fields[f], &fdof)); 1429566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldConstraintDof(section, p, fields[f], &fcdof)); 143b9eca265Ssarens pSubSize += fdof - fcdof; 144b9eca265Ssarens } 145b9eca265Ssarens subSize += pSubSize; 1463a544194SStefano Zampini if (pSubSize && bs != pSubSize) { 147b9eca265Ssarens /* Layout does not admit a pointwise block size */ 148b9eca265Ssarens bs = 1; 1494d9407bcSMatthew G. Knepley } 1504d9407bcSMatthew G. Knepley } 1514d9407bcSMatthew G. Knepley } 152b9eca265Ssarens /* Must have same blocksize on all procs (some might have no points) */ 1539371c9d4SSatish Balay bsLocal[0] = bs < 0 ? PETSC_MAX_INT : bs; 1549371c9d4SSatish Balay bsLocal[1] = bs; 1559566063dSJacob Faibussowitsch PetscCall(PetscGlobalMinMaxInt(PetscObjectComm((PetscObject)dm), bsLocal, bsMinMax)); 1569371c9d4SSatish Balay if (bsMinMax[0] != bsMinMax[1]) { 1579371c9d4SSatish Balay bs = 1; 1589371c9d4SSatish Balay } else { 1599371c9d4SSatish Balay bs = bsMinMax[0]; 1609371c9d4SSatish Balay } 1619566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(subSize, &subIndices)); 1624d9407bcSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 1634d9407bcSMatthew G. Knepley PetscInt gdof, goff; 1644d9407bcSMatthew G. Knepley 1659566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(sectionGlobal, p, &gdof)); 1664d9407bcSMatthew G. Knepley if (gdof > 0) { 1679566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(sectionGlobal, p, &goff)); 1684d9407bcSMatthew G. Knepley for (f = 0; f < numFields; ++f) { 1694d9407bcSMatthew G. Knepley PetscInt fdof, fcdof, fc, f2, poff = 0; 1704d9407bcSMatthew G. Knepley 1714d9407bcSMatthew G. Knepley /* Can get rid of this loop by storing field information in the global section */ 1724d9407bcSMatthew G. Knepley for (f2 = 0; f2 < fields[f]; ++f2) { 1739566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(section, p, f2, &fdof)); 1749566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldConstraintDof(section, p, f2, &fcdof)); 1754d9407bcSMatthew G. Knepley poff += fdof - fcdof; 1764d9407bcSMatthew G. Knepley } 1779566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(section, p, fields[f], &fdof)); 1789566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldConstraintDof(section, p, fields[f], &fcdof)); 179ad540459SPierre Jolivet for (fc = 0; fc < fdof - fcdof; ++fc, ++subOff) subIndices[subOff] = goff + poff + fc; 1804d9407bcSMatthew G. Knepley } 1814d9407bcSMatthew G. Knepley } 1824d9407bcSMatthew G. Knepley } 1839566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), subSize, subIndices, PETSC_OWN_POINTER, is)); 184627349f5SMatthew G. Knepley if (bs > 1) { 185627349f5SMatthew G. Knepley /* We need to check that the block size does not come from non-contiguous fields */ 18688206212SAlexis Marboeuf PetscInt i, j, set = 1, rset = 1; 187627349f5SMatthew G. Knepley for (i = 0; i < subSize; i += bs) { 188627349f5SMatthew G. Knepley for (j = 0; j < bs; ++j) { 1899371c9d4SSatish Balay if (subIndices[i + j] != subIndices[i] + j) { 1909371c9d4SSatish Balay set = 0; 1919371c9d4SSatish Balay break; 1929371c9d4SSatish Balay } 193627349f5SMatthew G. Knepley } 194627349f5SMatthew G. Knepley } 195712fec58SPierre Jolivet PetscCall(MPIU_Allreduce(&set, &rset, 1, MPIU_INT, MPI_PROD, PetscObjectComm((PetscObject)dm))); 19688206212SAlexis Marboeuf if (rset) PetscCall(ISSetBlockSize(*is, bs)); 197627349f5SMatthew G. Knepley } 1984d9407bcSMatthew G. Knepley } 1994d9407bcSMatthew G. Knepley if (subdm) { 2004d9407bcSMatthew G. Knepley PetscSection subsection; 2014d9407bcSMatthew G. Knepley PetscBool haveNull = PETSC_FALSE; 2028cda7954SMatthew G. Knepley PetscInt f, nf = 0, of = 0; 2034d9407bcSMatthew G. Knepley 2049566063dSJacob Faibussowitsch PetscCall(PetscSectionCreateSubsection(section, numFields, fields, &subsection)); 2059566063dSJacob Faibussowitsch PetscCall(DMSetLocalSection(*subdm, subsection)); 2069566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&subsection)); 2074d9407bcSMatthew G. Knepley for (f = 0; f < numFields; ++f) { 2084d9407bcSMatthew G. Knepley (*subdm)->nullspaceConstructors[f] = dm->nullspaceConstructors[fields[f]]; 2094d9407bcSMatthew G. Knepley if ((*subdm)->nullspaceConstructors[f]) { 2104d9407bcSMatthew G. Knepley haveNull = PETSC_TRUE; 2114d9407bcSMatthew G. Knepley nf = f; 2128cda7954SMatthew G. Knepley of = fields[f]; 2134d9407bcSMatthew G. Knepley } 2144d9407bcSMatthew G. Knepley } 215e5e52638SMatthew G. Knepley if (dm->probs) { 2169566063dSJacob Faibussowitsch PetscCall(DMSetNumFields(*subdm, numFields)); 2174d9407bcSMatthew G. Knepley for (f = 0; f < numFields; ++f) { 2180f21e855SMatthew G. Knepley PetscObject disc; 2190f21e855SMatthew G. Knepley 2209566063dSJacob Faibussowitsch PetscCall(DMGetField(dm, fields[f], NULL, &disc)); 2219566063dSJacob Faibussowitsch PetscCall(DMSetField(*subdm, f, NULL, disc)); 2224d9407bcSMatthew G. Knepley } 2239566063dSJacob Faibussowitsch PetscCall(DMCreateDS(*subdm)); 224f646a522SMatthew G. Knepley if (numFields == 1 && is) { 2250f21e855SMatthew G. Knepley PetscObject disc, space, pmat; 2264d9407bcSMatthew G. Knepley 2279566063dSJacob Faibussowitsch PetscCall(DMGetField(*subdm, 0, NULL, &disc)); 2289566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(disc, "nullspace", &space)); 2299566063dSJacob Faibussowitsch if (space) PetscCall(PetscObjectCompose((PetscObject)*is, "nullspace", space)); 2309566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(disc, "nearnullspace", &space)); 2319566063dSJacob Faibussowitsch if (space) PetscCall(PetscObjectCompose((PetscObject)*is, "nearnullspace", space)); 2329566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(disc, "pmat", &pmat)); 2339566063dSJacob Faibussowitsch if (pmat) PetscCall(PetscObjectCompose((PetscObject)*is, "pmat", pmat)); 2344d9407bcSMatthew G. Knepley } 2359310035eSMatthew G. Knepley /* Check if DSes record their DM fields */ 236e21d936aSMatthew G. Knepley if (dm->probs[0].fields) { 2379310035eSMatthew G. Knepley PetscInt d, e; 2389310035eSMatthew G. Knepley 2399310035eSMatthew G. Knepley for (d = 0, e = 0; d < dm->Nds && e < (*subdm)->Nds; ++d) { 2409310035eSMatthew G. Knepley const PetscInt Nf = dm->probs[d].ds->Nf; 2419310035eSMatthew G. Knepley const PetscInt *fld; 2429310035eSMatthew G. Knepley PetscInt f, g; 2439310035eSMatthew G. Knepley 2449566063dSJacob Faibussowitsch PetscCall(ISGetIndices(dm->probs[d].fields, &fld)); 2459310035eSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 2469371c9d4SSatish Balay for (g = 0; g < numFields; ++g) 2479371c9d4SSatish Balay if (fld[f] == fields[g]) break; 2489310035eSMatthew G. Knepley if (g < numFields) break; 2499310035eSMatthew G. Knepley } 2509566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(dm->probs[d].fields, &fld)); 2519310035eSMatthew G. Knepley if (f == Nf) continue; 2529566063dSJacob Faibussowitsch PetscCall(PetscDSCopyConstants(dm->probs[d].ds, (*subdm)->probs[e].ds)); 2539566063dSJacob Faibussowitsch PetscCall(PetscDSCopyBoundary(dm->probs[d].ds, numFields, fields, (*subdm)->probs[e].ds)); 2549310035eSMatthew G. Knepley /* Translate DM fields to DS fields */ 2559310035eSMatthew G. Knepley { 256e21d936aSMatthew G. Knepley IS infields, dsfields; 25774a61aaaSMatthew G. Knepley const PetscInt *fld, *ofld; 25874a61aaaSMatthew G. Knepley PetscInt *fidx; 25974a61aaaSMatthew G. Knepley PetscInt onf, nf, f, g; 260e21d936aSMatthew G. Knepley 2619566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, numFields, fields, PETSC_USE_POINTER, &infields)); 2629566063dSJacob Faibussowitsch PetscCall(ISIntersect(infields, dm->probs[d].fields, &dsfields)); 2639566063dSJacob Faibussowitsch PetscCall(ISDestroy(&infields)); 2649566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(dsfields, &nf)); 2657a8be351SBarry Smith PetscCheck(nf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "DS cannot be supported on 0 fields"); 2669566063dSJacob Faibussowitsch PetscCall(ISGetIndices(dsfields, &fld)); 2679566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(dm->probs[d].fields, &onf)); 2689566063dSJacob Faibussowitsch PetscCall(ISGetIndices(dm->probs[d].fields, &ofld)); 2699566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nf, &fidx)); 2703268a0aaSMatthew G. Knepley for (f = 0, g = 0; f < onf && g < nf; ++f) { 27174a61aaaSMatthew G. Knepley if (ofld[f] == fld[g]) fidx[g++] = f; 27274a61aaaSMatthew G. Knepley } 2739566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(dm->probs[d].fields, &ofld)); 2749566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(dsfields, &fld)); 2759566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dsfields)); 2769566063dSJacob Faibussowitsch PetscCall(PetscDSSelectDiscretizations(dm->probs[0].ds, nf, fidx, (*subdm)->probs[0].ds)); 2779566063dSJacob Faibussowitsch PetscCall(PetscDSSelectEquations(dm->probs[0].ds, nf, fidx, (*subdm)->probs[0].ds)); 2789566063dSJacob Faibussowitsch PetscCall(PetscFree(fidx)); 2799310035eSMatthew G. Knepley } 2809310035eSMatthew G. Knepley ++e; 2819310035eSMatthew G. Knepley } 282e21d936aSMatthew G. Knepley } else { 2839566063dSJacob Faibussowitsch PetscCall(PetscDSCopyConstants(dm->probs[0].ds, (*subdm)->probs[0].ds)); 2849566063dSJacob Faibussowitsch PetscCall(PetscDSCopyBoundary(dm->probs[0].ds, PETSC_DETERMINE, NULL, (*subdm)->probs[0].ds)); 2859566063dSJacob Faibussowitsch PetscCall(PetscDSSelectDiscretizations(dm->probs[0].ds, numFields, fields, (*subdm)->probs[0].ds)); 2869566063dSJacob Faibussowitsch PetscCall(PetscDSSelectEquations(dm->probs[0].ds, numFields, fields, (*subdm)->probs[0].ds)); 287d5af271aSMatthew G. Knepley } 288e21d936aSMatthew G. Knepley } 2898cda7954SMatthew G. Knepley if (haveNull && is) { 2908cda7954SMatthew G. Knepley MatNullSpace nullSpace; 2918cda7954SMatthew G. Knepley 2929566063dSJacob Faibussowitsch PetscCall((*(*subdm)->nullspaceConstructors[nf])(*subdm, of, nf, &nullSpace)); 2939566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)*is, "nullspace", (PetscObject)nullSpace)); 2949566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&nullSpace)); 2958cda7954SMatthew G. Knepley } 29648a46eb9SPierre Jolivet if (dm->coarseMesh) PetscCall(DMCreateSubDM(dm->coarseMesh, numFields, fields, NULL, &(*subdm)->coarseMesh)); 2974d9407bcSMatthew G. Knepley } 2983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2994d9407bcSMatthew G. Knepley } 3002adcc780SMatthew G. Knepley 301792b654fSMatthew G. Knepley /*@C 302*0b3275a6SBarry Smith 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` 303792b654fSMatthew G. Knepley 30420f4b53cSBarry Smith Not Collective 305792b654fSMatthew G. Knepley 306d8d19677SJose E. Roman Input Parameters: 307*0b3275a6SBarry Smith + dms - The `DM` objects, the must all have the same topology; for example obtained with `DMClone()` 308*0b3275a6SBarry Smith - len - The number of `DM` in `dms` 309792b654fSMatthew G. Knepley 310792b654fSMatthew G. Knepley Output Parameters: 31120f4b53cSBarry Smith + is - The global indices for the subproblem, or `NULL` 312*0b3275a6SBarry Smith - superdm - The `DM` for the superproblem, which must already have be cloned and contain the same topology as the `dms` 313792b654fSMatthew G. Knepley 314792b654fSMatthew G. Knepley Level: intermediate 315792b654fSMatthew G. Knepley 316a4e35b19SJacob Faibussowitsch .seealso: `DMCreateSuperDM()`, `DMGetLocalSection()`, `DMPlexSetMigrationSF()`, `DMView()` 317792b654fSMatthew G. Knepley @*/ 318d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateSectionSuperDM(DM dms[], PetscInt len, IS **is, DM *superdm) 319d71ae5a4SJacob Faibussowitsch { 3209796d432SMatthew G. Knepley MPI_Comm comm; 3219796d432SMatthew G. Knepley PetscSection supersection, *sections, *sectionGlobals; 3228cda7954SMatthew G. Knepley PetscInt *Nfs, Nf = 0, f, supf, oldf = -1, nullf = -1, i; 3239796d432SMatthew G. Knepley PetscBool haveNull = PETSC_FALSE; 3242adcc780SMatthew G. Knepley 3252adcc780SMatthew G. Knepley PetscFunctionBegin; 3269566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dms[0], &comm)); 3279796d432SMatthew G. Knepley /* Pull out local and global sections */ 3289566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(len, &Nfs, len, §ions, len, §ionGlobals)); 3292adcc780SMatthew G. Knepley for (i = 0; i < len; ++i) { 3309566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dms[i], §ions[i])); 3319566063dSJacob Faibussowitsch PetscCall(DMGetGlobalSection(dms[i], §ionGlobals[i])); 3327a8be351SBarry Smith PetscCheck(sections[i], comm, PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields"); 3337a8be351SBarry Smith PetscCheck(sectionGlobals[i], comm, PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields"); 3349566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(sections[i], &Nfs[i])); 3352adcc780SMatthew G. Knepley Nf += Nfs[i]; 3362adcc780SMatthew G. Knepley } 3379796d432SMatthew G. Knepley /* Create the supersection */ 3389566063dSJacob Faibussowitsch PetscCall(PetscSectionCreateSupersection(sections, len, &supersection)); 3399566063dSJacob Faibussowitsch PetscCall(DMSetLocalSection(*superdm, supersection)); 3409796d432SMatthew G. Knepley /* Create ISes */ 3412adcc780SMatthew G. Knepley if (is) { 3429796d432SMatthew G. Knepley PetscSection supersectionGlobal; 3439796d432SMatthew G. Knepley PetscInt bs = -1, startf = 0; 3442adcc780SMatthew G. Knepley 3459566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(len, is)); 3469566063dSJacob Faibussowitsch PetscCall(DMGetGlobalSection(*superdm, &supersectionGlobal)); 3479796d432SMatthew G. Knepley for (i = 0; i < len; startf += Nfs[i], ++i) { 3489796d432SMatthew G. Knepley PetscInt *subIndices; 349ec4c761aSMatthew G. Knepley PetscInt subSize, subOff, pStart, pEnd, p, start, end, dummy; 3502adcc780SMatthew G. Knepley 3519566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(sectionGlobals[i], &pStart, &pEnd)); 3529566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstrainedStorageSize(sectionGlobals[i], &subSize)); 3539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(subSize, &subIndices)); 3549796d432SMatthew G. Knepley for (p = pStart, subOff = 0; p < pEnd; ++p) { 3559796d432SMatthew G. Knepley PetscInt gdof, gcdof, gtdof, d; 3562adcc780SMatthew G. Knepley 3579566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(sectionGlobals[i], p, &gdof)); 3589566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(sections[i], p, &gcdof)); 3592adcc780SMatthew G. Knepley gtdof = gdof - gcdof; 3602adcc780SMatthew G. Knepley if (gdof > 0 && gtdof) { 3619371c9d4SSatish Balay if (bs < 0) { 3629371c9d4SSatish Balay bs = gtdof; 3639371c9d4SSatish Balay } else if (bs != gtdof) { 3649371c9d4SSatish Balay bs = 1; 3659371c9d4SSatish Balay } 3669566063dSJacob Faibussowitsch PetscCall(DMGetGlobalFieldOffset_Private(*superdm, p, startf, &start, &dummy)); 3679566063dSJacob Faibussowitsch PetscCall(DMGetGlobalFieldOffset_Private(*superdm, p, startf + Nfs[i] - 1, &dummy, &end)); 36863a3b9bcSJacob Faibussowitsch 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); 3699796d432SMatthew G. Knepley for (d = start; d < end; ++d, ++subOff) subIndices[subOff] = d; 3702adcc780SMatthew G. Knepley } 3712adcc780SMatthew G. Knepley } 3729566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, subSize, subIndices, PETSC_OWN_POINTER, &(*is)[i])); 3732adcc780SMatthew G. Knepley /* Must have same blocksize on all procs (some might have no points) */ 3749796d432SMatthew G. Knepley { 3759796d432SMatthew G. Knepley PetscInt bs = -1, bsLocal[2], bsMinMax[2]; 3769796d432SMatthew G. Knepley 3779371c9d4SSatish Balay bsLocal[0] = bs < 0 ? PETSC_MAX_INT : bs; 3789371c9d4SSatish Balay bsLocal[1] = bs; 3799566063dSJacob Faibussowitsch PetscCall(PetscGlobalMinMaxInt(comm, bsLocal, bsMinMax)); 3809371c9d4SSatish Balay if (bsMinMax[0] != bsMinMax[1]) { 3819371c9d4SSatish Balay bs = 1; 3829371c9d4SSatish Balay } else { 3839371c9d4SSatish Balay bs = bsMinMax[0]; 3849371c9d4SSatish Balay } 3859566063dSJacob Faibussowitsch PetscCall(ISSetBlockSize((*is)[i], bs)); 3862adcc780SMatthew G. Knepley } 3872adcc780SMatthew G. Knepley } 3882adcc780SMatthew G. Knepley } 3899796d432SMatthew G. Knepley /* Preserve discretizations */ 390e5e52638SMatthew G. Knepley if (len && dms[0]->probs) { 3919566063dSJacob Faibussowitsch PetscCall(DMSetNumFields(*superdm, Nf)); 3929796d432SMatthew G. Knepley for (i = 0, supf = 0; i < len; ++i) { 3939796d432SMatthew G. Knepley for (f = 0; f < Nfs[i]; ++f, ++supf) { 3942adcc780SMatthew G. Knepley PetscObject disc; 3952adcc780SMatthew G. Knepley 3969566063dSJacob Faibussowitsch PetscCall(DMGetField(dms[i], f, NULL, &disc)); 3979566063dSJacob Faibussowitsch PetscCall(DMSetField(*superdm, supf, NULL, disc)); 3982adcc780SMatthew G. Knepley } 3992adcc780SMatthew G. Knepley } 4009566063dSJacob Faibussowitsch PetscCall(DMCreateDS(*superdm)); 4012adcc780SMatthew G. Knepley } 4029796d432SMatthew G. Knepley /* Preserve nullspaces */ 4039796d432SMatthew G. Knepley for (i = 0, supf = 0; i < len; ++i) { 4049796d432SMatthew G. Knepley for (f = 0; f < Nfs[i]; ++f, ++supf) { 4059796d432SMatthew G. Knepley (*superdm)->nullspaceConstructors[supf] = dms[i]->nullspaceConstructors[f]; 4069796d432SMatthew G. Knepley if ((*superdm)->nullspaceConstructors[supf]) { 4079796d432SMatthew G. Knepley haveNull = PETSC_TRUE; 4089796d432SMatthew G. Knepley nullf = supf; 4098cda7954SMatthew G. Knepley oldf = f; 4102adcc780SMatthew G. Knepley } 4119796d432SMatthew G. Knepley } 4129796d432SMatthew G. Knepley } 4139796d432SMatthew G. Knepley /* Attach nullspace to IS */ 4149796d432SMatthew G. Knepley if (haveNull && is) { 4159796d432SMatthew G. Knepley MatNullSpace nullSpace; 4169796d432SMatthew G. Knepley 4179566063dSJacob Faibussowitsch PetscCall((*(*superdm)->nullspaceConstructors[nullf])(*superdm, oldf, nullf, &nullSpace)); 4189566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)(*is)[nullf], "nullspace", (PetscObject)nullSpace)); 4199566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&nullSpace)); 4209796d432SMatthew G. Knepley } 4219566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&supersection)); 4229566063dSJacob Faibussowitsch PetscCall(PetscFree3(Nfs, sections, sectionGlobals)); 4233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4242adcc780SMatthew G. Knepley } 425