xref: /petsc/src/dm/interface/dmi.c (revision dd072f5f68c35d6a8c33edba0c786d5785b10330)
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 
88*dd072f5fSMatthew G. Knepley static PetscErrorCode PetscSectionSelectFields_Private(PetscSection s, PetscSection gs, PetscInt numFields, const PetscInt fields[], const PetscInt numComps[], const PetscInt comps[], IS *is)
89d71ae5a4SJacob Faibussowitsch {
904d9407bcSMatthew G. Knepley   PetscInt *subIndices;
91*dd072f5fSMatthew G. Knepley   PetscInt  bs = 0, bsLocal[2], bsMinMax[2];
92*dd072f5fSMatthew G. Knepley   PetscInt  pStart, pEnd, Nc, subSize = 0, subOff = 0;
934d9407bcSMatthew G. Knepley 
944d9407bcSMatthew G. Knepley   PetscFunctionBegin;
959e47a1b0SMatthew G. Knepley   PetscCall(PetscSectionGetChart(gs, &pStart, &pEnd));
96*dd072f5fSMatthew G. Knepley   if (numComps) {
97*dd072f5fSMatthew G. Knepley     for (PetscInt f = 0, off = 0; f < numFields; ++f) {
983a544194SStefano Zampini       PetscInt Nc;
993a544194SStefano Zampini 
100*dd072f5fSMatthew G. Knepley       if (numComps[f] < 0) continue;
101*dd072f5fSMatthew G. Knepley       PetscCall(PetscSectionGetFieldComponents(s, f, &Nc));
102*dd072f5fSMatthew 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);
103*dd072f5fSMatthew 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);
104*dd072f5fSMatthew G. Knepley       bs += numComps[f];
1059e47a1b0SMatthew G. Knepley     }
106*dd072f5fSMatthew G. Knepley   } else {
1079e47a1b0SMatthew G. Knepley     for (PetscInt f = 0; f < numFields; ++f) {
1089e47a1b0SMatthew G. Knepley       PetscInt Nc;
1099e47a1b0SMatthew G. Knepley 
1109e47a1b0SMatthew G. Knepley       PetscCall(PetscSectionGetFieldComponents(s, fields[f], &Nc));
1113a544194SStefano Zampini       bs += Nc;
1123a544194SStefano Zampini     }
113*dd072f5fSMatthew G. Knepley   }
1149e47a1b0SMatthew G. Knepley   for (PetscInt p = pStart; p < pEnd; ++p) {
115b9eca265Ssarens     PetscInt gdof, pSubSize = 0;
1164d9407bcSMatthew G. Knepley 
1179e47a1b0SMatthew G. Knepley     PetscCall(PetscSectionGetDof(gs, p, &gdof));
1184d9407bcSMatthew G. Knepley     if (gdof > 0) {
119*dd072f5fSMatthew G. Knepley       PetscInt off = 0;
1204d9407bcSMatthew G. Knepley 
121*dd072f5fSMatthew G. Knepley       for (PetscInt f = 0; f < numFields; ++f) {
122*dd072f5fSMatthew G. Knepley         PetscInt fdof, fcdof, sfdof, sfcdof = 0;
123*dd072f5fSMatthew G. Knepley 
124*dd072f5fSMatthew G. Knepley         PetscCall(PetscSectionGetFieldComponents(s, f, &Nc));
1259e47a1b0SMatthew G. Knepley         PetscCall(PetscSectionGetFieldDof(s, p, fields[f], &fdof));
1269e47a1b0SMatthew G. Knepley         PetscCall(PetscSectionGetFieldConstraintDof(s, p, fields[f], &fcdof));
127*dd072f5fSMatthew G. Knepley         if (numComps && numComps[f] >= 0) {
128*dd072f5fSMatthew G. Knepley           const PetscInt *ind;
129*dd072f5fSMatthew G. Knepley 
130*dd072f5fSMatthew G. Knepley           // Assume sets of dofs on points are of size Nc
131*dd072f5fSMatthew 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, p);
132*dd072f5fSMatthew G. Knepley           sfdof = (fdof / Nc) * numComps[f];
133*dd072f5fSMatthew G. Knepley           PetscCall(PetscSectionGetFieldConstraintIndices(s, p, fields[f], &ind));
134*dd072f5fSMatthew G. Knepley           for (PetscInt i = 0; i < (fdof / Nc); ++i) {
135*dd072f5fSMatthew G. Knepley             for (PetscInt c = 0, fcc = 0; c < Nc; ++c) {
136*dd072f5fSMatthew G. Knepley               if (c == comps[off + fcc]) {
137*dd072f5fSMatthew G. Knepley                 ++fcc;
138*dd072f5fSMatthew G. Knepley                 ++sfcdof;
139*dd072f5fSMatthew G. Knepley               }
140*dd072f5fSMatthew G. Knepley             }
141*dd072f5fSMatthew G. Knepley           }
142*dd072f5fSMatthew G. Knepley           pSubSize += sfdof - sfcdof;
143*dd072f5fSMatthew G. Knepley           off += numComps[f];
144*dd072f5fSMatthew G. Knepley         } else {
145b9eca265Ssarens           pSubSize += fdof - fcdof;
146b9eca265Ssarens         }
147*dd072f5fSMatthew G. Knepley       }
148b9eca265Ssarens       subSize += pSubSize;
1493a544194SStefano Zampini       if (pSubSize && bs != pSubSize) {
1509e47a1b0SMatthew G. Knepley         // Layout does not admit a pointwise block size
151b9eca265Ssarens         bs = 1;
1524d9407bcSMatthew G. Knepley       }
1534d9407bcSMatthew G. Knepley     }
1544d9407bcSMatthew G. Knepley   }
1559e47a1b0SMatthew G. Knepley   // Must have same blocksize on all procs (some might have no points)
1569371c9d4SSatish Balay   bsLocal[0] = bs < 0 ? PETSC_MAX_INT : bs;
1579371c9d4SSatish Balay   bsLocal[1] = bs;
1589e47a1b0SMatthew G. Knepley   PetscCall(PetscGlobalMinMaxInt(PetscObjectComm((PetscObject)gs), bsLocal, bsMinMax));
1599371c9d4SSatish Balay   if (bsMinMax[0] != bsMinMax[1]) {
1609371c9d4SSatish Balay     bs = 1;
1619371c9d4SSatish Balay   } else {
1629371c9d4SSatish Balay     bs = bsMinMax[0];
1639371c9d4SSatish Balay   }
1649566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(subSize, &subIndices));
1659e47a1b0SMatthew G. Knepley   for (PetscInt p = pStart; p < pEnd; ++p) {
1664d9407bcSMatthew G. Knepley     PetscInt gdof, goff;
1674d9407bcSMatthew G. Knepley 
1689e47a1b0SMatthew G. Knepley     PetscCall(PetscSectionGetDof(gs, p, &gdof));
1694d9407bcSMatthew G. Knepley     if (gdof > 0) {
170*dd072f5fSMatthew G. Knepley       PetscInt off = 0;
171*dd072f5fSMatthew G. Knepley 
1729e47a1b0SMatthew G. Knepley       PetscCall(PetscSectionGetOffset(gs, p, &goff));
1739e47a1b0SMatthew G. Knepley       for (PetscInt f = 0; f < numFields; ++f) {
1749e47a1b0SMatthew G. Knepley         PetscInt fdof, fcdof, poff = 0;
1754d9407bcSMatthew G. Knepley 
1764d9407bcSMatthew G. Knepley         /* Can get rid of this loop by storing field information in the global section */
1779e47a1b0SMatthew G. Knepley         for (PetscInt f2 = 0; f2 < fields[f]; ++f2) {
1789e47a1b0SMatthew G. Knepley           PetscCall(PetscSectionGetFieldDof(s, p, f2, &fdof));
1799e47a1b0SMatthew G. Knepley           PetscCall(PetscSectionGetFieldConstraintDof(s, p, f2, &fcdof));
1804d9407bcSMatthew G. Knepley           poff += fdof - fcdof;
1814d9407bcSMatthew G. Knepley         }
1829e47a1b0SMatthew G. Knepley         PetscCall(PetscSectionGetFieldDof(s, p, fields[f], &fdof));
1839e47a1b0SMatthew G. Knepley         PetscCall(PetscSectionGetFieldConstraintDof(s, p, fields[f], &fcdof));
184*dd072f5fSMatthew G. Knepley 
185*dd072f5fSMatthew G. Knepley         if (numComps && numComps[f] >= 0) {
186*dd072f5fSMatthew G. Knepley           const PetscInt *ind;
187*dd072f5fSMatthew G. Knepley 
188*dd072f5fSMatthew G. Knepley           // Assume sets of dofs on points are of size Nc
189*dd072f5fSMatthew G. Knepley           PetscCall(PetscSectionGetFieldConstraintIndices(s, p, fields[f], &ind));
190*dd072f5fSMatthew G. Knepley           for (PetscInt i = 0, fcoff = 0, pfoff = 0; i < (fdof / Nc); ++i) {
191*dd072f5fSMatthew G. Knepley             for (PetscInt c = 0, fcc = 0; c < Nc; ++c) {
192*dd072f5fSMatthew G. Knepley               const PetscInt k = i * Nc + c;
193*dd072f5fSMatthew G. Knepley 
194*dd072f5fSMatthew G. Knepley               if (ind[fcoff] == k) {
195*dd072f5fSMatthew G. Knepley                 ++fcoff;
196*dd072f5fSMatthew G. Knepley                 continue;
197*dd072f5fSMatthew G. Knepley               }
198*dd072f5fSMatthew G. Knepley               if (c == comps[off + fcc]) {
199*dd072f5fSMatthew G. Knepley                 ++fcc;
200*dd072f5fSMatthew G. Knepley                 subIndices[subOff++] = goff + poff + pfoff;
201*dd072f5fSMatthew G. Knepley               }
202*dd072f5fSMatthew G. Knepley               ++pfoff;
203*dd072f5fSMatthew G. Knepley             }
204*dd072f5fSMatthew G. Knepley           }
205*dd072f5fSMatthew G. Knepley           off += numComps[f];
206*dd072f5fSMatthew G. Knepley         } else {
2079e47a1b0SMatthew G. Knepley           for (PetscInt fc = 0; fc < fdof - fcdof; ++fc, ++subOff) subIndices[subOff] = goff + poff + fc;
2084d9407bcSMatthew G. Knepley         }
2094d9407bcSMatthew G. Knepley       }
2104d9407bcSMatthew G. Knepley     }
211*dd072f5fSMatthew G. Knepley   }
2129e47a1b0SMatthew 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);
2139e47a1b0SMatthew G. Knepley   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)gs), subSize, subIndices, PETSC_OWN_POINTER, is));
214627349f5SMatthew G. Knepley   if (bs > 1) {
2159e47a1b0SMatthew G. Knepley     // We need to check that the block size does not come from non-contiguous fields
2169e47a1b0SMatthew G. Knepley     PetscInt set = 1, rset = 1;
2179e47a1b0SMatthew G. Knepley     for (PetscInt i = 0; i < subSize; i += bs) {
2189e47a1b0SMatthew G. Knepley       for (PetscInt j = 0; j < bs; ++j) {
2199371c9d4SSatish Balay         if (subIndices[i + j] != subIndices[i] + j) {
2209371c9d4SSatish Balay           set = 0;
2219371c9d4SSatish Balay           break;
2229371c9d4SSatish Balay         }
223627349f5SMatthew G. Knepley       }
224627349f5SMatthew G. Knepley     }
2259e47a1b0SMatthew G. Knepley     PetscCall(MPIU_Allreduce(&set, &rset, 1, MPIU_INT, MPI_PROD, PetscObjectComm((PetscObject)gs)));
22688206212SAlexis Marboeuf     if (rset) PetscCall(ISSetBlockSize(*is, bs));
227627349f5SMatthew G. Knepley   }
2289e47a1b0SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
2294d9407bcSMatthew G. Knepley }
2309e47a1b0SMatthew G. Knepley 
231*dd072f5fSMatthew 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)
2329e47a1b0SMatthew G. Knepley {
2334d9407bcSMatthew G. Knepley   PetscSection subsection;
2344d9407bcSMatthew G. Knepley   PetscBool    haveNull = PETSC_FALSE;
2359e47a1b0SMatthew G. Knepley   PetscInt     nf = 0, of = 0;
2364d9407bcSMatthew G. Knepley 
2379e47a1b0SMatthew G. Knepley   PetscFunctionBegin;
238*dd072f5fSMatthew G. Knepley   if (numComps) {
239*dd072f5fSMatthew G. Knepley     const PetscInt field = fields[0];
240*dd072f5fSMatthew G. Knepley 
241*dd072f5fSMatthew G. Knepley     PetscCheck(numFields == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "We only support a single field for component selection right now");
242*dd072f5fSMatthew G. Knepley     PetscCall(PetscSectionCreateComponentSubsection(section, numComps[field], comps, &subsection));
2439566063dSJacob Faibussowitsch     PetscCall(DMSetLocalSection(*subdm, subsection));
2449566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&subsection));
245*dd072f5fSMatthew G. Knepley     (*subdm)->nullspaceConstructors[field] = dm->nullspaceConstructors[field];
246e5e52638SMatthew G. Knepley     if (dm->probs) {
2479e47a1b0SMatthew G. Knepley       PetscFV  fv, fvNew;
248*dd072f5fSMatthew G. Knepley       PetscInt fnum[1] = {field};
2490f21e855SMatthew G. Knepley 
2509e47a1b0SMatthew G. Knepley       PetscCall(DMSetNumFields(*subdm, 1));
251*dd072f5fSMatthew G. Knepley       PetscCall(DMGetField(dm, field, NULL, (PetscObject *)&fv));
2529e47a1b0SMatthew G. Knepley       PetscCall(PetscFVClone(fv, &fvNew));
253*dd072f5fSMatthew G. Knepley       PetscCall(PetscFVSetNumComponents(fvNew, numComps[0]));
2549e47a1b0SMatthew G. Knepley       PetscCall(DMSetField(*subdm, 0, NULL, (PetscObject)fvNew));
2559e47a1b0SMatthew G. Knepley       PetscCall(PetscFVDestroy(&fvNew));
2569566063dSJacob Faibussowitsch       PetscCall(DMCreateDS(*subdm));
257*dd072f5fSMatthew G. Knepley       if (numComps[0] == 1 && is) {
2580f21e855SMatthew G. Knepley         PetscObject disc, space, pmat;
2594d9407bcSMatthew G. Knepley 
260*dd072f5fSMatthew G. Knepley         PetscCall(DMGetField(*subdm, field, NULL, &disc));
2619566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "nullspace", &space));
2629566063dSJacob Faibussowitsch         if (space) PetscCall(PetscObjectCompose((PetscObject)*is, "nullspace", space));
2639566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "nearnullspace", &space));
2649566063dSJacob Faibussowitsch         if (space) PetscCall(PetscObjectCompose((PetscObject)*is, "nearnullspace", space));
2659566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "pmat", &pmat));
2669566063dSJacob Faibussowitsch         if (pmat) PetscCall(PetscObjectCompose((PetscObject)*is, "pmat", pmat));
2674d9407bcSMatthew G. Knepley       }
268*dd072f5fSMatthew G. Knepley       PetscCall(PetscDSCopyConstants(dm->probs[field].ds, (*subdm)->probs[0].ds));
269*dd072f5fSMatthew G. Knepley       PetscCall(PetscDSCopyBoundary(dm->probs[field].ds, 1, fnum, (*subdm)->probs[0].ds));
270*dd072f5fSMatthew G. Knepley       PetscCall(PetscDSSelectEquations(dm->probs[field].ds, 1, fnum, (*subdm)->probs[0].ds));
2719e47a1b0SMatthew G. Knepley     }
2729e47a1b0SMatthew G. Knepley     if ((*subdm)->nullspaceConstructors[0] && is) {
2739e47a1b0SMatthew G. Knepley       MatNullSpace nullSpace;
2749e47a1b0SMatthew G. Knepley 
2759e47a1b0SMatthew G. Knepley       PetscCall((*(*subdm)->nullspaceConstructors[0])(*subdm, 0, 0, &nullSpace));
2769e47a1b0SMatthew G. Knepley       PetscCall(PetscObjectCompose((PetscObject)*is, "nullspace", (PetscObject)nullSpace));
2779e47a1b0SMatthew G. Knepley       PetscCall(MatNullSpaceDestroy(&nullSpace));
2789e47a1b0SMatthew G. Knepley     }
2799e47a1b0SMatthew G. Knepley     if (dm->coarseMesh) PetscCall(DMCreateSubDM(dm->coarseMesh, numFields, fields, NULL, &(*subdm)->coarseMesh));
2809e47a1b0SMatthew G. Knepley     PetscFunctionReturn(PETSC_SUCCESS);
2819e47a1b0SMatthew G. Knepley   }
2829e47a1b0SMatthew G. Knepley   PetscCall(PetscSectionCreateSubsection(section, numFields, fields, &subsection));
2839e47a1b0SMatthew G. Knepley   PetscCall(DMSetLocalSection(*subdm, subsection));
2849e47a1b0SMatthew G. Knepley   PetscCall(PetscSectionDestroy(&subsection));
2859e47a1b0SMatthew G. Knepley   for (PetscInt f = 0; f < numFields; ++f) {
2869e47a1b0SMatthew G. Knepley     (*subdm)->nullspaceConstructors[f] = dm->nullspaceConstructors[fields[f]];
2879e47a1b0SMatthew G. Knepley     if ((*subdm)->nullspaceConstructors[f]) {
2889e47a1b0SMatthew G. Knepley       haveNull = PETSC_TRUE;
2899e47a1b0SMatthew G. Knepley       nf       = f;
2909e47a1b0SMatthew G. Knepley       of       = fields[f];
2919e47a1b0SMatthew G. Knepley     }
2929e47a1b0SMatthew G. Knepley   }
2939e47a1b0SMatthew G. Knepley   if (dm->probs) {
2949e47a1b0SMatthew G. Knepley     PetscCall(DMSetNumFields(*subdm, numFields));
2959e47a1b0SMatthew G. Knepley     for (PetscInt f = 0; f < numFields; ++f) {
2969e47a1b0SMatthew G. Knepley       PetscObject disc;
2979e47a1b0SMatthew G. Knepley 
2989e47a1b0SMatthew G. Knepley       PetscCall(DMGetField(dm, fields[f], NULL, &disc));
2999e47a1b0SMatthew G. Knepley       PetscCall(DMSetField(*subdm, f, NULL, disc));
3009e47a1b0SMatthew G. Knepley     }
3019e47a1b0SMatthew G. Knepley     // TODO: if only FV, then cut down the components
3029e47a1b0SMatthew G. Knepley     PetscCall(DMCreateDS(*subdm));
3039e47a1b0SMatthew G. Knepley     if (numFields == 1 && is) {
3049e47a1b0SMatthew G. Knepley       PetscObject disc, space, pmat;
3059e47a1b0SMatthew G. Knepley 
3069e47a1b0SMatthew G. Knepley       PetscCall(DMGetField(*subdm, 0, NULL, &disc));
3079e47a1b0SMatthew G. Knepley       PetscCall(PetscObjectQuery(disc, "nullspace", &space));
3089e47a1b0SMatthew G. Knepley       if (space) PetscCall(PetscObjectCompose((PetscObject)*is, "nullspace", space));
3099e47a1b0SMatthew G. Knepley       PetscCall(PetscObjectQuery(disc, "nearnullspace", &space));
3109e47a1b0SMatthew G. Knepley       if (space) PetscCall(PetscObjectCompose((PetscObject)*is, "nearnullspace", space));
3119e47a1b0SMatthew G. Knepley       PetscCall(PetscObjectQuery(disc, "pmat", &pmat));
3129e47a1b0SMatthew G. Knepley       if (pmat) PetscCall(PetscObjectCompose((PetscObject)*is, "pmat", pmat));
3139e47a1b0SMatthew G. Knepley     }
3149e47a1b0SMatthew G. Knepley     // Check if DSes record their DM fields
315e21d936aSMatthew G. Knepley     if (dm->probs[0].fields) {
3169310035eSMatthew G. Knepley       PetscInt d, e;
3179310035eSMatthew G. Knepley 
3189310035eSMatthew G. Knepley       for (d = 0, e = 0; d < dm->Nds && e < (*subdm)->Nds; ++d) {
3199310035eSMatthew G. Knepley         const PetscInt  Nf = dm->probs[d].ds->Nf;
3209310035eSMatthew G. Knepley         const PetscInt *fld;
3219310035eSMatthew G. Knepley         PetscInt        f, g;
3229310035eSMatthew G. Knepley 
3239566063dSJacob Faibussowitsch         PetscCall(ISGetIndices(dm->probs[d].fields, &fld));
3249310035eSMatthew G. Knepley         for (f = 0; f < Nf; ++f) {
3259371c9d4SSatish Balay           for (g = 0; g < numFields; ++g)
3269371c9d4SSatish Balay             if (fld[f] == fields[g]) break;
3279310035eSMatthew G. Knepley           if (g < numFields) break;
3289310035eSMatthew G. Knepley         }
3299566063dSJacob Faibussowitsch         PetscCall(ISRestoreIndices(dm->probs[d].fields, &fld));
3309310035eSMatthew G. Knepley         if (f == Nf) continue;
3319566063dSJacob Faibussowitsch         PetscCall(PetscDSCopyConstants(dm->probs[d].ds, (*subdm)->probs[e].ds));
3329566063dSJacob Faibussowitsch         PetscCall(PetscDSCopyBoundary(dm->probs[d].ds, numFields, fields, (*subdm)->probs[e].ds));
3339e47a1b0SMatthew G. Knepley         // Translate DM fields to DS fields
3349310035eSMatthew G. Knepley         {
335e21d936aSMatthew G. Knepley           IS              infields, dsfields;
33674a61aaaSMatthew G. Knepley           const PetscInt *fld, *ofld;
33774a61aaaSMatthew G. Knepley           PetscInt       *fidx;
3389e47a1b0SMatthew G. Knepley           PetscInt        onf, nf;
339e21d936aSMatthew G. Knepley 
3409566063dSJacob Faibussowitsch           PetscCall(ISCreateGeneral(PETSC_COMM_SELF, numFields, fields, PETSC_USE_POINTER, &infields));
3419566063dSJacob Faibussowitsch           PetscCall(ISIntersect(infields, dm->probs[d].fields, &dsfields));
3429566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&infields));
3439566063dSJacob Faibussowitsch           PetscCall(ISGetLocalSize(dsfields, &nf));
3447a8be351SBarry Smith           PetscCheck(nf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "DS cannot be supported on 0 fields");
3459566063dSJacob Faibussowitsch           PetscCall(ISGetIndices(dsfields, &fld));
3469566063dSJacob Faibussowitsch           PetscCall(ISGetLocalSize(dm->probs[d].fields, &onf));
3479566063dSJacob Faibussowitsch           PetscCall(ISGetIndices(dm->probs[d].fields, &ofld));
3489566063dSJacob Faibussowitsch           PetscCall(PetscMalloc1(nf, &fidx));
3499e47a1b0SMatthew G. Knepley           for (PetscInt f = 0, g = 0; f < onf && g < nf; ++f) {
35074a61aaaSMatthew G. Knepley             if (ofld[f] == fld[g]) fidx[g++] = f;
35174a61aaaSMatthew G. Knepley           }
3529566063dSJacob Faibussowitsch           PetscCall(ISRestoreIndices(dm->probs[d].fields, &ofld));
3539566063dSJacob Faibussowitsch           PetscCall(ISRestoreIndices(dsfields, &fld));
3549566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&dsfields));
3559566063dSJacob Faibussowitsch           PetscCall(PetscDSSelectDiscretizations(dm->probs[0].ds, nf, fidx, (*subdm)->probs[0].ds));
3569566063dSJacob Faibussowitsch           PetscCall(PetscDSSelectEquations(dm->probs[0].ds, nf, fidx, (*subdm)->probs[0].ds));
3579566063dSJacob Faibussowitsch           PetscCall(PetscFree(fidx));
3589310035eSMatthew G. Knepley         }
3599310035eSMatthew G. Knepley         ++e;
3609310035eSMatthew G. Knepley       }
361e21d936aSMatthew G. Knepley     } else {
3629566063dSJacob Faibussowitsch       PetscCall(PetscDSCopyConstants(dm->probs[0].ds, (*subdm)->probs[0].ds));
3639566063dSJacob Faibussowitsch       PetscCall(PetscDSCopyBoundary(dm->probs[0].ds, PETSC_DETERMINE, NULL, (*subdm)->probs[0].ds));
3649566063dSJacob Faibussowitsch       PetscCall(PetscDSSelectDiscretizations(dm->probs[0].ds, numFields, fields, (*subdm)->probs[0].ds));
3659566063dSJacob Faibussowitsch       PetscCall(PetscDSSelectEquations(dm->probs[0].ds, numFields, fields, (*subdm)->probs[0].ds));
366d5af271aSMatthew G. Knepley     }
367e21d936aSMatthew G. Knepley   }
3688cda7954SMatthew G. Knepley   if (haveNull && is) {
3698cda7954SMatthew G. Knepley     MatNullSpace nullSpace;
3708cda7954SMatthew G. Knepley 
3719566063dSJacob Faibussowitsch     PetscCall((*(*subdm)->nullspaceConstructors[nf])(*subdm, of, nf, &nullSpace));
3729566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)*is, "nullspace", (PetscObject)nullSpace));
3739566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceDestroy(&nullSpace));
3748cda7954SMatthew G. Knepley   }
37548a46eb9SPierre Jolivet   if (dm->coarseMesh) PetscCall(DMCreateSubDM(dm->coarseMesh, numFields, fields, NULL, &(*subdm)->coarseMesh));
3769e47a1b0SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
3774d9407bcSMatthew G. Knepley }
3789e47a1b0SMatthew G. Knepley 
3799e47a1b0SMatthew G. Knepley /*@C
3809e47a1b0SMatthew 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`.
3819e47a1b0SMatthew G. Knepley 
3829e47a1b0SMatthew G. Knepley   Not Collective
3839e47a1b0SMatthew G. Knepley 
3849e47a1b0SMatthew G. Knepley   Input Parameters:
3859e47a1b0SMatthew G. Knepley + dm        - The `DM` object
3869e47a1b0SMatthew G. Knepley . numFields - The number of fields to incorporate into `subdm`
387*dd072f5fSMatthew G. Knepley . fields    - The field numbers of the selected fields
388*dd072f5fSMatthew G. Knepley . numComps  - The number of components from each field to incorporate into `subdm`, or PETSC_DECIDE for all components
389*dd072f5fSMatthew G. Knepley - comps     - The component numbers of the selected fields (omitted for PTESC_DECIDE fields)
3909e47a1b0SMatthew G. Knepley 
3919e47a1b0SMatthew G. Knepley   Output Parameters:
3929e47a1b0SMatthew G. Knepley + is    - The global indices for the subproblem or `NULL`
3939e47a1b0SMatthew G. Knepley - subdm - The `DM` for the subproblem, which must already have be cloned from `dm` or `NULL`
3949e47a1b0SMatthew G. Knepley 
3959e47a1b0SMatthew G. Knepley   Level: intermediate
3969e47a1b0SMatthew G. Knepley 
3979e47a1b0SMatthew G. Knepley   Notes:
3989e47a1b0SMatthew G. Knepley   If `is` and `subdm` are both `NULL` this does nothing
3999e47a1b0SMatthew G. Knepley 
4009e47a1b0SMatthew G. Knepley .seealso: `DMCreateSubDM()`, `DMGetLocalSection()`, `DMPlexSetMigrationSF()`, `DMView()`
4019e47a1b0SMatthew G. Knepley @*/
402*dd072f5fSMatthew G. Knepley PetscErrorCode DMCreateSectionSubDM(DM dm, PetscInt numFields, const PetscInt fields[], const PetscInt numComps[], const PetscInt comps[], IS *is, DM *subdm)
4039e47a1b0SMatthew G. Knepley {
4049e47a1b0SMatthew G. Knepley   PetscSection section, sectionGlobal;
4059e47a1b0SMatthew G. Knepley   PetscInt     Nf;
4069e47a1b0SMatthew G. Knepley 
4079e47a1b0SMatthew G. Knepley   PetscFunctionBegin;
4089e47a1b0SMatthew G. Knepley   if (!numFields) PetscFunctionReturn(PETSC_SUCCESS);
4099e47a1b0SMatthew G. Knepley   PetscCall(DMGetLocalSection(dm, &section));
4109e47a1b0SMatthew G. Knepley   PetscCall(DMGetGlobalSection(dm, &sectionGlobal));
4119e47a1b0SMatthew G. Knepley   PetscCheck(section, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields");
4129e47a1b0SMatthew G. Knepley   PetscCheck(sectionGlobal, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields");
4139e47a1b0SMatthew G. Knepley   PetscCall(PetscSectionGetNumFields(section, &Nf));
4149e47a1b0SMatthew 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);
4159e47a1b0SMatthew G. Knepley 
416*dd072f5fSMatthew G. Knepley   if (is) PetscCall(PetscSectionSelectFields_Private(section, sectionGlobal, numFields, fields, numComps, comps, is));
417*dd072f5fSMatthew G. Knepley   if (subdm) PetscCall(DMSelectFields_Private(dm, section, numFields, fields, numComps, comps, is, subdm));
4183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4194d9407bcSMatthew G. Knepley }
4202adcc780SMatthew G. Knepley 
421792b654fSMatthew G. Knepley /*@C
4220b3275a6SBarry 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`
423792b654fSMatthew G. Knepley 
42420f4b53cSBarry Smith   Not Collective
425792b654fSMatthew G. Knepley 
426d8d19677SJose E. Roman   Input Parameters:
4270b3275a6SBarry Smith + dms - The `DM` objects, the must all have the same topology; for example obtained with `DMClone()`
4280b3275a6SBarry Smith - len - The number of `DM` in `dms`
429792b654fSMatthew G. Knepley 
430792b654fSMatthew G. Knepley   Output Parameters:
43120f4b53cSBarry Smith + is      - The global indices for the subproblem, or `NULL`
4320b3275a6SBarry Smith - superdm - The `DM` for the superproblem, which must already have be cloned and contain the same topology as the `dms`
433792b654fSMatthew G. Knepley 
434792b654fSMatthew G. Knepley   Level: intermediate
435792b654fSMatthew G. Knepley 
436a4e35b19SJacob Faibussowitsch .seealso: `DMCreateSuperDM()`, `DMGetLocalSection()`, `DMPlexSetMigrationSF()`, `DMView()`
437792b654fSMatthew G. Knepley @*/
438d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateSectionSuperDM(DM dms[], PetscInt len, IS **is, DM *superdm)
439d71ae5a4SJacob Faibussowitsch {
4409796d432SMatthew G. Knepley   MPI_Comm     comm;
4419796d432SMatthew G. Knepley   PetscSection supersection, *sections, *sectionGlobals;
4428cda7954SMatthew G. Knepley   PetscInt    *Nfs, Nf = 0, f, supf, oldf = -1, nullf = -1, i;
4439796d432SMatthew G. Knepley   PetscBool    haveNull = PETSC_FALSE;
4442adcc780SMatthew G. Knepley 
4452adcc780SMatthew G. Knepley   PetscFunctionBegin;
4469566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dms[0], &comm));
4479796d432SMatthew G. Knepley   /* Pull out local and global sections */
4489566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(len, &Nfs, len, &sections, len, &sectionGlobals));
4492adcc780SMatthew G. Knepley   for (i = 0; i < len; ++i) {
4509566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dms[i], &sections[i]));
4519566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalSection(dms[i], &sectionGlobals[i]));
4527a8be351SBarry Smith     PetscCheck(sections[i], comm, PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields");
4537a8be351SBarry Smith     PetscCheck(sectionGlobals[i], comm, PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields");
4549566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetNumFields(sections[i], &Nfs[i]));
4552adcc780SMatthew G. Knepley     Nf += Nfs[i];
4562adcc780SMatthew G. Knepley   }
4579796d432SMatthew G. Knepley   /* Create the supersection */
4589566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreateSupersection(sections, len, &supersection));
4599566063dSJacob Faibussowitsch   PetscCall(DMSetLocalSection(*superdm, supersection));
4609796d432SMatthew G. Knepley   /* Create ISes */
4612adcc780SMatthew G. Knepley   if (is) {
4629796d432SMatthew G. Knepley     PetscSection supersectionGlobal;
4639796d432SMatthew G. Knepley     PetscInt     bs = -1, startf = 0;
4642adcc780SMatthew G. Knepley 
4659566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(len, is));
4669566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalSection(*superdm, &supersectionGlobal));
4679796d432SMatthew G. Knepley     for (i = 0; i < len; startf += Nfs[i], ++i) {
4689796d432SMatthew G. Knepley       PetscInt *subIndices;
469ec4c761aSMatthew G. Knepley       PetscInt  subSize, subOff, pStart, pEnd, p, start, end, dummy;
4702adcc780SMatthew G. Knepley 
4719566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetChart(sectionGlobals[i], &pStart, &pEnd));
4729566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetConstrainedStorageSize(sectionGlobals[i], &subSize));
4739566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(subSize, &subIndices));
4749796d432SMatthew G. Knepley       for (p = pStart, subOff = 0; p < pEnd; ++p) {
4759796d432SMatthew G. Knepley         PetscInt gdof, gcdof, gtdof, d;
4762adcc780SMatthew G. Knepley 
4779566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(sectionGlobals[i], p, &gdof));
4789566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetConstraintDof(sections[i], p, &gcdof));
4792adcc780SMatthew G. Knepley         gtdof = gdof - gcdof;
4802adcc780SMatthew G. Knepley         if (gdof > 0 && gtdof) {
4819371c9d4SSatish Balay           if (bs < 0) {
4829371c9d4SSatish Balay             bs = gtdof;
4839371c9d4SSatish Balay           } else if (bs != gtdof) {
4849371c9d4SSatish Balay             bs = 1;
4859371c9d4SSatish Balay           }
4869566063dSJacob Faibussowitsch           PetscCall(DMGetGlobalFieldOffset_Private(*superdm, p, startf, &start, &dummy));
4879566063dSJacob Faibussowitsch           PetscCall(DMGetGlobalFieldOffset_Private(*superdm, p, startf + Nfs[i] - 1, &dummy, &end));
48863a3b9bcSJacob 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);
4899796d432SMatthew G. Knepley           for (d = start; d < end; ++d, ++subOff) subIndices[subOff] = d;
4902adcc780SMatthew G. Knepley         }
4912adcc780SMatthew G. Knepley       }
4929566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(comm, subSize, subIndices, PETSC_OWN_POINTER, &(*is)[i]));
4932adcc780SMatthew G. Knepley       /* Must have same blocksize on all procs (some might have no points) */
4949796d432SMatthew G. Knepley       {
4959796d432SMatthew G. Knepley         PetscInt bs = -1, bsLocal[2], bsMinMax[2];
4969796d432SMatthew G. Knepley 
4979371c9d4SSatish Balay         bsLocal[0] = bs < 0 ? PETSC_MAX_INT : bs;
4989371c9d4SSatish Balay         bsLocal[1] = bs;
4999566063dSJacob Faibussowitsch         PetscCall(PetscGlobalMinMaxInt(comm, bsLocal, bsMinMax));
5009371c9d4SSatish Balay         if (bsMinMax[0] != bsMinMax[1]) {
5019371c9d4SSatish Balay           bs = 1;
5029371c9d4SSatish Balay         } else {
5039371c9d4SSatish Balay           bs = bsMinMax[0];
5049371c9d4SSatish Balay         }
5059566063dSJacob Faibussowitsch         PetscCall(ISSetBlockSize((*is)[i], bs));
5062adcc780SMatthew G. Knepley       }
5072adcc780SMatthew G. Knepley     }
5082adcc780SMatthew G. Knepley   }
5099796d432SMatthew G. Knepley   /* Preserve discretizations */
510e5e52638SMatthew G. Knepley   if (len && dms[0]->probs) {
5119566063dSJacob Faibussowitsch     PetscCall(DMSetNumFields(*superdm, Nf));
5129796d432SMatthew G. Knepley     for (i = 0, supf = 0; i < len; ++i) {
5139796d432SMatthew G. Knepley       for (f = 0; f < Nfs[i]; ++f, ++supf) {
5142adcc780SMatthew G. Knepley         PetscObject disc;
5152adcc780SMatthew G. Knepley 
5169566063dSJacob Faibussowitsch         PetscCall(DMGetField(dms[i], f, NULL, &disc));
5179566063dSJacob Faibussowitsch         PetscCall(DMSetField(*superdm, supf, NULL, disc));
5182adcc780SMatthew G. Knepley       }
5192adcc780SMatthew G. Knepley     }
5209566063dSJacob Faibussowitsch     PetscCall(DMCreateDS(*superdm));
5212adcc780SMatthew G. Knepley   }
5229796d432SMatthew G. Knepley   /* Preserve nullspaces */
5239796d432SMatthew G. Knepley   for (i = 0, supf = 0; i < len; ++i) {
5249796d432SMatthew G. Knepley     for (f = 0; f < Nfs[i]; ++f, ++supf) {
5259796d432SMatthew G. Knepley       (*superdm)->nullspaceConstructors[supf] = dms[i]->nullspaceConstructors[f];
5269796d432SMatthew G. Knepley       if ((*superdm)->nullspaceConstructors[supf]) {
5279796d432SMatthew G. Knepley         haveNull = PETSC_TRUE;
5289796d432SMatthew G. Knepley         nullf    = supf;
5298cda7954SMatthew G. Knepley         oldf     = f;
5302adcc780SMatthew G. Knepley       }
5319796d432SMatthew G. Knepley     }
5329796d432SMatthew G. Knepley   }
5339796d432SMatthew G. Knepley   /* Attach nullspace to IS */
5349796d432SMatthew G. Knepley   if (haveNull && is) {
5359796d432SMatthew G. Knepley     MatNullSpace nullSpace;
5369796d432SMatthew G. Knepley 
5379566063dSJacob Faibussowitsch     PetscCall((*(*superdm)->nullspaceConstructors[nullf])(*superdm, oldf, nullf, &nullSpace));
5389566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)(*is)[nullf], "nullspace", (PetscObject)nullSpace));
5399566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceDestroy(&nullSpace));
5409796d432SMatthew G. Knepley   }
5419566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&supersection));
5429566063dSJacob Faibussowitsch   PetscCall(PetscFree3(Nfs, sections, sectionGlobals));
5433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5442adcc780SMatthew G. Knepley }
545