xref: /petsc/src/dm/interface/dmi.c (revision 4758e3ce3f4cc9adc2a49a02e2f05c2a0f943969)
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 
401690c2aeSBarry Smith   // You cannot negate PETSC_INT_MIN
411690c2aeSBarry Smith   in[0] = blockSize < 0 ? -PETSC_INT_MAX : -blockSize;
425d1c0279SHong Zhang   in[1] = blockSize;
43462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(in, out, 2, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm)));
445d1c0279SHong Zhang   /* -out[0] = min(blockSize), out[1] = max(blockSize) */
455d1c0279SHong Zhang   if (-out[0] == out[1]) {
465d1c0279SHong Zhang     bs = out[1];
475d1c0279SHong Zhang   } else bs = 1;
485d1c0279SHong Zhang 
495d1c0279SHong Zhang   if (bs < 0) { /* Everyone was empty */
505d1c0279SHong Zhang     blockSize = 1;
515d1c0279SHong Zhang     bs        = 1;
525d1c0279SHong Zhang   }
535d1c0279SHong Zhang 
549566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetConstrainedStorageSize(gSection, &localSize));
557a8be351SBarry 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);
569566063dSJacob Faibussowitsch   PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), vec));
579566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(*vec, localSize, PETSC_DETERMINE));
589566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(*vec, bs));
599566063dSJacob Faibussowitsch   PetscCall(VecSetType(*vec, dm->vectype));
609566063dSJacob Faibussowitsch   PetscCall(VecSetDM(*vec, dm));
619566063dSJacob Faibussowitsch   /* PetscCall(VecSetLocalToGlobalMapping(*vec, dm->ltogmap)); */
623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6311689aebSJed Brown }
6411689aebSJed Brown 
65d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateLocalVector_Section_Private(DM dm, Vec *vec)
66d71ae5a4SJacob Faibussowitsch {
67fad22124SMatthew G Knepley   PetscSection section;
6811689aebSJed Brown   PetscInt     localSize, blockSize = -1, pStart, pEnd, p;
6911689aebSJed Brown 
7011689aebSJed Brown   PetscFunctionBegin;
719566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
729566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
7311689aebSJed Brown   for (p = pStart; p < pEnd; ++p) {
7411689aebSJed Brown     PetscInt dof;
7511689aebSJed Brown 
769566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
7711689aebSJed Brown     if ((blockSize < 0) && (dof > 0)) blockSize = dof;
785b8243e1SJed Brown     if (dof > 0) blockSize = PetscGCD(dof, blockSize);
7911689aebSJed Brown   }
809566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(section, &localSize));
819566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF, vec));
829566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(*vec, localSize, localSize));
8358b7e2c1SStefano Zampini   PetscCall(VecSetBlockSize(*vec, PetscAbs(blockSize)));
849566063dSJacob Faibussowitsch   PetscCall(VecSetType(*vec, dm->vectype));
859566063dSJacob Faibussowitsch   PetscCall(VecSetDM(*vec, dm));
863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8711689aebSJed Brown }
884d9407bcSMatthew G. Knepley 
89dd072f5fSMatthew G. Knepley static PetscErrorCode PetscSectionSelectFields_Private(PetscSection s, PetscSection gs, PetscInt numFields, const PetscInt fields[], const PetscInt numComps[], const PetscInt comps[], IS *is)
90d71ae5a4SJacob Faibussowitsch {
9155f9ef1dSMatthew G. Knepley   IS              permutation;
9255f9ef1dSMatthew G. Knepley   const PetscInt *perm = NULL;
934d9407bcSMatthew G. Knepley   PetscInt       *subIndices;
94a5e5a98fSStefano Zampini   PetscInt        mbs, bs = 0, bsLocal[2], bsMinMax[2];
95dd072f5fSMatthew G. Knepley   PetscInt        pStart, pEnd, Nc, subSize = 0, subOff = 0;
964d9407bcSMatthew G. Knepley 
974d9407bcSMatthew G. Knepley   PetscFunctionBegin;
989e47a1b0SMatthew G. Knepley   PetscCall(PetscSectionGetChart(gs, &pStart, &pEnd));
9955f9ef1dSMatthew G. Knepley   PetscCall(PetscSectionGetPermutation(s, &permutation));
10055f9ef1dSMatthew G. Knepley   if (permutation) PetscCall(ISGetIndices(permutation, &perm));
101dd072f5fSMatthew G. Knepley   if (numComps) {
102dd072f5fSMatthew G. Knepley     for (PetscInt f = 0, off = 0; f < numFields; ++f) {
1033a544194SStefano Zampini       PetscInt Nc;
1043a544194SStefano Zampini 
105dd072f5fSMatthew G. Knepley       if (numComps[f] < 0) continue;
106dd072f5fSMatthew G. Knepley       PetscCall(PetscSectionGetFieldComponents(s, f, &Nc));
107dd072f5fSMatthew G. Knepley       PetscCheck(numComps[f] <= Nc, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field %" PetscInt_FMT ": Number of selected components %" PetscInt_FMT " > %" PetscInt_FMT " number of field components", f, numComps[f], Nc);
108dd072f5fSMatthew G. Knepley       for (PetscInt c = 0; c < numComps[f]; ++c, ++off) PetscCheck(comps[off] < Nc, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field %" PetscInt_FMT ": component %" PetscInt_FMT " not in [0, %" PetscInt_FMT ")", f, comps[off], Nc);
109dd072f5fSMatthew G. Knepley       bs += numComps[f];
1109e47a1b0SMatthew G. Knepley     }
111dd072f5fSMatthew G. Knepley   } else {
1129e47a1b0SMatthew G. Knepley     for (PetscInt f = 0; f < numFields; ++f) {
1139e47a1b0SMatthew G. Knepley       PetscInt Nc;
1149e47a1b0SMatthew G. Knepley 
1159e47a1b0SMatthew G. Knepley       PetscCall(PetscSectionGetFieldComponents(s, fields[f], &Nc));
1163a544194SStefano Zampini       bs += Nc;
1173a544194SStefano Zampini     }
118dd072f5fSMatthew G. Knepley   }
119a5e5a98fSStefano Zampini   mbs = -1; /* multiple of block size not set */
1209e47a1b0SMatthew G. Knepley   for (PetscInt p = pStart; p < pEnd; ++p) {
12155f9ef1dSMatthew G. Knepley     const PetscInt point = perm ? perm[p - pStart] : p;
122b9eca265Ssarens     PetscInt       gdof, pSubSize = 0;
1234d9407bcSMatthew G. Knepley 
12455f9ef1dSMatthew G. Knepley     PetscCall(PetscSectionGetDof(gs, point, &gdof));
1254d9407bcSMatthew G. Knepley     if (gdof > 0) {
126dd072f5fSMatthew G. Knepley       PetscInt off = 0;
1274d9407bcSMatthew G. Knepley 
128dd072f5fSMatthew G. Knepley       for (PetscInt f = 0; f < numFields; ++f) {
129dd072f5fSMatthew G. Knepley         PetscInt fdof, fcdof, sfdof, sfcdof = 0;
130dd072f5fSMatthew G. Knepley 
131dd072f5fSMatthew G. Knepley         PetscCall(PetscSectionGetFieldComponents(s, f, &Nc));
13255f9ef1dSMatthew G. Knepley         PetscCall(PetscSectionGetFieldDof(s, point, fields[f], &fdof));
13355f9ef1dSMatthew G. Knepley         PetscCall(PetscSectionGetFieldConstraintDof(s, point, fields[f], &fcdof));
134dd072f5fSMatthew G. Knepley         if (numComps && numComps[f] >= 0) {
135dd072f5fSMatthew G. Knepley           const PetscInt *ind;
136dd072f5fSMatthew G. Knepley 
137dd072f5fSMatthew G. Knepley           // Assume sets of dofs on points are of size Nc
13855f9ef1dSMatthew G. Knepley           PetscCheck(!(fdof % Nc), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components %" PetscInt_FMT " should evenly divide the dofs %" PetscInt_FMT " on point %" PetscInt_FMT, Nc, fdof, point);
139dd072f5fSMatthew G. Knepley           sfdof = (fdof / Nc) * numComps[f];
14055f9ef1dSMatthew G. Knepley           PetscCall(PetscSectionGetFieldConstraintIndices(s, point, fields[f], &ind));
141dd072f5fSMatthew G. Knepley           for (PetscInt i = 0; i < (fdof / Nc); ++i) {
142dd072f5fSMatthew G. Knepley             for (PetscInt c = 0, fcc = 0; c < Nc; ++c) {
143dd072f5fSMatthew G. Knepley               if (c == comps[off + fcc]) {
144dd072f5fSMatthew G. Knepley                 ++fcc;
145dd072f5fSMatthew G. Knepley                 ++sfcdof;
146dd072f5fSMatthew G. Knepley               }
147dd072f5fSMatthew G. Knepley             }
148dd072f5fSMatthew G. Knepley           }
149dd072f5fSMatthew G. Knepley           pSubSize += sfdof - sfcdof;
150dd072f5fSMatthew G. Knepley           off += numComps[f];
151dd072f5fSMatthew G. Knepley         } else {
152b9eca265Ssarens           pSubSize += fdof - fcdof;
153b9eca265Ssarens         }
154dd072f5fSMatthew G. Knepley       }
155b9eca265Ssarens       subSize += pSubSize;
156a5e5a98fSStefano Zampini       if (pSubSize && pSubSize % bs) {
157a5e5a98fSStefano Zampini         // Layout does not admit a pointwise block size -> set mbs to 0
158a5e5a98fSStefano Zampini         mbs = 0;
159a5e5a98fSStefano Zampini       } else if (pSubSize) {
160a5e5a98fSStefano Zampini         if (mbs == -1) mbs = pSubSize / bs;
161a5e5a98fSStefano Zampini         else mbs = PetscMin(mbs, pSubSize / bs);
1624d9407bcSMatthew G. Knepley       }
1634d9407bcSMatthew G. Knepley     }
1644d9407bcSMatthew G. Knepley   }
165a5e5a98fSStefano Zampini 
1669e47a1b0SMatthew G. Knepley   // Must have same blocksize on all procs (some might have no points)
167a5e5a98fSStefano Zampini   bsLocal[0] = mbs < 0 ? PETSC_INT_MAX : mbs;
168a5e5a98fSStefano Zampini   bsLocal[1] = mbs;
1699e47a1b0SMatthew G. Knepley   PetscCall(PetscGlobalMinMaxInt(PetscObjectComm((PetscObject)gs), bsLocal, bsMinMax));
170a5e5a98fSStefano Zampini   if (bsMinMax[0] != bsMinMax[1]) { /* different multiple of block size -> set bs to 1 */
1719371c9d4SSatish Balay     bs = 1;
172a5e5a98fSStefano Zampini   } else { /* same multiple */
173a5e5a98fSStefano Zampini     mbs = bsMinMax[0];
174a5e5a98fSStefano Zampini     bs *= mbs;
1759371c9d4SSatish Balay   }
1769566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(subSize, &subIndices));
1779e47a1b0SMatthew G. Knepley   for (PetscInt p = pStart; p < pEnd; ++p) {
17855f9ef1dSMatthew G. Knepley     const PetscInt point = perm ? perm[p - pStart] : p;
1794d9407bcSMatthew G. Knepley     PetscInt       gdof, goff;
1804d9407bcSMatthew G. Knepley 
18155f9ef1dSMatthew G. Knepley     PetscCall(PetscSectionGetDof(gs, point, &gdof));
1824d9407bcSMatthew G. Knepley     if (gdof > 0) {
183dd072f5fSMatthew G. Knepley       PetscInt off = 0;
184dd072f5fSMatthew G. Knepley 
18555f9ef1dSMatthew G. Knepley       PetscCall(PetscSectionGetOffset(gs, point, &goff));
1869e47a1b0SMatthew G. Knepley       for (PetscInt f = 0; f < numFields; ++f) {
1879e47a1b0SMatthew G. Knepley         PetscInt fdof, fcdof, poff = 0;
1884d9407bcSMatthew G. Knepley 
1894d9407bcSMatthew G. Knepley         /* Can get rid of this loop by storing field information in the global section */
1909e47a1b0SMatthew G. Knepley         for (PetscInt f2 = 0; f2 < fields[f]; ++f2) {
19155f9ef1dSMatthew G. Knepley           PetscCall(PetscSectionGetFieldDof(s, point, f2, &fdof));
19255f9ef1dSMatthew G. Knepley           PetscCall(PetscSectionGetFieldConstraintDof(s, point, f2, &fcdof));
1934d9407bcSMatthew G. Knepley           poff += fdof - fcdof;
1944d9407bcSMatthew G. Knepley         }
19555f9ef1dSMatthew G. Knepley         PetscCall(PetscSectionGetFieldDof(s, point, fields[f], &fdof));
19655f9ef1dSMatthew G. Knepley         PetscCall(PetscSectionGetFieldConstraintDof(s, point, fields[f], &fcdof));
197dd072f5fSMatthew G. Knepley 
198dd072f5fSMatthew G. Knepley         if (numComps && numComps[f] >= 0) {
199dd072f5fSMatthew G. Knepley           const PetscInt *ind;
200dd072f5fSMatthew G. Knepley 
201dd072f5fSMatthew G. Knepley           // Assume sets of dofs on points are of size Nc
20255f9ef1dSMatthew G. Knepley           PetscCall(PetscSectionGetFieldConstraintIndices(s, point, fields[f], &ind));
203dd072f5fSMatthew G. Knepley           for (PetscInt i = 0, fcoff = 0, pfoff = 0; i < (fdof / Nc); ++i) {
204dd072f5fSMatthew G. Knepley             for (PetscInt c = 0, fcc = 0; c < Nc; ++c) {
205dd072f5fSMatthew G. Knepley               const PetscInt k = i * Nc + c;
206dd072f5fSMatthew G. Knepley 
207dd072f5fSMatthew G. Knepley               if (ind[fcoff] == k) {
208dd072f5fSMatthew G. Knepley                 ++fcoff;
209dd072f5fSMatthew G. Knepley                 continue;
210dd072f5fSMatthew G. Knepley               }
211dd072f5fSMatthew G. Knepley               if (c == comps[off + fcc]) {
212dd072f5fSMatthew G. Knepley                 ++fcc;
213dd072f5fSMatthew G. Knepley                 subIndices[subOff++] = goff + poff + pfoff;
214dd072f5fSMatthew G. Knepley               }
215dd072f5fSMatthew G. Knepley               ++pfoff;
216dd072f5fSMatthew G. Knepley             }
217dd072f5fSMatthew G. Knepley           }
218dd072f5fSMatthew G. Knepley           off += numComps[f];
219dd072f5fSMatthew G. Knepley         } else {
2209e47a1b0SMatthew G. Knepley           for (PetscInt fc = 0; fc < fdof - fcdof; ++fc, ++subOff) subIndices[subOff] = goff + poff + fc;
2214d9407bcSMatthew G. Knepley         }
2224d9407bcSMatthew G. Knepley       }
2234d9407bcSMatthew G. Knepley     }
224dd072f5fSMatthew G. Knepley   }
22555f9ef1dSMatthew G. Knepley   if (permutation) PetscCall(ISRestoreIndices(permutation, &perm));
2269e47a1b0SMatthew G. Knepley   PetscCheck(subSize == subOff, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "The offset array size %" PetscInt_FMT " != %" PetscInt_FMT " the number of indices", subSize, subOff);
2279e47a1b0SMatthew G. Knepley   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)gs), subSize, subIndices, PETSC_OWN_POINTER, is));
228627349f5SMatthew G. Knepley   if (bs > 1) {
2299e47a1b0SMatthew G. Knepley     // We need to check that the block size does not come from non-contiguous fields
2309e47a1b0SMatthew G. Knepley     PetscInt set = 1, rset = 1;
2319e47a1b0SMatthew G. Knepley     for (PetscInt i = 0; i < subSize; i += bs) {
2329e47a1b0SMatthew G. Knepley       for (PetscInt j = 0; j < bs; ++j) {
2339371c9d4SSatish Balay         if (subIndices[i + j] != subIndices[i] + j) {
2349371c9d4SSatish Balay           set = 0;
2359371c9d4SSatish Balay           break;
2369371c9d4SSatish Balay         }
237627349f5SMatthew G. Knepley       }
238627349f5SMatthew G. Knepley     }
239462c564dSBarry Smith     PetscCallMPI(MPIU_Allreduce(&set, &rset, 1, MPIU_INT, MPI_PROD, PetscObjectComm((PetscObject)gs)));
24088206212SAlexis Marboeuf     if (rset) PetscCall(ISSetBlockSize(*is, bs));
241627349f5SMatthew G. Knepley   }
2429e47a1b0SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
2434d9407bcSMatthew G. Knepley }
2449e47a1b0SMatthew G. Knepley 
245dd072f5fSMatthew G. Knepley static PetscErrorCode DMSelectFields_Private(DM dm, PetscSection section, PetscInt numFields, const PetscInt fields[], const PetscInt numComps[], const PetscInt comps[], IS *is, DM *subdm)
2469e47a1b0SMatthew G. Knepley {
2474d9407bcSMatthew G. Knepley   PetscSection subsection;
2484d9407bcSMatthew G. Knepley   PetscBool    haveNull = PETSC_FALSE;
2499e47a1b0SMatthew G. Knepley   PetscInt     nf = 0, of = 0;
2504d9407bcSMatthew G. Knepley 
2519e47a1b0SMatthew G. Knepley   PetscFunctionBegin;
252*4758e3ceSMatthew G. Knepley   // Create nullspace constructor slots
253*4758e3ceSMatthew G. Knepley   if (dm->nullspaceConstructors) {
254*4758e3ceSMatthew G. Knepley     PetscCall(PetscFree2((*subdm)->nullspaceConstructors, (*subdm)->nearnullspaceConstructors));
255*4758e3ceSMatthew G. Knepley     PetscCall(PetscCalloc2(numFields, &(*subdm)->nullspaceConstructors, numFields, &(*subdm)->nearnullspaceConstructors));
256*4758e3ceSMatthew G. Knepley   }
257dd072f5fSMatthew G. Knepley   if (numComps) {
258dd072f5fSMatthew G. Knepley     const PetscInt field = fields[0];
259dd072f5fSMatthew G. Knepley 
260dd072f5fSMatthew G. Knepley     PetscCheck(numFields == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "We only support a single field for component selection right now");
261dd072f5fSMatthew G. Knepley     PetscCall(PetscSectionCreateComponentSubsection(section, numComps[field], comps, &subsection));
2629566063dSJacob Faibussowitsch     PetscCall(DMSetLocalSection(*subdm, subsection));
2639566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&subsection));
264*4758e3ceSMatthew G. Knepley     if (dm->nullspaceConstructors) (*subdm)->nullspaceConstructors[field] = dm->nullspaceConstructors[field];
265e5e52638SMatthew G. Knepley     if (dm->probs) {
2669e47a1b0SMatthew G. Knepley       PetscFV  fv, fvNew;
267dd072f5fSMatthew G. Knepley       PetscInt fnum[1] = {field};
2680f21e855SMatthew G. Knepley 
2699e47a1b0SMatthew G. Knepley       PetscCall(DMSetNumFields(*subdm, 1));
270dd072f5fSMatthew G. Knepley       PetscCall(DMGetField(dm, field, NULL, (PetscObject *)&fv));
2719e47a1b0SMatthew G. Knepley       PetscCall(PetscFVClone(fv, &fvNew));
272dd072f5fSMatthew G. Knepley       PetscCall(PetscFVSetNumComponents(fvNew, numComps[0]));
2739e47a1b0SMatthew G. Knepley       PetscCall(DMSetField(*subdm, 0, NULL, (PetscObject)fvNew));
2749e47a1b0SMatthew G. Knepley       PetscCall(PetscFVDestroy(&fvNew));
2759566063dSJacob Faibussowitsch       PetscCall(DMCreateDS(*subdm));
276dd072f5fSMatthew G. Knepley       if (numComps[0] == 1 && is) {
2770f21e855SMatthew G. Knepley         PetscObject disc, space, pmat;
2784d9407bcSMatthew G. Knepley 
279dd072f5fSMatthew G. Knepley         PetscCall(DMGetField(*subdm, field, NULL, &disc));
2809566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "nullspace", &space));
2819566063dSJacob Faibussowitsch         if (space) PetscCall(PetscObjectCompose((PetscObject)*is, "nullspace", space));
2829566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "nearnullspace", &space));
2839566063dSJacob Faibussowitsch         if (space) PetscCall(PetscObjectCompose((PetscObject)*is, "nearnullspace", space));
2849566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "pmat", &pmat));
2859566063dSJacob Faibussowitsch         if (pmat) PetscCall(PetscObjectCompose((PetscObject)*is, "pmat", pmat));
2864d9407bcSMatthew G. Knepley       }
287dd072f5fSMatthew G. Knepley       PetscCall(PetscDSCopyConstants(dm->probs[field].ds, (*subdm)->probs[0].ds));
288dd072f5fSMatthew G. Knepley       PetscCall(PetscDSCopyBoundary(dm->probs[field].ds, 1, fnum, (*subdm)->probs[0].ds));
289dd072f5fSMatthew G. Knepley       PetscCall(PetscDSSelectEquations(dm->probs[field].ds, 1, fnum, (*subdm)->probs[0].ds));
2909e47a1b0SMatthew G. Knepley     }
291*4758e3ceSMatthew G. Knepley     if ((*subdm)->nullspaceConstructors && (*subdm)->nullspaceConstructors[0] && is) {
2929e47a1b0SMatthew G. Knepley       MatNullSpace nullSpace;
2939e47a1b0SMatthew G. Knepley 
2949e47a1b0SMatthew G. Knepley       PetscCall((*(*subdm)->nullspaceConstructors[0])(*subdm, 0, 0, &nullSpace));
2959e47a1b0SMatthew G. Knepley       PetscCall(PetscObjectCompose((PetscObject)*is, "nullspace", (PetscObject)nullSpace));
2969e47a1b0SMatthew G. Knepley       PetscCall(MatNullSpaceDestroy(&nullSpace));
2979e47a1b0SMatthew G. Knepley     }
2989e47a1b0SMatthew G. Knepley     if (dm->coarseMesh) PetscCall(DMCreateSubDM(dm->coarseMesh, numFields, fields, NULL, &(*subdm)->coarseMesh));
2999e47a1b0SMatthew G. Knepley     PetscFunctionReturn(PETSC_SUCCESS);
3009e47a1b0SMatthew G. Knepley   }
3019e47a1b0SMatthew G. Knepley   PetscCall(PetscSectionCreateSubsection(section, numFields, fields, &subsection));
3029e47a1b0SMatthew G. Knepley   PetscCall(DMSetLocalSection(*subdm, subsection));
3039e47a1b0SMatthew G. Knepley   PetscCall(PetscSectionDestroy(&subsection));
3049e47a1b0SMatthew G. Knepley   if (dm->probs) {
3059e47a1b0SMatthew G. Knepley     PetscCall(DMSetNumFields(*subdm, numFields));
3069e47a1b0SMatthew G. Knepley     for (PetscInt f = 0; f < numFields; ++f) {
3079e47a1b0SMatthew G. Knepley       PetscObject disc;
3089e47a1b0SMatthew G. Knepley 
3099e47a1b0SMatthew G. Knepley       PetscCall(DMGetField(dm, fields[f], NULL, &disc));
3109e47a1b0SMatthew G. Knepley       PetscCall(DMSetField(*subdm, f, NULL, disc));
3119e47a1b0SMatthew G. Knepley     }
3129e47a1b0SMatthew G. Knepley     // TODO: if only FV, then cut down the components
3139e47a1b0SMatthew G. Knepley     PetscCall(DMCreateDS(*subdm));
3149e47a1b0SMatthew G. Knepley     if (numFields == 1 && is) {
3159e47a1b0SMatthew G. Knepley       PetscObject disc, space, pmat;
3169e47a1b0SMatthew G. Knepley 
3179e47a1b0SMatthew G. Knepley       PetscCall(DMGetField(*subdm, 0, NULL, &disc));
3189e47a1b0SMatthew G. Knepley       PetscCall(PetscObjectQuery(disc, "nullspace", &space));
3199e47a1b0SMatthew G. Knepley       if (space) PetscCall(PetscObjectCompose((PetscObject)*is, "nullspace", space));
3209e47a1b0SMatthew G. Knepley       PetscCall(PetscObjectQuery(disc, "nearnullspace", &space));
3219e47a1b0SMatthew G. Knepley       if (space) PetscCall(PetscObjectCompose((PetscObject)*is, "nearnullspace", space));
3229e47a1b0SMatthew G. Knepley       PetscCall(PetscObjectQuery(disc, "pmat", &pmat));
3239e47a1b0SMatthew G. Knepley       if (pmat) PetscCall(PetscObjectCompose((PetscObject)*is, "pmat", pmat));
3249e47a1b0SMatthew G. Knepley     }
3259e47a1b0SMatthew G. Knepley     // Check if DSes record their DM fields
326e21d936aSMatthew G. Knepley     if (dm->probs[0].fields) {
3279310035eSMatthew G. Knepley       PetscInt d, e;
3289310035eSMatthew G. Knepley 
3299310035eSMatthew G. Knepley       for (d = 0, e = 0; d < dm->Nds && e < (*subdm)->Nds; ++d) {
3309310035eSMatthew G. Knepley         const PetscInt  Nf = dm->probs[d].ds->Nf;
3319310035eSMatthew G. Knepley         const PetscInt *fld;
3329310035eSMatthew G. Knepley         PetscInt        f, g;
3339310035eSMatthew G. Knepley 
3349566063dSJacob Faibussowitsch         PetscCall(ISGetIndices(dm->probs[d].fields, &fld));
3359310035eSMatthew G. Knepley         for (f = 0; f < Nf; ++f) {
3369371c9d4SSatish Balay           for (g = 0; g < numFields; ++g)
3379371c9d4SSatish Balay             if (fld[f] == fields[g]) break;
3389310035eSMatthew G. Knepley           if (g < numFields) break;
3399310035eSMatthew G. Knepley         }
3409566063dSJacob Faibussowitsch         PetscCall(ISRestoreIndices(dm->probs[d].fields, &fld));
3419310035eSMatthew G. Knepley         if (f == Nf) continue;
3429566063dSJacob Faibussowitsch         PetscCall(PetscDSCopyConstants(dm->probs[d].ds, (*subdm)->probs[e].ds));
3439566063dSJacob Faibussowitsch         PetscCall(PetscDSCopyBoundary(dm->probs[d].ds, numFields, fields, (*subdm)->probs[e].ds));
3449e47a1b0SMatthew G. Knepley         // Translate DM fields to DS fields
3459310035eSMatthew G. Knepley         {
346e21d936aSMatthew G. Knepley           IS              infields, dsfields;
34774a61aaaSMatthew G. Knepley           const PetscInt *fld, *ofld;
34874a61aaaSMatthew G. Knepley           PetscInt       *fidx;
3499e47a1b0SMatthew G. Knepley           PetscInt        onf, nf;
350e21d936aSMatthew G. Knepley 
3519566063dSJacob Faibussowitsch           PetscCall(ISCreateGeneral(PETSC_COMM_SELF, numFields, fields, PETSC_USE_POINTER, &infields));
3529566063dSJacob Faibussowitsch           PetscCall(ISIntersect(infields, dm->probs[d].fields, &dsfields));
3539566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&infields));
3549566063dSJacob Faibussowitsch           PetscCall(ISGetLocalSize(dsfields, &nf));
3557a8be351SBarry Smith           PetscCheck(nf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "DS cannot be supported on 0 fields");
3569566063dSJacob Faibussowitsch           PetscCall(ISGetIndices(dsfields, &fld));
3579566063dSJacob Faibussowitsch           PetscCall(ISGetLocalSize(dm->probs[d].fields, &onf));
3589566063dSJacob Faibussowitsch           PetscCall(ISGetIndices(dm->probs[d].fields, &ofld));
3599566063dSJacob Faibussowitsch           PetscCall(PetscMalloc1(nf, &fidx));
3609e47a1b0SMatthew G. Knepley           for (PetscInt f = 0, g = 0; f < onf && g < nf; ++f) {
36174a61aaaSMatthew G. Knepley             if (ofld[f] == fld[g]) fidx[g++] = f;
36274a61aaaSMatthew G. Knepley           }
3639566063dSJacob Faibussowitsch           PetscCall(ISRestoreIndices(dm->probs[d].fields, &ofld));
3649566063dSJacob Faibussowitsch           PetscCall(ISRestoreIndices(dsfields, &fld));
3659566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&dsfields));
366bb4b53efSMatthew G. Knepley           PetscCall(PetscDSSelectDiscretizations(dm->probs[0].ds, nf, fidx, PETSC_DETERMINE, PETSC_DETERMINE, (*subdm)->probs[0].ds));
3679566063dSJacob Faibussowitsch           PetscCall(PetscDSSelectEquations(dm->probs[0].ds, nf, fidx, (*subdm)->probs[0].ds));
3689566063dSJacob Faibussowitsch           PetscCall(PetscFree(fidx));
3699310035eSMatthew G. Knepley         }
3709310035eSMatthew G. Knepley         ++e;
3719310035eSMatthew G. Knepley       }
372e21d936aSMatthew G. Knepley     } else {
3739566063dSJacob Faibussowitsch       PetscCall(PetscDSCopyConstants(dm->probs[0].ds, (*subdm)->probs[0].ds));
3749566063dSJacob Faibussowitsch       PetscCall(PetscDSCopyBoundary(dm->probs[0].ds, PETSC_DETERMINE, NULL, (*subdm)->probs[0].ds));
375bb4b53efSMatthew G. Knepley       PetscCall(PetscDSSelectDiscretizations(dm->probs[0].ds, numFields, fields, PETSC_DETERMINE, PETSC_DETERMINE, (*subdm)->probs[0].ds));
3769566063dSJacob Faibussowitsch       PetscCall(PetscDSSelectEquations(dm->probs[0].ds, numFields, fields, (*subdm)->probs[0].ds));
377d5af271aSMatthew G. Knepley     }
378e21d936aSMatthew G. Knepley   }
379*4758e3ceSMatthew G. Knepley   for (PetscInt f = 0; f < numFields; ++f) {
380*4758e3ceSMatthew G. Knepley     if (dm->nullspaceConstructors) {
381*4758e3ceSMatthew G. Knepley       (*subdm)->nullspaceConstructors[f] = dm->nullspaceConstructors[fields[f]];
382*4758e3ceSMatthew G. Knepley       if ((*subdm)->nullspaceConstructors[f]) {
383*4758e3ceSMatthew G. Knepley         haveNull = PETSC_TRUE;
384*4758e3ceSMatthew G. Knepley         nf       = f;
385*4758e3ceSMatthew G. Knepley         of       = fields[f];
386*4758e3ceSMatthew G. Knepley       }
387*4758e3ceSMatthew G. Knepley     }
388*4758e3ceSMatthew G. Knepley   }
3898cda7954SMatthew G. Knepley   if (haveNull && is) {
3908cda7954SMatthew G. Knepley     MatNullSpace nullSpace;
3918cda7954SMatthew G. Knepley 
3929566063dSJacob Faibussowitsch     PetscCall((*(*subdm)->nullspaceConstructors[nf])(*subdm, of, nf, &nullSpace));
3939566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)*is, "nullspace", (PetscObject)nullSpace));
3949566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceDestroy(&nullSpace));
3958cda7954SMatthew G. Knepley   }
39648a46eb9SPierre Jolivet   if (dm->coarseMesh) PetscCall(DMCreateSubDM(dm->coarseMesh, numFields, fields, NULL, &(*subdm)->coarseMesh));
3979e47a1b0SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
3984d9407bcSMatthew G. Knepley }
3999e47a1b0SMatthew G. Knepley 
4005d83a8b1SBarry Smith /*@
4019e47a1b0SMatthew G. Knepley   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`.
4029e47a1b0SMatthew G. Knepley 
4039e47a1b0SMatthew G. Knepley   Not Collective
4049e47a1b0SMatthew G. Knepley 
4059e47a1b0SMatthew G. Knepley   Input Parameters:
4069e47a1b0SMatthew G. Knepley + dm        - The `DM` object
4079e47a1b0SMatthew G. Knepley . numFields - The number of fields to incorporate into `subdm`
408dd072f5fSMatthew G. Knepley . fields    - The field numbers of the selected fields
409dd072f5fSMatthew G. Knepley . numComps  - The number of components from each field to incorporate into `subdm`, or PETSC_DECIDE for all components
410dd072f5fSMatthew G. Knepley - comps     - The component numbers of the selected fields (omitted for PTESC_DECIDE fields)
4119e47a1b0SMatthew G. Knepley 
4129e47a1b0SMatthew G. Knepley   Output Parameters:
4139e47a1b0SMatthew G. Knepley + is    - The global indices for the subproblem or `NULL`
4149e47a1b0SMatthew G. Knepley - subdm - The `DM` for the subproblem, which must already have be cloned from `dm` or `NULL`
4159e47a1b0SMatthew G. Knepley 
4169e47a1b0SMatthew G. Knepley   Level: intermediate
4179e47a1b0SMatthew G. Knepley 
4189e47a1b0SMatthew G. Knepley   Notes:
4199e47a1b0SMatthew G. Knepley   If `is` and `subdm` are both `NULL` this does nothing
4209e47a1b0SMatthew G. Knepley 
4219e47a1b0SMatthew G. Knepley .seealso: `DMCreateSubDM()`, `DMGetLocalSection()`, `DMPlexSetMigrationSF()`, `DMView()`
4229e47a1b0SMatthew G. Knepley @*/
423dd072f5fSMatthew G. Knepley PetscErrorCode DMCreateSectionSubDM(DM dm, PetscInt numFields, const PetscInt fields[], const PetscInt numComps[], const PetscInt comps[], IS *is, DM *subdm)
4249e47a1b0SMatthew G. Knepley {
4259e47a1b0SMatthew G. Knepley   PetscSection section, sectionGlobal;
4269e47a1b0SMatthew G. Knepley   PetscInt     Nf;
4279e47a1b0SMatthew G. Knepley 
4289e47a1b0SMatthew G. Knepley   PetscFunctionBegin;
4299e47a1b0SMatthew G. Knepley   if (!numFields) PetscFunctionReturn(PETSC_SUCCESS);
4309e47a1b0SMatthew G. Knepley   PetscCall(DMGetLocalSection(dm, &section));
4319e47a1b0SMatthew G. Knepley   PetscCall(DMGetGlobalSection(dm, &sectionGlobal));
4329e47a1b0SMatthew G. Knepley   PetscCheck(section, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields");
4339e47a1b0SMatthew G. Knepley   PetscCheck(sectionGlobal, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields");
4349e47a1b0SMatthew G. Knepley   PetscCall(PetscSectionGetNumFields(section, &Nf));
4359e47a1b0SMatthew G. Knepley   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);
4369e47a1b0SMatthew G. Knepley 
437dd072f5fSMatthew G. Knepley   if (is) PetscCall(PetscSectionSelectFields_Private(section, sectionGlobal, numFields, fields, numComps, comps, is));
438dd072f5fSMatthew G. Knepley   if (subdm) PetscCall(DMSelectFields_Private(dm, section, numFields, fields, numComps, comps, is, subdm));
4393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4404d9407bcSMatthew G. Knepley }
4412adcc780SMatthew G. Knepley 
442792b654fSMatthew G. Knepley /*@C
4430b3275a6SBarry 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`
444792b654fSMatthew G. Knepley 
44520f4b53cSBarry Smith   Not Collective
446792b654fSMatthew G. Knepley 
447d8d19677SJose E. Roman   Input Parameters:
4480b3275a6SBarry Smith + dms - The `DM` objects, the must all have the same topology; for example obtained with `DMClone()`
4490b3275a6SBarry Smith - len - The number of `DM` in `dms`
450792b654fSMatthew G. Knepley 
451792b654fSMatthew G. Knepley   Output Parameters:
45220f4b53cSBarry Smith + is      - The global indices for the subproblem, or `NULL`
4530b3275a6SBarry Smith - superdm - The `DM` for the superproblem, which must already have be cloned and contain the same topology as the `dms`
454792b654fSMatthew G. Knepley 
455792b654fSMatthew G. Knepley   Level: intermediate
456792b654fSMatthew G. Knepley 
457a4e35b19SJacob Faibussowitsch .seealso: `DMCreateSuperDM()`, `DMGetLocalSection()`, `DMPlexSetMigrationSF()`, `DMView()`
458792b654fSMatthew G. Knepley @*/
4595d83a8b1SBarry Smith PetscErrorCode DMCreateSectionSuperDM(DM dms[], PetscInt len, IS *is[], DM *superdm)
460d71ae5a4SJacob Faibussowitsch {
4619796d432SMatthew G. Knepley   MPI_Comm     comm;
4629796d432SMatthew G. Knepley   PetscSection supersection, *sections, *sectionGlobals;
4638cda7954SMatthew G. Knepley   PetscInt    *Nfs, Nf = 0, f, supf, oldf = -1, nullf = -1, i;
4649796d432SMatthew G. Knepley   PetscBool    haveNull = PETSC_FALSE;
4652adcc780SMatthew G. Knepley 
4662adcc780SMatthew G. Knepley   PetscFunctionBegin;
4679566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dms[0], &comm));
4689796d432SMatthew G. Knepley   /* Pull out local and global sections */
4699566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(len, &Nfs, len, &sections, len, &sectionGlobals));
4702adcc780SMatthew G. Knepley   for (i = 0; i < len; ++i) {
4719566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dms[i], &sections[i]));
4729566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalSection(dms[i], &sectionGlobals[i]));
4737a8be351SBarry Smith     PetscCheck(sections[i], comm, PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields");
4747a8be351SBarry Smith     PetscCheck(sectionGlobals[i], comm, PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields");
4759566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetNumFields(sections[i], &Nfs[i]));
4762adcc780SMatthew G. Knepley     Nf += Nfs[i];
4772adcc780SMatthew G. Knepley   }
4789796d432SMatthew G. Knepley   /* Create the supersection */
4799566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreateSupersection(sections, len, &supersection));
4809566063dSJacob Faibussowitsch   PetscCall(DMSetLocalSection(*superdm, supersection));
4819796d432SMatthew G. Knepley   /* Create ISes */
4822adcc780SMatthew G. Knepley   if (is) {
4839796d432SMatthew G. Knepley     PetscSection supersectionGlobal;
4849796d432SMatthew G. Knepley     PetscInt     bs = -1, startf = 0;
4852adcc780SMatthew G. Knepley 
4869566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(len, is));
4879566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalSection(*superdm, &supersectionGlobal));
4889796d432SMatthew G. Knepley     for (i = 0; i < len; startf += Nfs[i], ++i) {
4899796d432SMatthew G. Knepley       PetscInt *subIndices;
490ec4c761aSMatthew G. Knepley       PetscInt  subSize, subOff, pStart, pEnd, p, start, end, dummy;
4912adcc780SMatthew G. Knepley 
4929566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetChart(sectionGlobals[i], &pStart, &pEnd));
4939566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetConstrainedStorageSize(sectionGlobals[i], &subSize));
4949566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(subSize, &subIndices));
4959796d432SMatthew G. Knepley       for (p = pStart, subOff = 0; p < pEnd; ++p) {
4969796d432SMatthew G. Knepley         PetscInt gdof, gcdof, gtdof, d;
4972adcc780SMatthew G. Knepley 
4989566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(sectionGlobals[i], p, &gdof));
4999566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetConstraintDof(sections[i], p, &gcdof));
5002adcc780SMatthew G. Knepley         gtdof = gdof - gcdof;
5012adcc780SMatthew G. Knepley         if (gdof > 0 && gtdof) {
5029371c9d4SSatish Balay           if (bs < 0) {
5039371c9d4SSatish Balay             bs = gtdof;
5049371c9d4SSatish Balay           } else if (bs != gtdof) {
5059371c9d4SSatish Balay             bs = 1;
5069371c9d4SSatish Balay           }
5079566063dSJacob Faibussowitsch           PetscCall(DMGetGlobalFieldOffset_Private(*superdm, p, startf, &start, &dummy));
5089566063dSJacob Faibussowitsch           PetscCall(DMGetGlobalFieldOffset_Private(*superdm, p, startf + Nfs[i] - 1, &dummy, &end));
50963a3b9bcSJacob 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);
5109796d432SMatthew G. Knepley           for (d = start; d < end; ++d, ++subOff) subIndices[subOff] = d;
5112adcc780SMatthew G. Knepley         }
5122adcc780SMatthew G. Knepley       }
5139566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(comm, subSize, subIndices, PETSC_OWN_POINTER, &(*is)[i]));
5142adcc780SMatthew G. Knepley       /* Must have same blocksize on all procs (some might have no points) */
5159796d432SMatthew G. Knepley       {
5169796d432SMatthew G. Knepley         PetscInt bs = -1, bsLocal[2], bsMinMax[2];
5179796d432SMatthew G. Knepley 
5181690c2aeSBarry Smith         bsLocal[0] = bs < 0 ? PETSC_INT_MAX : bs;
5199371c9d4SSatish Balay         bsLocal[1] = bs;
5209566063dSJacob Faibussowitsch         PetscCall(PetscGlobalMinMaxInt(comm, bsLocal, bsMinMax));
5219371c9d4SSatish Balay         if (bsMinMax[0] != bsMinMax[1]) {
5229371c9d4SSatish Balay           bs = 1;
5239371c9d4SSatish Balay         } else {
5249371c9d4SSatish Balay           bs = bsMinMax[0];
5259371c9d4SSatish Balay         }
5269566063dSJacob Faibussowitsch         PetscCall(ISSetBlockSize((*is)[i], bs));
5272adcc780SMatthew G. Knepley       }
5282adcc780SMatthew G. Knepley     }
5292adcc780SMatthew G. Knepley   }
5309796d432SMatthew G. Knepley   /* Preserve discretizations */
531e5e52638SMatthew G. Knepley   if (len && dms[0]->probs) {
5329566063dSJacob Faibussowitsch     PetscCall(DMSetNumFields(*superdm, Nf));
5339796d432SMatthew G. Knepley     for (i = 0, supf = 0; i < len; ++i) {
5349796d432SMatthew G. Knepley       for (f = 0; f < Nfs[i]; ++f, ++supf) {
5352adcc780SMatthew G. Knepley         PetscObject disc;
5362adcc780SMatthew G. Knepley 
5379566063dSJacob Faibussowitsch         PetscCall(DMGetField(dms[i], f, NULL, &disc));
5389566063dSJacob Faibussowitsch         PetscCall(DMSetField(*superdm, supf, NULL, disc));
5392adcc780SMatthew G. Knepley       }
5402adcc780SMatthew G. Knepley     }
5419566063dSJacob Faibussowitsch     PetscCall(DMCreateDS(*superdm));
5422adcc780SMatthew G. Knepley   }
543*4758e3ceSMatthew G. Knepley   // Create nullspace constructor slots
544*4758e3ceSMatthew G. Knepley   PetscCall(PetscFree2((*superdm)->nullspaceConstructors, (*superdm)->nearnullspaceConstructors));
545*4758e3ceSMatthew G. Knepley   PetscCall(PetscCalloc2(Nf, &(*superdm)->nullspaceConstructors, Nf, &(*superdm)->nearnullspaceConstructors));
5469796d432SMatthew G. Knepley   /* Preserve nullspaces */
5479796d432SMatthew G. Knepley   for (i = 0, supf = 0; i < len; ++i) {
5489796d432SMatthew G. Knepley     for (f = 0; f < Nfs[i]; ++f, ++supf) {
549*4758e3ceSMatthew G. Knepley       if (dms[i]->nullspaceConstructors) {
5509796d432SMatthew G. Knepley         (*superdm)->nullspaceConstructors[supf] = dms[i]->nullspaceConstructors[f];
5519796d432SMatthew G. Knepley         if ((*superdm)->nullspaceConstructors[supf]) {
5529796d432SMatthew G. Knepley           haveNull = PETSC_TRUE;
5539796d432SMatthew G. Knepley           nullf    = supf;
5548cda7954SMatthew G. Knepley           oldf     = f;
5552adcc780SMatthew G. Knepley         }
5569796d432SMatthew G. Knepley       }
5579796d432SMatthew G. Knepley     }
558*4758e3ceSMatthew G. Knepley   }
5599796d432SMatthew G. Knepley   /* Attach nullspace to IS */
5609796d432SMatthew G. Knepley   if (haveNull && is) {
5619796d432SMatthew G. Knepley     MatNullSpace nullSpace;
5629796d432SMatthew G. Knepley 
5639566063dSJacob Faibussowitsch     PetscCall((*(*superdm)->nullspaceConstructors[nullf])(*superdm, oldf, nullf, &nullSpace));
5649566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)(*is)[nullf], "nullspace", (PetscObject)nullSpace));
5659566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceDestroy(&nullSpace));
5669796d432SMatthew G. Knepley   }
5679566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&supersection));
5689566063dSJacob Faibussowitsch   PetscCall(PetscFree3(Nfs, sections, sectionGlobals));
5693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5702adcc780SMatthew G. Knepley }
571