xref: /petsc/src/dm/interface/dmi.c (revision a4e35b1925eceef64945ea472b84f2bf06a67b5e)
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, &section));
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
8920f4b53cSBarry Smith   DMCreateSectionSubDM - Returns an `IS` and subDM+subSection encapsulating a subproblem defined by 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
95792b654fSMatthew G. Knepley . numFields - The number of fields in this subproblem
96792b654fSMatthew G. Knepley - fields    - The field numbers of the selected fields
97792b654fSMatthew G. Knepley 
98792b654fSMatthew G. Knepley   Output Parameters:
99792b654fSMatthew G. Knepley + is    - The global indices for the subproblem
10020f4b53cSBarry Smith - subdm - The `DM` for the subproblem, which must already have be cloned from `dm`
101792b654fSMatthew G. Knepley 
102792b654fSMatthew G. Knepley   Level: intermediate
103792b654fSMatthew G. Knepley 
10420f4b53cSBarry Smith   Note:
10520f4b53cSBarry Smith   This handles all information in the `DM` class and the `PetscSection`. This is used as the basis for creating subDMs in specialized classes,
10620f4b53cSBarry Smith   such as `DMPLEX` and `DMFOREST`
10720f4b53cSBarry Smith 
108*a4e35b19SJacob Faibussowitsch .seealso: `DMCreateSubDM()`, `DMGetLocalSection()`, `DMPlexSetMigrationSF()`, `DMView()`
109792b654fSMatthew G. Knepley @*/
110d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateSectionSubDM(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm)
111d71ae5a4SJacob Faibussowitsch {
1124d9407bcSMatthew G. Knepley   PetscSection section, sectionGlobal;
1134d9407bcSMatthew G. Knepley   PetscInt    *subIndices;
114696188eaSMatthew G. Knepley   PetscInt     subSize = 0, subOff = 0, Nf, f, pStart, pEnd, p;
1154d9407bcSMatthew G. Knepley 
1164d9407bcSMatthew G. Knepley   PetscFunctionBegin;
1173ba16761SJacob Faibussowitsch   if (!numFields) PetscFunctionReturn(PETSC_SUCCESS);
1189566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
1199566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalSection(dm, &sectionGlobal));
1207a8be351SBarry Smith   PetscCheck(section, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields");
1217a8be351SBarry Smith   PetscCheck(sectionGlobal, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields");
1229566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetNumFields(section, &Nf));
1237a8be351SBarry 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);
1244d9407bcSMatthew G. Knepley   if (is) {
1253a544194SStefano Zampini     PetscInt bs, bsLocal[2], bsMinMax[2];
1263a544194SStefano Zampini 
1273a544194SStefano Zampini     for (f = 0, bs = 0; f < numFields; ++f) {
1283a544194SStefano Zampini       PetscInt Nc;
1293a544194SStefano Zampini 
1309566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldComponents(section, fields[f], &Nc));
1313a544194SStefano Zampini       bs += Nc;
1323a544194SStefano Zampini     }
1339566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetChart(sectionGlobal, &pStart, &pEnd));
1344d9407bcSMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
135b9eca265Ssarens       PetscInt gdof, pSubSize = 0;
1364d9407bcSMatthew G. Knepley 
1379566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(sectionGlobal, p, &gdof));
1384d9407bcSMatthew G. Knepley       if (gdof > 0) {
1394d9407bcSMatthew G. Knepley         for (f = 0; f < numFields; ++f) {
1404d9407bcSMatthew G. Knepley           PetscInt fdof, fcdof;
1414d9407bcSMatthew G. Knepley 
1429566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetFieldDof(section, p, fields[f], &fdof));
1439566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetFieldConstraintDof(section, p, fields[f], &fcdof));
144b9eca265Ssarens           pSubSize += fdof - fcdof;
145b9eca265Ssarens         }
146b9eca265Ssarens         subSize += pSubSize;
1473a544194SStefano Zampini         if (pSubSize && bs != pSubSize) {
148b9eca265Ssarens           /* Layout does not admit a pointwise block size */
149b9eca265Ssarens           bs = 1;
1504d9407bcSMatthew G. Knepley         }
1514d9407bcSMatthew G. Knepley       }
1524d9407bcSMatthew G. Knepley     }
153b9eca265Ssarens     /* Must have same blocksize on all procs (some might have no points) */
1549371c9d4SSatish Balay     bsLocal[0] = bs < 0 ? PETSC_MAX_INT : bs;
1559371c9d4SSatish Balay     bsLocal[1] = bs;
1569566063dSJacob Faibussowitsch     PetscCall(PetscGlobalMinMaxInt(PetscObjectComm((PetscObject)dm), bsLocal, bsMinMax));
1579371c9d4SSatish Balay     if (bsMinMax[0] != bsMinMax[1]) {
1589371c9d4SSatish Balay       bs = 1;
1599371c9d4SSatish Balay     } else {
1609371c9d4SSatish Balay       bs = bsMinMax[0];
1619371c9d4SSatish Balay     }
1629566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(subSize, &subIndices));
1634d9407bcSMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
1644d9407bcSMatthew G. Knepley       PetscInt gdof, goff;
1654d9407bcSMatthew G. Knepley 
1669566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(sectionGlobal, p, &gdof));
1674d9407bcSMatthew G. Knepley       if (gdof > 0) {
1689566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetOffset(sectionGlobal, p, &goff));
1694d9407bcSMatthew G. Knepley         for (f = 0; f < numFields; ++f) {
1704d9407bcSMatthew G. Knepley           PetscInt fdof, fcdof, fc, f2, poff = 0;
1714d9407bcSMatthew G. Knepley 
1724d9407bcSMatthew G. Knepley           /* Can get rid of this loop by storing field information in the global section */
1734d9407bcSMatthew G. Knepley           for (f2 = 0; f2 < fields[f]; ++f2) {
1749566063dSJacob Faibussowitsch             PetscCall(PetscSectionGetFieldDof(section, p, f2, &fdof));
1759566063dSJacob Faibussowitsch             PetscCall(PetscSectionGetFieldConstraintDof(section, p, f2, &fcdof));
1764d9407bcSMatthew G. Knepley             poff += fdof - fcdof;
1774d9407bcSMatthew G. Knepley           }
1789566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetFieldDof(section, p, fields[f], &fdof));
1799566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetFieldConstraintDof(section, p, fields[f], &fcdof));
180ad540459SPierre Jolivet           for (fc = 0; fc < fdof - fcdof; ++fc, ++subOff) subIndices[subOff] = goff + poff + fc;
1814d9407bcSMatthew G. Knepley         }
1824d9407bcSMatthew G. Knepley       }
1834d9407bcSMatthew G. Knepley     }
1849566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), subSize, subIndices, PETSC_OWN_POINTER, is));
185627349f5SMatthew G. Knepley     if (bs > 1) {
186627349f5SMatthew G. Knepley       /* We need to check that the block size does not come from non-contiguous fields */
18788206212SAlexis Marboeuf       PetscInt i, j, set = 1, rset = 1;
188627349f5SMatthew G. Knepley       for (i = 0; i < subSize; i += bs) {
189627349f5SMatthew G. Knepley         for (j = 0; j < bs; ++j) {
1909371c9d4SSatish Balay           if (subIndices[i + j] != subIndices[i] + j) {
1919371c9d4SSatish Balay             set = 0;
1929371c9d4SSatish Balay             break;
1939371c9d4SSatish Balay           }
194627349f5SMatthew G. Knepley         }
195627349f5SMatthew G. Knepley       }
196712fec58SPierre Jolivet       PetscCall(MPIU_Allreduce(&set, &rset, 1, MPIU_INT, MPI_PROD, PetscObjectComm((PetscObject)dm)));
19788206212SAlexis Marboeuf       if (rset) PetscCall(ISSetBlockSize(*is, bs));
198627349f5SMatthew G. Knepley     }
1994d9407bcSMatthew G. Knepley   }
2004d9407bcSMatthew G. Knepley   if (subdm) {
2014d9407bcSMatthew G. Knepley     PetscSection subsection;
2024d9407bcSMatthew G. Knepley     PetscBool    haveNull = PETSC_FALSE;
2038cda7954SMatthew G. Knepley     PetscInt     f, nf = 0, of = 0;
2044d9407bcSMatthew G. Knepley 
2059566063dSJacob Faibussowitsch     PetscCall(PetscSectionCreateSubsection(section, numFields, fields, &subsection));
2069566063dSJacob Faibussowitsch     PetscCall(DMSetLocalSection(*subdm, subsection));
2079566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&subsection));
2084d9407bcSMatthew G. Knepley     for (f = 0; f < numFields; ++f) {
2094d9407bcSMatthew G. Knepley       (*subdm)->nullspaceConstructors[f] = dm->nullspaceConstructors[fields[f]];
2104d9407bcSMatthew G. Knepley       if ((*subdm)->nullspaceConstructors[f]) {
2114d9407bcSMatthew G. Knepley         haveNull = PETSC_TRUE;
2124d9407bcSMatthew G. Knepley         nf       = f;
2138cda7954SMatthew G. Knepley         of       = fields[f];
2144d9407bcSMatthew G. Knepley       }
2154d9407bcSMatthew G. Knepley     }
216e5e52638SMatthew G. Knepley     if (dm->probs) {
2179566063dSJacob Faibussowitsch       PetscCall(DMSetNumFields(*subdm, numFields));
2184d9407bcSMatthew G. Knepley       for (f = 0; f < numFields; ++f) {
2190f21e855SMatthew G. Knepley         PetscObject disc;
2200f21e855SMatthew G. Knepley 
2219566063dSJacob Faibussowitsch         PetscCall(DMGetField(dm, fields[f], NULL, &disc));
2229566063dSJacob Faibussowitsch         PetscCall(DMSetField(*subdm, f, NULL, disc));
2234d9407bcSMatthew G. Knepley       }
2249566063dSJacob Faibussowitsch       PetscCall(DMCreateDS(*subdm));
225f646a522SMatthew G. Knepley       if (numFields == 1 && is) {
2260f21e855SMatthew G. Knepley         PetscObject disc, space, pmat;
2274d9407bcSMatthew G. Knepley 
2289566063dSJacob Faibussowitsch         PetscCall(DMGetField(*subdm, 0, NULL, &disc));
2299566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "nullspace", &space));
2309566063dSJacob Faibussowitsch         if (space) PetscCall(PetscObjectCompose((PetscObject)*is, "nullspace", space));
2319566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "nearnullspace", &space));
2329566063dSJacob Faibussowitsch         if (space) PetscCall(PetscObjectCompose((PetscObject)*is, "nearnullspace", space));
2339566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "pmat", &pmat));
2349566063dSJacob Faibussowitsch         if (pmat) PetscCall(PetscObjectCompose((PetscObject)*is, "pmat", pmat));
2354d9407bcSMatthew G. Knepley       }
2369310035eSMatthew G. Knepley       /* Check if DSes record their DM fields */
237e21d936aSMatthew G. Knepley       if (dm->probs[0].fields) {
2389310035eSMatthew G. Knepley         PetscInt d, e;
2399310035eSMatthew G. Knepley 
2409310035eSMatthew G. Knepley         for (d = 0, e = 0; d < dm->Nds && e < (*subdm)->Nds; ++d) {
2419310035eSMatthew G. Knepley           const PetscInt  Nf = dm->probs[d].ds->Nf;
2429310035eSMatthew G. Knepley           const PetscInt *fld;
2439310035eSMatthew G. Knepley           PetscInt        f, g;
2449310035eSMatthew G. Knepley 
2459566063dSJacob Faibussowitsch           PetscCall(ISGetIndices(dm->probs[d].fields, &fld));
2469310035eSMatthew G. Knepley           for (f = 0; f < Nf; ++f) {
2479371c9d4SSatish Balay             for (g = 0; g < numFields; ++g)
2489371c9d4SSatish Balay               if (fld[f] == fields[g]) break;
2499310035eSMatthew G. Knepley             if (g < numFields) break;
2509310035eSMatthew G. Knepley           }
2519566063dSJacob Faibussowitsch           PetscCall(ISRestoreIndices(dm->probs[d].fields, &fld));
2529310035eSMatthew G. Knepley           if (f == Nf) continue;
2539566063dSJacob Faibussowitsch           PetscCall(PetscDSCopyConstants(dm->probs[d].ds, (*subdm)->probs[e].ds));
2549566063dSJacob Faibussowitsch           PetscCall(PetscDSCopyBoundary(dm->probs[d].ds, numFields, fields, (*subdm)->probs[e].ds));
2559310035eSMatthew G. Knepley           /* Translate DM fields to DS fields */
2569310035eSMatthew G. Knepley           {
257e21d936aSMatthew G. Knepley             IS              infields, dsfields;
25874a61aaaSMatthew G. Knepley             const PetscInt *fld, *ofld;
25974a61aaaSMatthew G. Knepley             PetscInt       *fidx;
26074a61aaaSMatthew G. Knepley             PetscInt        onf, nf, f, g;
261e21d936aSMatthew G. Knepley 
2629566063dSJacob Faibussowitsch             PetscCall(ISCreateGeneral(PETSC_COMM_SELF, numFields, fields, PETSC_USE_POINTER, &infields));
2639566063dSJacob Faibussowitsch             PetscCall(ISIntersect(infields, dm->probs[d].fields, &dsfields));
2649566063dSJacob Faibussowitsch             PetscCall(ISDestroy(&infields));
2659566063dSJacob Faibussowitsch             PetscCall(ISGetLocalSize(dsfields, &nf));
2667a8be351SBarry Smith             PetscCheck(nf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "DS cannot be supported on 0 fields");
2679566063dSJacob Faibussowitsch             PetscCall(ISGetIndices(dsfields, &fld));
2689566063dSJacob Faibussowitsch             PetscCall(ISGetLocalSize(dm->probs[d].fields, &onf));
2699566063dSJacob Faibussowitsch             PetscCall(ISGetIndices(dm->probs[d].fields, &ofld));
2709566063dSJacob Faibussowitsch             PetscCall(PetscMalloc1(nf, &fidx));
2713268a0aaSMatthew G. Knepley             for (f = 0, g = 0; f < onf && g < nf; ++f) {
27274a61aaaSMatthew G. Knepley               if (ofld[f] == fld[g]) fidx[g++] = f;
27374a61aaaSMatthew G. Knepley             }
2749566063dSJacob Faibussowitsch             PetscCall(ISRestoreIndices(dm->probs[d].fields, &ofld));
2759566063dSJacob Faibussowitsch             PetscCall(ISRestoreIndices(dsfields, &fld));
2769566063dSJacob Faibussowitsch             PetscCall(ISDestroy(&dsfields));
2779566063dSJacob Faibussowitsch             PetscCall(PetscDSSelectDiscretizations(dm->probs[0].ds, nf, fidx, (*subdm)->probs[0].ds));
2789566063dSJacob Faibussowitsch             PetscCall(PetscDSSelectEquations(dm->probs[0].ds, nf, fidx, (*subdm)->probs[0].ds));
2799566063dSJacob Faibussowitsch             PetscCall(PetscFree(fidx));
2809310035eSMatthew G. Knepley           }
2819310035eSMatthew G. Knepley           ++e;
2829310035eSMatthew G. Knepley         }
283e21d936aSMatthew G. Knepley       } else {
2849566063dSJacob Faibussowitsch         PetscCall(PetscDSCopyConstants(dm->probs[0].ds, (*subdm)->probs[0].ds));
2859566063dSJacob Faibussowitsch         PetscCall(PetscDSCopyBoundary(dm->probs[0].ds, PETSC_DETERMINE, NULL, (*subdm)->probs[0].ds));
2869566063dSJacob Faibussowitsch         PetscCall(PetscDSSelectDiscretizations(dm->probs[0].ds, numFields, fields, (*subdm)->probs[0].ds));
2879566063dSJacob Faibussowitsch         PetscCall(PetscDSSelectEquations(dm->probs[0].ds, numFields, fields, (*subdm)->probs[0].ds));
288d5af271aSMatthew G. Knepley       }
289e21d936aSMatthew G. Knepley     }
2908cda7954SMatthew G. Knepley     if (haveNull && is) {
2918cda7954SMatthew G. Knepley       MatNullSpace nullSpace;
2928cda7954SMatthew G. Knepley 
2939566063dSJacob Faibussowitsch       PetscCall((*(*subdm)->nullspaceConstructors[nf])(*subdm, of, nf, &nullSpace));
2949566063dSJacob Faibussowitsch       PetscCall(PetscObjectCompose((PetscObject)*is, "nullspace", (PetscObject)nullSpace));
2959566063dSJacob Faibussowitsch       PetscCall(MatNullSpaceDestroy(&nullSpace));
2968cda7954SMatthew G. Knepley     }
29748a46eb9SPierre Jolivet     if (dm->coarseMesh) PetscCall(DMCreateSubDM(dm->coarseMesh, numFields, fields, NULL, &(*subdm)->coarseMesh));
2984d9407bcSMatthew G. Knepley   }
2993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3004d9407bcSMatthew G. Knepley }
3012adcc780SMatthew G. Knepley 
302792b654fSMatthew G. Knepley /*@C
30320f4b53cSBarry Smith   DMCreateSectionSuperDM - Returns an arrays of `IS` and DM+Section encapsulating a superproblem defined by the DM+Sections passed in.
304792b654fSMatthew G. Knepley 
30520f4b53cSBarry Smith   Not Collective
306792b654fSMatthew G. Knepley 
307d8d19677SJose E. Roman   Input Parameters:
30820f4b53cSBarry Smith + dms - The `DM` objects
30920f4b53cSBarry Smith - len - The number of `DM`s
310792b654fSMatthew G. Knepley 
311792b654fSMatthew G. Knepley   Output Parameters:
31220f4b53cSBarry Smith + is      - The global indices for the subproblem, or `NULL`
31320f4b53cSBarry Smith - superdm - The `DM` for the superproblem, which must already have be cloned
314792b654fSMatthew G. Knepley 
315792b654fSMatthew G. Knepley   Level: intermediate
316792b654fSMatthew G. Knepley 
31720f4b53cSBarry Smith   Note:
31820f4b53cSBarry Smith   This handles all information in the `DM` class and the `PetscSection`. This is used as the basis for creating subDMs in specialized classes,
31920f4b53cSBarry Smith   such as `DMPLEX` and `DMFOREST`
32020f4b53cSBarry Smith 
321*a4e35b19SJacob Faibussowitsch .seealso: `DMCreateSuperDM()`, `DMGetLocalSection()`, `DMPlexSetMigrationSF()`, `DMView()`
322792b654fSMatthew G. Knepley @*/
323d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateSectionSuperDM(DM dms[], PetscInt len, IS **is, DM *superdm)
324d71ae5a4SJacob Faibussowitsch {
3259796d432SMatthew G. Knepley   MPI_Comm     comm;
3269796d432SMatthew G. Knepley   PetscSection supersection, *sections, *sectionGlobals;
3278cda7954SMatthew G. Knepley   PetscInt    *Nfs, Nf = 0, f, supf, oldf = -1, nullf = -1, i;
3289796d432SMatthew G. Knepley   PetscBool    haveNull = PETSC_FALSE;
3292adcc780SMatthew G. Knepley 
3302adcc780SMatthew G. Knepley   PetscFunctionBegin;
3319566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dms[0], &comm));
3329796d432SMatthew G. Knepley   /* Pull out local and global sections */
3339566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(len, &Nfs, len, &sections, len, &sectionGlobals));
3342adcc780SMatthew G. Knepley   for (i = 0; i < len; ++i) {
3359566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dms[i], &sections[i]));
3369566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalSection(dms[i], &sectionGlobals[i]));
3377a8be351SBarry Smith     PetscCheck(sections[i], comm, PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields");
3387a8be351SBarry Smith     PetscCheck(sectionGlobals[i], comm, PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields");
3399566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetNumFields(sections[i], &Nfs[i]));
3402adcc780SMatthew G. Knepley     Nf += Nfs[i];
3412adcc780SMatthew G. Knepley   }
3429796d432SMatthew G. Knepley   /* Create the supersection */
3439566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreateSupersection(sections, len, &supersection));
3449566063dSJacob Faibussowitsch   PetscCall(DMSetLocalSection(*superdm, supersection));
3459796d432SMatthew G. Knepley   /* Create ISes */
3462adcc780SMatthew G. Knepley   if (is) {
3479796d432SMatthew G. Knepley     PetscSection supersectionGlobal;
3489796d432SMatthew G. Knepley     PetscInt     bs = -1, startf = 0;
3492adcc780SMatthew G. Knepley 
3509566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(len, is));
3519566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalSection(*superdm, &supersectionGlobal));
3529796d432SMatthew G. Knepley     for (i = 0; i < len; startf += Nfs[i], ++i) {
3539796d432SMatthew G. Knepley       PetscInt *subIndices;
354ec4c761aSMatthew G. Knepley       PetscInt  subSize, subOff, pStart, pEnd, p, start, end, dummy;
3552adcc780SMatthew G. Knepley 
3569566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetChart(sectionGlobals[i], &pStart, &pEnd));
3579566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetConstrainedStorageSize(sectionGlobals[i], &subSize));
3589566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(subSize, &subIndices));
3599796d432SMatthew G. Knepley       for (p = pStart, subOff = 0; p < pEnd; ++p) {
3609796d432SMatthew G. Knepley         PetscInt gdof, gcdof, gtdof, d;
3612adcc780SMatthew G. Knepley 
3629566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(sectionGlobals[i], p, &gdof));
3639566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetConstraintDof(sections[i], p, &gcdof));
3642adcc780SMatthew G. Knepley         gtdof = gdof - gcdof;
3652adcc780SMatthew G. Knepley         if (gdof > 0 && gtdof) {
3669371c9d4SSatish Balay           if (bs < 0) {
3679371c9d4SSatish Balay             bs = gtdof;
3689371c9d4SSatish Balay           } else if (bs != gtdof) {
3699371c9d4SSatish Balay             bs = 1;
3709371c9d4SSatish Balay           }
3719566063dSJacob Faibussowitsch           PetscCall(DMGetGlobalFieldOffset_Private(*superdm, p, startf, &start, &dummy));
3729566063dSJacob Faibussowitsch           PetscCall(DMGetGlobalFieldOffset_Private(*superdm, p, startf + Nfs[i] - 1, &dummy, &end));
37363a3b9bcSJacob 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);
3749796d432SMatthew G. Knepley           for (d = start; d < end; ++d, ++subOff) subIndices[subOff] = d;
3752adcc780SMatthew G. Knepley         }
3762adcc780SMatthew G. Knepley       }
3779566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(comm, subSize, subIndices, PETSC_OWN_POINTER, &(*is)[i]));
3782adcc780SMatthew G. Knepley       /* Must have same blocksize on all procs (some might have no points) */
3799796d432SMatthew G. Knepley       {
3809796d432SMatthew G. Knepley         PetscInt bs = -1, bsLocal[2], bsMinMax[2];
3819796d432SMatthew G. Knepley 
3829371c9d4SSatish Balay         bsLocal[0] = bs < 0 ? PETSC_MAX_INT : bs;
3839371c9d4SSatish Balay         bsLocal[1] = bs;
3849566063dSJacob Faibussowitsch         PetscCall(PetscGlobalMinMaxInt(comm, bsLocal, bsMinMax));
3859371c9d4SSatish Balay         if (bsMinMax[0] != bsMinMax[1]) {
3869371c9d4SSatish Balay           bs = 1;
3879371c9d4SSatish Balay         } else {
3889371c9d4SSatish Balay           bs = bsMinMax[0];
3899371c9d4SSatish Balay         }
3909566063dSJacob Faibussowitsch         PetscCall(ISSetBlockSize((*is)[i], bs));
3912adcc780SMatthew G. Knepley       }
3922adcc780SMatthew G. Knepley     }
3932adcc780SMatthew G. Knepley   }
3949796d432SMatthew G. Knepley   /* Preserve discretizations */
395e5e52638SMatthew G. Knepley   if (len && dms[0]->probs) {
3969566063dSJacob Faibussowitsch     PetscCall(DMSetNumFields(*superdm, Nf));
3979796d432SMatthew G. Knepley     for (i = 0, supf = 0; i < len; ++i) {
3989796d432SMatthew G. Knepley       for (f = 0; f < Nfs[i]; ++f, ++supf) {
3992adcc780SMatthew G. Knepley         PetscObject disc;
4002adcc780SMatthew G. Knepley 
4019566063dSJacob Faibussowitsch         PetscCall(DMGetField(dms[i], f, NULL, &disc));
4029566063dSJacob Faibussowitsch         PetscCall(DMSetField(*superdm, supf, NULL, disc));
4032adcc780SMatthew G. Knepley       }
4042adcc780SMatthew G. Knepley     }
4059566063dSJacob Faibussowitsch     PetscCall(DMCreateDS(*superdm));
4062adcc780SMatthew G. Knepley   }
4079796d432SMatthew G. Knepley   /* Preserve nullspaces */
4089796d432SMatthew G. Knepley   for (i = 0, supf = 0; i < len; ++i) {
4099796d432SMatthew G. Knepley     for (f = 0; f < Nfs[i]; ++f, ++supf) {
4109796d432SMatthew G. Knepley       (*superdm)->nullspaceConstructors[supf] = dms[i]->nullspaceConstructors[f];
4119796d432SMatthew G. Knepley       if ((*superdm)->nullspaceConstructors[supf]) {
4129796d432SMatthew G. Knepley         haveNull = PETSC_TRUE;
4139796d432SMatthew G. Knepley         nullf    = supf;
4148cda7954SMatthew G. Knepley         oldf     = f;
4152adcc780SMatthew G. Knepley       }
4169796d432SMatthew G. Knepley     }
4179796d432SMatthew G. Knepley   }
4189796d432SMatthew G. Knepley   /* Attach nullspace to IS */
4199796d432SMatthew G. Knepley   if (haveNull && is) {
4209796d432SMatthew G. Knepley     MatNullSpace nullSpace;
4219796d432SMatthew G. Knepley 
4229566063dSJacob Faibussowitsch     PetscCall((*(*superdm)->nullspaceConstructors[nullf])(*superdm, oldf, nullf, &nullSpace));
4239566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)(*is)[nullf], "nullspace", (PetscObject)nullSpace));
4249566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceDestroy(&nullSpace));
4259796d432SMatthew G. Knepley   }
4269566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&supersection));
4279566063dSJacob Faibussowitsch   PetscCall(PetscFree3(Nfs, sections, sectionGlobals));
4283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4292adcc780SMatthew G. Knepley }
430