xref: /petsc/src/dm/interface/dmi.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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
55b8243e1SJed Brown PetscInt PetscGCD(PetscInt a, PetscInt b) {
65b8243e1SJed Brown   while (b != 0) {
75b8243e1SJed Brown     PetscInt tmp = b;
85b8243e1SJed Brown     b            = a % b;
95b8243e1SJed Brown     a            = tmp;
105b8243e1SJed Brown   }
115b8243e1SJed Brown   return a;
125b8243e1SJed Brown }
135b8243e1SJed Brown 
149371c9d4SSatish Balay PetscErrorCode DMCreateGlobalVector_Section_Private(DM dm, Vec *vec) {
1511689aebSJed Brown   PetscSection gSection;
16cc30b04dSMatthew G Knepley   PetscInt     localSize, bs, blockSize = -1, pStart, pEnd, p;
175d1c0279SHong Zhang   PetscInt     in[2], out[2];
1811689aebSJed Brown 
1911689aebSJed Brown   PetscFunctionBegin;
209566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalSection(dm, &gSection));
219566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(gSection, &pStart, &pEnd));
2211689aebSJed Brown   for (p = pStart; p < pEnd; ++p) {
2311689aebSJed Brown     PetscInt dof, cdof;
2411689aebSJed Brown 
259566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(gSection, p, &dof));
269566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(gSection, p, &cdof));
270680ce0aSHong Zhang 
285b8243e1SJed Brown     if (dof - cdof > 0) {
295b8243e1SJed Brown       if (blockSize < 0) {
300680ce0aSHong Zhang         /* set blockSize */
310680ce0aSHong Zhang         blockSize = dof - cdof;
325b8243e1SJed Brown       } else {
335b8243e1SJed Brown         blockSize = PetscGCD(dof - cdof, blockSize);
3411689aebSJed Brown       }
3511689aebSJed Brown     }
360680ce0aSHong Zhang   }
375d1c0279SHong Zhang 
382c8be454SMatthew G. Knepley   in[0] = blockSize < 0 ? PETSC_MIN_INT : -blockSize;
395d1c0279SHong Zhang   in[1] = blockSize;
401c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(in, out, 2, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm)));
415d1c0279SHong Zhang   /* -out[0] = min(blockSize), out[1] = max(blockSize) */
425d1c0279SHong Zhang   if (-out[0] == out[1]) {
435d1c0279SHong Zhang     bs = out[1];
445d1c0279SHong Zhang   } else bs = 1;
455d1c0279SHong Zhang 
465d1c0279SHong Zhang   if (bs < 0) { /* Everyone was empty */
475d1c0279SHong Zhang     blockSize = 1;
485d1c0279SHong Zhang     bs        = 1;
495d1c0279SHong Zhang   }
505d1c0279SHong Zhang 
519566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetConstrainedStorageSize(gSection, &localSize));
527a8be351SBarry 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);
539566063dSJacob Faibussowitsch   PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), vec));
549566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(*vec, localSize, PETSC_DETERMINE));
559566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(*vec, bs));
569566063dSJacob Faibussowitsch   PetscCall(VecSetType(*vec, dm->vectype));
579566063dSJacob Faibussowitsch   PetscCall(VecSetDM(*vec, dm));
589566063dSJacob Faibussowitsch   /* PetscCall(VecSetLocalToGlobalMapping(*vec, dm->ltogmap)); */
5911689aebSJed Brown   PetscFunctionReturn(0);
6011689aebSJed Brown }
6111689aebSJed Brown 
629371c9d4SSatish Balay PetscErrorCode DMCreateLocalVector_Section_Private(DM dm, Vec *vec) {
63fad22124SMatthew G Knepley   PetscSection section;
6411689aebSJed Brown   PetscInt     localSize, blockSize = -1, pStart, pEnd, p;
6511689aebSJed Brown 
6611689aebSJed Brown   PetscFunctionBegin;
679566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
689566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
6911689aebSJed Brown   for (p = pStart; p < pEnd; ++p) {
7011689aebSJed Brown     PetscInt dof;
7111689aebSJed Brown 
729566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
7311689aebSJed Brown     if ((blockSize < 0) && (dof > 0)) blockSize = dof;
745b8243e1SJed Brown     if (dof > 0) blockSize = PetscGCD(dof, blockSize);
7511689aebSJed Brown   }
769566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(section, &localSize));
779566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF, vec));
789566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(*vec, localSize, localSize));
799566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(*vec, blockSize));
809566063dSJacob Faibussowitsch   PetscCall(VecSetType(*vec, dm->vectype));
819566063dSJacob Faibussowitsch   PetscCall(VecSetDM(*vec, dm));
8211689aebSJed Brown   PetscFunctionReturn(0);
8311689aebSJed Brown }
844d9407bcSMatthew G. Knepley 
85792b654fSMatthew G. Knepley /*@C
86792b654fSMatthew G. Knepley   DMCreateSectionSubDM - Returns an IS and subDM+subSection encapsulating a subproblem defined by the fields in a PetscSection in the DM.
87792b654fSMatthew G. Knepley 
88792b654fSMatthew G. Knepley   Not collective
89792b654fSMatthew G. Knepley 
90792b654fSMatthew G. Knepley   Input Parameters:
91792b654fSMatthew G. Knepley + dm        - The DM object
92792b654fSMatthew G. Knepley . numFields - The number of fields in this subproblem
93792b654fSMatthew G. Knepley - fields    - The field numbers of the selected fields
94792b654fSMatthew G. Knepley 
95792b654fSMatthew G. Knepley   Output Parameters:
96792b654fSMatthew G. Knepley + is - The global indices for the subproblem
97792b654fSMatthew G. Knepley - subdm - The DM for the subproblem, which must already have be cloned from dm
98792b654fSMatthew G. Knepley 
99792b654fSMatthew G. Knepley   Note: This handles all information in the DM class and the PetscSection. This is used as the basis for creating subDMs in specialized classes,
100792b654fSMatthew G. Knepley   such as Plex and Forest.
101792b654fSMatthew G. Knepley 
102792b654fSMatthew G. Knepley   Level: intermediate
103792b654fSMatthew G. Knepley 
104db781477SPatrick Sanan .seealso `DMCreateSubDM()`, `DMGetLocalSection()`, `DMPlexSetMigrationSF()`, `DMView()`
105792b654fSMatthew G. Knepley @*/
1069371c9d4SSatish Balay PetscErrorCode DMCreateSectionSubDM(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm) {
1074d9407bcSMatthew G. Knepley   PetscSection section, sectionGlobal;
1084d9407bcSMatthew G. Knepley   PetscInt    *subIndices;
109696188eaSMatthew G. Knepley   PetscInt     subSize = 0, subOff = 0, Nf, f, pStart, pEnd, p;
1104d9407bcSMatthew G. Knepley 
1114d9407bcSMatthew G. Knepley   PetscFunctionBegin;
1124d9407bcSMatthew G. Knepley   if (!numFields) PetscFunctionReturn(0);
1139566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
1149566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalSection(dm, &sectionGlobal));
1157a8be351SBarry Smith   PetscCheck(section, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields");
1167a8be351SBarry Smith   PetscCheck(sectionGlobal, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields");
1179566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetNumFields(section, &Nf));
1187a8be351SBarry 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);
1194d9407bcSMatthew G. Knepley   if (is) {
1203a544194SStefano Zampini     PetscInt bs, bsLocal[2], bsMinMax[2];
1213a544194SStefano Zampini 
1223a544194SStefano Zampini     for (f = 0, bs = 0; f < numFields; ++f) {
1233a544194SStefano Zampini       PetscInt Nc;
1243a544194SStefano Zampini 
1259566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldComponents(section, fields[f], &Nc));
1263a544194SStefano Zampini       bs += Nc;
1273a544194SStefano Zampini     }
1289566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetChart(sectionGlobal, &pStart, &pEnd));
1294d9407bcSMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
130b9eca265Ssarens       PetscInt gdof, pSubSize = 0;
1314d9407bcSMatthew G. Knepley 
1329566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(sectionGlobal, p, &gdof));
1334d9407bcSMatthew G. Knepley       if (gdof > 0) {
1344d9407bcSMatthew G. Knepley         for (f = 0; f < numFields; ++f) {
1354d9407bcSMatthew G. Knepley           PetscInt fdof, fcdof;
1364d9407bcSMatthew G. Knepley 
1379566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetFieldDof(section, p, fields[f], &fdof));
1389566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetFieldConstraintDof(section, p, fields[f], &fcdof));
139b9eca265Ssarens           pSubSize += fdof - fcdof;
140b9eca265Ssarens         }
141b9eca265Ssarens         subSize += pSubSize;
1423a544194SStefano Zampini         if (pSubSize && bs != pSubSize) {
143b9eca265Ssarens           /* Layout does not admit a pointwise block size */
144b9eca265Ssarens           bs = 1;
1454d9407bcSMatthew G. Knepley         }
1464d9407bcSMatthew G. Knepley       }
1474d9407bcSMatthew G. Knepley     }
148b9eca265Ssarens     /* Must have same blocksize on all procs (some might have no points) */
1499371c9d4SSatish Balay     bsLocal[0] = bs < 0 ? PETSC_MAX_INT : bs;
1509371c9d4SSatish Balay     bsLocal[1] = bs;
1519566063dSJacob Faibussowitsch     PetscCall(PetscGlobalMinMaxInt(PetscObjectComm((PetscObject)dm), bsLocal, bsMinMax));
1529371c9d4SSatish Balay     if (bsMinMax[0] != bsMinMax[1]) {
1539371c9d4SSatish Balay       bs = 1;
1549371c9d4SSatish Balay     } else {
1559371c9d4SSatish Balay       bs = bsMinMax[0];
1569371c9d4SSatish Balay     }
1579566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(subSize, &subIndices));
1584d9407bcSMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
1594d9407bcSMatthew G. Knepley       PetscInt gdof, goff;
1604d9407bcSMatthew G. Knepley 
1619566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(sectionGlobal, p, &gdof));
1624d9407bcSMatthew G. Knepley       if (gdof > 0) {
1639566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetOffset(sectionGlobal, p, &goff));
1644d9407bcSMatthew G. Knepley         for (f = 0; f < numFields; ++f) {
1654d9407bcSMatthew G. Knepley           PetscInt fdof, fcdof, fc, f2, poff = 0;
1664d9407bcSMatthew G. Knepley 
1674d9407bcSMatthew G. Knepley           /* Can get rid of this loop by storing field information in the global section */
1684d9407bcSMatthew G. Knepley           for (f2 = 0; f2 < fields[f]; ++f2) {
1699566063dSJacob Faibussowitsch             PetscCall(PetscSectionGetFieldDof(section, p, f2, &fdof));
1709566063dSJacob Faibussowitsch             PetscCall(PetscSectionGetFieldConstraintDof(section, p, f2, &fcdof));
1714d9407bcSMatthew G. Knepley             poff += fdof - fcdof;
1724d9407bcSMatthew G. Knepley           }
1739566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetFieldDof(section, p, fields[f], &fdof));
1749566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetFieldConstraintDof(section, p, fields[f], &fcdof));
1759371c9d4SSatish Balay           for (fc = 0; fc < fdof - fcdof; ++fc, ++subOff) { subIndices[subOff] = goff + poff + fc; }
1764d9407bcSMatthew G. Knepley         }
1774d9407bcSMatthew G. Knepley       }
1784d9407bcSMatthew G. Knepley     }
1799566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), subSize, subIndices, PETSC_OWN_POINTER, is));
180627349f5SMatthew G. Knepley     if (bs > 1) {
181627349f5SMatthew G. Knepley       /* We need to check that the block size does not come from non-contiguous fields */
18288206212SAlexis Marboeuf       PetscInt i, j, set = 1, rset = 1;
183627349f5SMatthew G. Knepley       for (i = 0; i < subSize; i += bs) {
184627349f5SMatthew G. Knepley         for (j = 0; j < bs; ++j) {
1859371c9d4SSatish Balay           if (subIndices[i + j] != subIndices[i] + j) {
1869371c9d4SSatish Balay             set = 0;
1879371c9d4SSatish Balay             break;
1889371c9d4SSatish Balay           }
189627349f5SMatthew G. Knepley         }
190627349f5SMatthew G. Knepley       }
19188206212SAlexis Marboeuf       PetscCallMPI(MPI_Allreduce(&set, &rset, 1, MPIU_INT, MPI_PROD, PetscObjectComm((PetscObject)dm)));
19288206212SAlexis Marboeuf       if (rset) PetscCall(ISSetBlockSize(*is, bs));
193627349f5SMatthew G. Knepley     }
1944d9407bcSMatthew G. Knepley   }
1954d9407bcSMatthew G. Knepley   if (subdm) {
1964d9407bcSMatthew G. Knepley     PetscSection subsection;
1974d9407bcSMatthew G. Knepley     PetscBool    haveNull = PETSC_FALSE;
1988cda7954SMatthew G. Knepley     PetscInt     f, nf = 0, of = 0;
1994d9407bcSMatthew G. Knepley 
2009566063dSJacob Faibussowitsch     PetscCall(PetscSectionCreateSubsection(section, numFields, fields, &subsection));
2019566063dSJacob Faibussowitsch     PetscCall(DMSetLocalSection(*subdm, subsection));
2029566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&subsection));
2034d9407bcSMatthew G. Knepley     for (f = 0; f < numFields; ++f) {
2044d9407bcSMatthew G. Knepley       (*subdm)->nullspaceConstructors[f] = dm->nullspaceConstructors[fields[f]];
2054d9407bcSMatthew G. Knepley       if ((*subdm)->nullspaceConstructors[f]) {
2064d9407bcSMatthew G. Knepley         haveNull = PETSC_TRUE;
2074d9407bcSMatthew G. Knepley         nf       = f;
2088cda7954SMatthew G. Knepley         of       = fields[f];
2094d9407bcSMatthew G. Knepley       }
2104d9407bcSMatthew G. Knepley     }
211e5e52638SMatthew G. Knepley     if (dm->probs) {
2129566063dSJacob Faibussowitsch       PetscCall(DMSetNumFields(*subdm, numFields));
2134d9407bcSMatthew G. Knepley       for (f = 0; f < numFields; ++f) {
2140f21e855SMatthew G. Knepley         PetscObject disc;
2150f21e855SMatthew G. Knepley 
2169566063dSJacob Faibussowitsch         PetscCall(DMGetField(dm, fields[f], NULL, &disc));
2179566063dSJacob Faibussowitsch         PetscCall(DMSetField(*subdm, f, NULL, disc));
2184d9407bcSMatthew G. Knepley       }
2199566063dSJacob Faibussowitsch       PetscCall(DMCreateDS(*subdm));
220f646a522SMatthew G. Knepley       if (numFields == 1 && is) {
2210f21e855SMatthew G. Knepley         PetscObject disc, space, pmat;
2224d9407bcSMatthew G. Knepley 
2239566063dSJacob Faibussowitsch         PetscCall(DMGetField(*subdm, 0, NULL, &disc));
2249566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "nullspace", &space));
2259566063dSJacob Faibussowitsch         if (space) PetscCall(PetscObjectCompose((PetscObject)*is, "nullspace", space));
2269566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "nearnullspace", &space));
2279566063dSJacob Faibussowitsch         if (space) PetscCall(PetscObjectCompose((PetscObject)*is, "nearnullspace", space));
2289566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "pmat", &pmat));
2299566063dSJacob Faibussowitsch         if (pmat) PetscCall(PetscObjectCompose((PetscObject)*is, "pmat", pmat));
2304d9407bcSMatthew G. Knepley       }
2319310035eSMatthew G. Knepley       /* Check if DSes record their DM fields */
232e21d936aSMatthew G. Knepley       if (dm->probs[0].fields) {
2339310035eSMatthew G. Knepley         PetscInt d, e;
2349310035eSMatthew G. Knepley 
2359310035eSMatthew G. Knepley         for (d = 0, e = 0; d < dm->Nds && e < (*subdm)->Nds; ++d) {
2369310035eSMatthew G. Knepley           const PetscInt  Nf = dm->probs[d].ds->Nf;
2379310035eSMatthew G. Knepley           const PetscInt *fld;
2389310035eSMatthew G. Knepley           PetscInt        f, g;
2399310035eSMatthew G. Knepley 
2409566063dSJacob Faibussowitsch           PetscCall(ISGetIndices(dm->probs[d].fields, &fld));
2419310035eSMatthew G. Knepley           for (f = 0; f < Nf; ++f) {
2429371c9d4SSatish Balay             for (g = 0; g < numFields; ++g)
2439371c9d4SSatish Balay               if (fld[f] == fields[g]) break;
2449310035eSMatthew G. Knepley             if (g < numFields) break;
2459310035eSMatthew G. Knepley           }
2469566063dSJacob Faibussowitsch           PetscCall(ISRestoreIndices(dm->probs[d].fields, &fld));
2479310035eSMatthew G. Knepley           if (f == Nf) continue;
2489566063dSJacob Faibussowitsch           PetscCall(PetscDSCopyConstants(dm->probs[d].ds, (*subdm)->probs[e].ds));
2499566063dSJacob Faibussowitsch           PetscCall(PetscDSCopyBoundary(dm->probs[d].ds, numFields, fields, (*subdm)->probs[e].ds));
2509310035eSMatthew G. Knepley           /* Translate DM fields to DS fields */
2519310035eSMatthew G. Knepley           {
252e21d936aSMatthew G. Knepley             IS              infields, dsfields;
25374a61aaaSMatthew G. Knepley             const PetscInt *fld, *ofld;
25474a61aaaSMatthew G. Knepley             PetscInt       *fidx;
25574a61aaaSMatthew G. Knepley             PetscInt        onf, nf, f, g;
256e21d936aSMatthew G. Knepley 
2579566063dSJacob Faibussowitsch             PetscCall(ISCreateGeneral(PETSC_COMM_SELF, numFields, fields, PETSC_USE_POINTER, &infields));
2589566063dSJacob Faibussowitsch             PetscCall(ISIntersect(infields, dm->probs[d].fields, &dsfields));
2599566063dSJacob Faibussowitsch             PetscCall(ISDestroy(&infields));
2609566063dSJacob Faibussowitsch             PetscCall(ISGetLocalSize(dsfields, &nf));
2617a8be351SBarry Smith             PetscCheck(nf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "DS cannot be supported on 0 fields");
2629566063dSJacob Faibussowitsch             PetscCall(ISGetIndices(dsfields, &fld));
2639566063dSJacob Faibussowitsch             PetscCall(ISGetLocalSize(dm->probs[d].fields, &onf));
2649566063dSJacob Faibussowitsch             PetscCall(ISGetIndices(dm->probs[d].fields, &ofld));
2659566063dSJacob Faibussowitsch             PetscCall(PetscMalloc1(nf, &fidx));
2663268a0aaSMatthew G. Knepley             for (f = 0, g = 0; f < onf && g < nf; ++f) {
26774a61aaaSMatthew G. Knepley               if (ofld[f] == fld[g]) fidx[g++] = f;
26874a61aaaSMatthew G. Knepley             }
2699566063dSJacob Faibussowitsch             PetscCall(ISRestoreIndices(dm->probs[d].fields, &ofld));
2709566063dSJacob Faibussowitsch             PetscCall(ISRestoreIndices(dsfields, &fld));
2719566063dSJacob Faibussowitsch             PetscCall(ISDestroy(&dsfields));
2729566063dSJacob Faibussowitsch             PetscCall(PetscDSSelectDiscretizations(dm->probs[0].ds, nf, fidx, (*subdm)->probs[0].ds));
2739566063dSJacob Faibussowitsch             PetscCall(PetscDSSelectEquations(dm->probs[0].ds, nf, fidx, (*subdm)->probs[0].ds));
2749566063dSJacob Faibussowitsch             PetscCall(PetscFree(fidx));
2759310035eSMatthew G. Knepley           }
2769310035eSMatthew G. Knepley           ++e;
2779310035eSMatthew G. Knepley         }
278e21d936aSMatthew G. Knepley       } else {
2799566063dSJacob Faibussowitsch         PetscCall(PetscDSCopyConstants(dm->probs[0].ds, (*subdm)->probs[0].ds));
2809566063dSJacob Faibussowitsch         PetscCall(PetscDSCopyBoundary(dm->probs[0].ds, PETSC_DETERMINE, NULL, (*subdm)->probs[0].ds));
2819566063dSJacob Faibussowitsch         PetscCall(PetscDSSelectDiscretizations(dm->probs[0].ds, numFields, fields, (*subdm)->probs[0].ds));
2829566063dSJacob Faibussowitsch         PetscCall(PetscDSSelectEquations(dm->probs[0].ds, numFields, fields, (*subdm)->probs[0].ds));
283d5af271aSMatthew G. Knepley       }
284e21d936aSMatthew G. Knepley     }
2858cda7954SMatthew G. Knepley     if (haveNull && is) {
2868cda7954SMatthew G. Knepley       MatNullSpace nullSpace;
2878cda7954SMatthew G. Knepley 
2889566063dSJacob Faibussowitsch       PetscCall((*(*subdm)->nullspaceConstructors[nf])(*subdm, of, nf, &nullSpace));
2899566063dSJacob Faibussowitsch       PetscCall(PetscObjectCompose((PetscObject)*is, "nullspace", (PetscObject)nullSpace));
2909566063dSJacob Faibussowitsch       PetscCall(MatNullSpaceDestroy(&nullSpace));
2918cda7954SMatthew G. Knepley     }
292*48a46eb9SPierre Jolivet     if (dm->coarseMesh) PetscCall(DMCreateSubDM(dm->coarseMesh, numFields, fields, NULL, &(*subdm)->coarseMesh));
2934d9407bcSMatthew G. Knepley   }
2944d9407bcSMatthew G. Knepley   PetscFunctionReturn(0);
2954d9407bcSMatthew G. Knepley }
2962adcc780SMatthew G. Knepley 
297792b654fSMatthew G. Knepley /*@C
298792b654fSMatthew G. Knepley   DMCreateSectionSuperDM - Returns an arrays of ISes and DM+Section encapsulating a superproblem defined by the DM+Sections passed in.
299792b654fSMatthew G. Knepley 
300792b654fSMatthew G. Knepley   Not collective
301792b654fSMatthew G. Knepley 
302d8d19677SJose E. Roman   Input Parameters:
303792b654fSMatthew G. Knepley + dms - The DM objects
304792b654fSMatthew G. Knepley - len - The number of DMs
305792b654fSMatthew G. Knepley 
306792b654fSMatthew G. Knepley   Output Parameters:
307792b654fSMatthew G. Knepley + is - The global indices for the subproblem, or NULL
308792b654fSMatthew G. Knepley - superdm - The DM for the superproblem, which must already have be cloned
309792b654fSMatthew G. Knepley 
310792b654fSMatthew G. Knepley   Note: This handles all information in the DM class and the PetscSection. This is used as the basis for creating subDMs in specialized classes,
311792b654fSMatthew G. Knepley   such as Plex and Forest.
312792b654fSMatthew G. Knepley 
313792b654fSMatthew G. Knepley   Level: intermediate
314792b654fSMatthew G. Knepley 
315db781477SPatrick Sanan .seealso `DMCreateSuperDM()`, `DMGetLocalSection()`, `DMPlexSetMigrationSF()`, `DMView()`
316792b654fSMatthew G. Knepley @*/
3179371c9d4SSatish Balay PetscErrorCode DMCreateSectionSuperDM(DM dms[], PetscInt len, IS **is, DM *superdm) {
3189796d432SMatthew G. Knepley   MPI_Comm     comm;
3199796d432SMatthew G. Knepley   PetscSection supersection, *sections, *sectionGlobals;
3208cda7954SMatthew G. Knepley   PetscInt    *Nfs, Nf = 0, f, supf, oldf = -1, nullf = -1, i;
3219796d432SMatthew G. Knepley   PetscBool    haveNull = PETSC_FALSE;
3222adcc780SMatthew G. Knepley 
3232adcc780SMatthew G. Knepley   PetscFunctionBegin;
3249566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dms[0], &comm));
3259796d432SMatthew G. Knepley   /* Pull out local and global sections */
3269566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(len, &Nfs, len, &sections, len, &sectionGlobals));
3272adcc780SMatthew G. Knepley   for (i = 0; i < len; ++i) {
3289566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dms[i], &sections[i]));
3299566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalSection(dms[i], &sectionGlobals[i]));
3307a8be351SBarry Smith     PetscCheck(sections[i], comm, PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields");
3317a8be351SBarry Smith     PetscCheck(sectionGlobals[i], comm, PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields");
3329566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetNumFields(sections[i], &Nfs[i]));
3332adcc780SMatthew G. Knepley     Nf += Nfs[i];
3342adcc780SMatthew G. Knepley   }
3359796d432SMatthew G. Knepley   /* Create the supersection */
3369566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreateSupersection(sections, len, &supersection));
3379566063dSJacob Faibussowitsch   PetscCall(DMSetLocalSection(*superdm, supersection));
3389796d432SMatthew G. Knepley   /* Create ISes */
3392adcc780SMatthew G. Knepley   if (is) {
3409796d432SMatthew G. Knepley     PetscSection supersectionGlobal;
3419796d432SMatthew G. Knepley     PetscInt     bs = -1, startf = 0;
3422adcc780SMatthew G. Knepley 
3439566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(len, is));
3449566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalSection(*superdm, &supersectionGlobal));
3459796d432SMatthew G. Knepley     for (i = 0; i < len; startf += Nfs[i], ++i) {
3469796d432SMatthew G. Knepley       PetscInt *subIndices;
347ec4c761aSMatthew G. Knepley       PetscInt  subSize, subOff, pStart, pEnd, p, start, end, dummy;
3482adcc780SMatthew G. Knepley 
3499566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetChart(sectionGlobals[i], &pStart, &pEnd));
3509566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetConstrainedStorageSize(sectionGlobals[i], &subSize));
3519566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(subSize, &subIndices));
3529796d432SMatthew G. Knepley       for (p = pStart, subOff = 0; p < pEnd; ++p) {
3539796d432SMatthew G. Knepley         PetscInt gdof, gcdof, gtdof, d;
3542adcc780SMatthew G. Knepley 
3559566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(sectionGlobals[i], p, &gdof));
3569566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetConstraintDof(sections[i], p, &gcdof));
3572adcc780SMatthew G. Knepley         gtdof = gdof - gcdof;
3582adcc780SMatthew G. Knepley         if (gdof > 0 && gtdof) {
3599371c9d4SSatish Balay           if (bs < 0) {
3609371c9d4SSatish Balay             bs = gtdof;
3619371c9d4SSatish Balay           } else if (bs != gtdof) {
3629371c9d4SSatish Balay             bs = 1;
3639371c9d4SSatish Balay           }
3649566063dSJacob Faibussowitsch           PetscCall(DMGetGlobalFieldOffset_Private(*superdm, p, startf, &start, &dummy));
3659566063dSJacob Faibussowitsch           PetscCall(DMGetGlobalFieldOffset_Private(*superdm, p, startf + Nfs[i] - 1, &dummy, &end));
36663a3b9bcSJacob 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);
3679796d432SMatthew G. Knepley           for (d = start; d < end; ++d, ++subOff) subIndices[subOff] = d;
3682adcc780SMatthew G. Knepley         }
3692adcc780SMatthew G. Knepley       }
3709566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(comm, subSize, subIndices, PETSC_OWN_POINTER, &(*is)[i]));
3712adcc780SMatthew G. Knepley       /* Must have same blocksize on all procs (some might have no points) */
3729796d432SMatthew G. Knepley       {
3739796d432SMatthew G. Knepley         PetscInt bs = -1, bsLocal[2], bsMinMax[2];
3749796d432SMatthew G. Knepley 
3759371c9d4SSatish Balay         bsLocal[0] = bs < 0 ? PETSC_MAX_INT : bs;
3769371c9d4SSatish Balay         bsLocal[1] = bs;
3779566063dSJacob Faibussowitsch         PetscCall(PetscGlobalMinMaxInt(comm, bsLocal, bsMinMax));
3789371c9d4SSatish Balay         if (bsMinMax[0] != bsMinMax[1]) {
3799371c9d4SSatish Balay           bs = 1;
3809371c9d4SSatish Balay         } else {
3819371c9d4SSatish Balay           bs = bsMinMax[0];
3829371c9d4SSatish Balay         }
3839566063dSJacob Faibussowitsch         PetscCall(ISSetBlockSize((*is)[i], bs));
3842adcc780SMatthew G. Knepley       }
3852adcc780SMatthew G. Knepley     }
3862adcc780SMatthew G. Knepley   }
3879796d432SMatthew G. Knepley   /* Preserve discretizations */
388e5e52638SMatthew G. Knepley   if (len && dms[0]->probs) {
3899566063dSJacob Faibussowitsch     PetscCall(DMSetNumFields(*superdm, Nf));
3909796d432SMatthew G. Knepley     for (i = 0, supf = 0; i < len; ++i) {
3919796d432SMatthew G. Knepley       for (f = 0; f < Nfs[i]; ++f, ++supf) {
3922adcc780SMatthew G. Knepley         PetscObject disc;
3932adcc780SMatthew G. Knepley 
3949566063dSJacob Faibussowitsch         PetscCall(DMGetField(dms[i], f, NULL, &disc));
3959566063dSJacob Faibussowitsch         PetscCall(DMSetField(*superdm, supf, NULL, disc));
3962adcc780SMatthew G. Knepley       }
3972adcc780SMatthew G. Knepley     }
3989566063dSJacob Faibussowitsch     PetscCall(DMCreateDS(*superdm));
3992adcc780SMatthew G. Knepley   }
4009796d432SMatthew G. Knepley   /* Preserve nullspaces */
4019796d432SMatthew G. Knepley   for (i = 0, supf = 0; i < len; ++i) {
4029796d432SMatthew G. Knepley     for (f = 0; f < Nfs[i]; ++f, ++supf) {
4039796d432SMatthew G. Knepley       (*superdm)->nullspaceConstructors[supf] = dms[i]->nullspaceConstructors[f];
4049796d432SMatthew G. Knepley       if ((*superdm)->nullspaceConstructors[supf]) {
4059796d432SMatthew G. Knepley         haveNull = PETSC_TRUE;
4069796d432SMatthew G. Knepley         nullf    = supf;
4078cda7954SMatthew G. Knepley         oldf     = f;
4082adcc780SMatthew G. Knepley       }
4099796d432SMatthew G. Knepley     }
4109796d432SMatthew G. Knepley   }
4119796d432SMatthew G. Knepley   /* Attach nullspace to IS */
4129796d432SMatthew G. Knepley   if (haveNull && is) {
4139796d432SMatthew G. Knepley     MatNullSpace nullSpace;
4149796d432SMatthew G. Knepley 
4159566063dSJacob Faibussowitsch     PetscCall((*(*superdm)->nullspaceConstructors[nullf])(*superdm, oldf, nullf, &nullSpace));
4169566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)(*is)[nullf], "nullspace", (PetscObject)nullSpace));
4179566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceDestroy(&nullSpace));
4189796d432SMatthew G. Knepley   }
4199566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&supersection));
4209566063dSJacob Faibussowitsch   PetscCall(PetscFree3(Nfs, sections, sectionGlobals));
4212adcc780SMatthew G. Knepley   PetscFunctionReturn(0);
4222adcc780SMatthew G. Knepley }
423