xref: /petsc/src/dm/interface/dmi.c (revision bb4b53ef092968f72b740b90dbab8a2b6700db0d)
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 
400b697e01SMatthew G. Knepley   // You cannot negate PETSC_MIN_INT
410b697e01SMatthew G. Knepley   in[0] = blockSize < 0 ? -PETSC_MAX_INT : -blockSize;
425d1c0279SHong Zhang   in[1] = blockSize;
431c2dc1cbSBarry Smith   PetscCall(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));
839566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(*vec, 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 {
914d9407bcSMatthew G. Knepley   PetscInt *subIndices;
92dd072f5fSMatthew G. Knepley   PetscInt  bs = 0, bsLocal[2], bsMinMax[2];
93dd072f5fSMatthew G. Knepley   PetscInt  pStart, pEnd, Nc, subSize = 0, subOff = 0;
944d9407bcSMatthew G. Knepley 
954d9407bcSMatthew G. Knepley   PetscFunctionBegin;
969e47a1b0SMatthew G. Knepley   PetscCall(PetscSectionGetChart(gs, &pStart, &pEnd));
97dd072f5fSMatthew G. Knepley   if (numComps) {
98dd072f5fSMatthew G. Knepley     for (PetscInt f = 0, off = 0; f < numFields; ++f) {
993a544194SStefano Zampini       PetscInt Nc;
1003a544194SStefano Zampini 
101dd072f5fSMatthew G. Knepley       if (numComps[f] < 0) continue;
102dd072f5fSMatthew G. Knepley       PetscCall(PetscSectionGetFieldComponents(s, f, &Nc));
103dd072f5fSMatthew 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);
104dd072f5fSMatthew 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);
105dd072f5fSMatthew G. Knepley       bs += numComps[f];
1069e47a1b0SMatthew G. Knepley     }
107dd072f5fSMatthew G. Knepley   } else {
1089e47a1b0SMatthew G. Knepley     for (PetscInt f = 0; f < numFields; ++f) {
1099e47a1b0SMatthew G. Knepley       PetscInt Nc;
1109e47a1b0SMatthew G. Knepley 
1119e47a1b0SMatthew G. Knepley       PetscCall(PetscSectionGetFieldComponents(s, fields[f], &Nc));
1123a544194SStefano Zampini       bs += Nc;
1133a544194SStefano Zampini     }
114dd072f5fSMatthew G. Knepley   }
1159e47a1b0SMatthew G. Knepley   for (PetscInt p = pStart; p < pEnd; ++p) {
116b9eca265Ssarens     PetscInt gdof, pSubSize = 0;
1174d9407bcSMatthew G. Knepley 
1189e47a1b0SMatthew G. Knepley     PetscCall(PetscSectionGetDof(gs, p, &gdof));
1194d9407bcSMatthew G. Knepley     if (gdof > 0) {
120dd072f5fSMatthew G. Knepley       PetscInt off = 0;
1214d9407bcSMatthew G. Knepley 
122dd072f5fSMatthew G. Knepley       for (PetscInt f = 0; f < numFields; ++f) {
123dd072f5fSMatthew G. Knepley         PetscInt fdof, fcdof, sfdof, sfcdof = 0;
124dd072f5fSMatthew G. Knepley 
125dd072f5fSMatthew G. Knepley         PetscCall(PetscSectionGetFieldComponents(s, f, &Nc));
1269e47a1b0SMatthew G. Knepley         PetscCall(PetscSectionGetFieldDof(s, p, fields[f], &fdof));
1279e47a1b0SMatthew G. Knepley         PetscCall(PetscSectionGetFieldConstraintDof(s, p, fields[f], &fcdof));
128dd072f5fSMatthew G. Knepley         if (numComps && numComps[f] >= 0) {
129dd072f5fSMatthew G. Knepley           const PetscInt *ind;
130dd072f5fSMatthew G. Knepley 
131dd072f5fSMatthew G. Knepley           // Assume sets of dofs on points are of size Nc
132dd072f5fSMatthew 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);
133dd072f5fSMatthew G. Knepley           sfdof = (fdof / Nc) * numComps[f];
134dd072f5fSMatthew G. Knepley           PetscCall(PetscSectionGetFieldConstraintIndices(s, p, fields[f], &ind));
135dd072f5fSMatthew G. Knepley           for (PetscInt i = 0; i < (fdof / Nc); ++i) {
136dd072f5fSMatthew G. Knepley             for (PetscInt c = 0, fcc = 0; c < Nc; ++c) {
137dd072f5fSMatthew G. Knepley               if (c == comps[off + fcc]) {
138dd072f5fSMatthew G. Knepley                 ++fcc;
139dd072f5fSMatthew G. Knepley                 ++sfcdof;
140dd072f5fSMatthew G. Knepley               }
141dd072f5fSMatthew G. Knepley             }
142dd072f5fSMatthew G. Knepley           }
143dd072f5fSMatthew G. Knepley           pSubSize += sfdof - sfcdof;
144dd072f5fSMatthew G. Knepley           off += numComps[f];
145dd072f5fSMatthew G. Knepley         } else {
146b9eca265Ssarens           pSubSize += fdof - fcdof;
147b9eca265Ssarens         }
148dd072f5fSMatthew G. Knepley       }
149b9eca265Ssarens       subSize += pSubSize;
1503a544194SStefano Zampini       if (pSubSize && bs != pSubSize) {
1519e47a1b0SMatthew G. Knepley         // Layout does not admit a pointwise block size
152b9eca265Ssarens         bs = 1;
1534d9407bcSMatthew G. Knepley       }
1544d9407bcSMatthew G. Knepley     }
1554d9407bcSMatthew G. Knepley   }
1569e47a1b0SMatthew G. Knepley   // Must have same blocksize on all procs (some might have no points)
1579371c9d4SSatish Balay   bsLocal[0] = bs < 0 ? PETSC_MAX_INT : bs;
1589371c9d4SSatish Balay   bsLocal[1] = bs;
1599e47a1b0SMatthew G. Knepley   PetscCall(PetscGlobalMinMaxInt(PetscObjectComm((PetscObject)gs), bsLocal, bsMinMax));
1609371c9d4SSatish Balay   if (bsMinMax[0] != bsMinMax[1]) {
1619371c9d4SSatish Balay     bs = 1;
1629371c9d4SSatish Balay   } else {
1639371c9d4SSatish Balay     bs = bsMinMax[0];
1649371c9d4SSatish Balay   }
1659566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(subSize, &subIndices));
1669e47a1b0SMatthew G. Knepley   for (PetscInt p = pStart; p < pEnd; ++p) {
1674d9407bcSMatthew G. Knepley     PetscInt gdof, goff;
1684d9407bcSMatthew G. Knepley 
1699e47a1b0SMatthew G. Knepley     PetscCall(PetscSectionGetDof(gs, p, &gdof));
1704d9407bcSMatthew G. Knepley     if (gdof > 0) {
171dd072f5fSMatthew G. Knepley       PetscInt off = 0;
172dd072f5fSMatthew G. Knepley 
1739e47a1b0SMatthew G. Knepley       PetscCall(PetscSectionGetOffset(gs, p, &goff));
1749e47a1b0SMatthew G. Knepley       for (PetscInt f = 0; f < numFields; ++f) {
1759e47a1b0SMatthew G. Knepley         PetscInt fdof, fcdof, poff = 0;
1764d9407bcSMatthew G. Knepley 
1774d9407bcSMatthew G. Knepley         /* Can get rid of this loop by storing field information in the global section */
1789e47a1b0SMatthew G. Knepley         for (PetscInt f2 = 0; f2 < fields[f]; ++f2) {
1799e47a1b0SMatthew G. Knepley           PetscCall(PetscSectionGetFieldDof(s, p, f2, &fdof));
1809e47a1b0SMatthew G. Knepley           PetscCall(PetscSectionGetFieldConstraintDof(s, p, f2, &fcdof));
1814d9407bcSMatthew G. Knepley           poff += fdof - fcdof;
1824d9407bcSMatthew G. Knepley         }
1839e47a1b0SMatthew G. Knepley         PetscCall(PetscSectionGetFieldDof(s, p, fields[f], &fdof));
1849e47a1b0SMatthew G. Knepley         PetscCall(PetscSectionGetFieldConstraintDof(s, p, fields[f], &fcdof));
185dd072f5fSMatthew G. Knepley 
186dd072f5fSMatthew G. Knepley         if (numComps && numComps[f] >= 0) {
187dd072f5fSMatthew G. Knepley           const PetscInt *ind;
188dd072f5fSMatthew G. Knepley 
189dd072f5fSMatthew G. Knepley           // Assume sets of dofs on points are of size Nc
190dd072f5fSMatthew G. Knepley           PetscCall(PetscSectionGetFieldConstraintIndices(s, p, fields[f], &ind));
191dd072f5fSMatthew G. Knepley           for (PetscInt i = 0, fcoff = 0, pfoff = 0; i < (fdof / Nc); ++i) {
192dd072f5fSMatthew G. Knepley             for (PetscInt c = 0, fcc = 0; c < Nc; ++c) {
193dd072f5fSMatthew G. Knepley               const PetscInt k = i * Nc + c;
194dd072f5fSMatthew G. Knepley 
195dd072f5fSMatthew G. Knepley               if (ind[fcoff] == k) {
196dd072f5fSMatthew G. Knepley                 ++fcoff;
197dd072f5fSMatthew G. Knepley                 continue;
198dd072f5fSMatthew G. Knepley               }
199dd072f5fSMatthew G. Knepley               if (c == comps[off + fcc]) {
200dd072f5fSMatthew G. Knepley                 ++fcc;
201dd072f5fSMatthew G. Knepley                 subIndices[subOff++] = goff + poff + pfoff;
202dd072f5fSMatthew G. Knepley               }
203dd072f5fSMatthew G. Knepley               ++pfoff;
204dd072f5fSMatthew G. Knepley             }
205dd072f5fSMatthew G. Knepley           }
206dd072f5fSMatthew G. Knepley           off += numComps[f];
207dd072f5fSMatthew G. Knepley         } else {
2089e47a1b0SMatthew G. Knepley           for (PetscInt fc = 0; fc < fdof - fcdof; ++fc, ++subOff) subIndices[subOff] = goff + poff + fc;
2094d9407bcSMatthew G. Knepley         }
2104d9407bcSMatthew G. Knepley       }
2114d9407bcSMatthew G. Knepley     }
212dd072f5fSMatthew G. Knepley   }
2139e47a1b0SMatthew 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);
2149e47a1b0SMatthew G. Knepley   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)gs), subSize, subIndices, PETSC_OWN_POINTER, is));
215627349f5SMatthew G. Knepley   if (bs > 1) {
2169e47a1b0SMatthew G. Knepley     // We need to check that the block size does not come from non-contiguous fields
2179e47a1b0SMatthew G. Knepley     PetscInt set = 1, rset = 1;
2189e47a1b0SMatthew G. Knepley     for (PetscInt i = 0; i < subSize; i += bs) {
2199e47a1b0SMatthew G. Knepley       for (PetscInt j = 0; j < bs; ++j) {
2209371c9d4SSatish Balay         if (subIndices[i + j] != subIndices[i] + j) {
2219371c9d4SSatish Balay           set = 0;
2229371c9d4SSatish Balay           break;
2239371c9d4SSatish Balay         }
224627349f5SMatthew G. Knepley       }
225627349f5SMatthew G. Knepley     }
2269e47a1b0SMatthew G. Knepley     PetscCall(MPIU_Allreduce(&set, &rset, 1, MPIU_INT, MPI_PROD, PetscObjectComm((PetscObject)gs)));
22788206212SAlexis Marboeuf     if (rset) PetscCall(ISSetBlockSize(*is, bs));
228627349f5SMatthew G. Knepley   }
2299e47a1b0SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
2304d9407bcSMatthew G. Knepley }
2319e47a1b0SMatthew G. Knepley 
232dd072f5fSMatthew 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)
2339e47a1b0SMatthew G. Knepley {
2344d9407bcSMatthew G. Knepley   PetscSection subsection;
2354d9407bcSMatthew G. Knepley   PetscBool    haveNull = PETSC_FALSE;
2369e47a1b0SMatthew G. Knepley   PetscInt     nf = 0, of = 0;
2374d9407bcSMatthew G. Knepley 
2389e47a1b0SMatthew G. Knepley   PetscFunctionBegin;
239dd072f5fSMatthew G. Knepley   if (numComps) {
240dd072f5fSMatthew G. Knepley     const PetscInt field = fields[0];
241dd072f5fSMatthew G. Knepley 
242dd072f5fSMatthew G. Knepley     PetscCheck(numFields == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "We only support a single field for component selection right now");
243dd072f5fSMatthew G. Knepley     PetscCall(PetscSectionCreateComponentSubsection(section, numComps[field], comps, &subsection));
2449566063dSJacob Faibussowitsch     PetscCall(DMSetLocalSection(*subdm, subsection));
2459566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&subsection));
246dd072f5fSMatthew G. Knepley     (*subdm)->nullspaceConstructors[field] = dm->nullspaceConstructors[field];
247e5e52638SMatthew G. Knepley     if (dm->probs) {
2489e47a1b0SMatthew G. Knepley       PetscFV  fv, fvNew;
249dd072f5fSMatthew G. Knepley       PetscInt fnum[1] = {field};
2500f21e855SMatthew G. Knepley 
2519e47a1b0SMatthew G. Knepley       PetscCall(DMSetNumFields(*subdm, 1));
252dd072f5fSMatthew G. Knepley       PetscCall(DMGetField(dm, field, NULL, (PetscObject *)&fv));
2539e47a1b0SMatthew G. Knepley       PetscCall(PetscFVClone(fv, &fvNew));
254dd072f5fSMatthew G. Knepley       PetscCall(PetscFVSetNumComponents(fvNew, numComps[0]));
2559e47a1b0SMatthew G. Knepley       PetscCall(DMSetField(*subdm, 0, NULL, (PetscObject)fvNew));
2569e47a1b0SMatthew G. Knepley       PetscCall(PetscFVDestroy(&fvNew));
2579566063dSJacob Faibussowitsch       PetscCall(DMCreateDS(*subdm));
258dd072f5fSMatthew G. Knepley       if (numComps[0] == 1 && is) {
2590f21e855SMatthew G. Knepley         PetscObject disc, space, pmat;
2604d9407bcSMatthew G. Knepley 
261dd072f5fSMatthew G. Knepley         PetscCall(DMGetField(*subdm, field, NULL, &disc));
2629566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "nullspace", &space));
2639566063dSJacob Faibussowitsch         if (space) PetscCall(PetscObjectCompose((PetscObject)*is, "nullspace", space));
2649566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "nearnullspace", &space));
2659566063dSJacob Faibussowitsch         if (space) PetscCall(PetscObjectCompose((PetscObject)*is, "nearnullspace", space));
2669566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "pmat", &pmat));
2679566063dSJacob Faibussowitsch         if (pmat) PetscCall(PetscObjectCompose((PetscObject)*is, "pmat", pmat));
2684d9407bcSMatthew G. Knepley       }
269dd072f5fSMatthew G. Knepley       PetscCall(PetscDSCopyConstants(dm->probs[field].ds, (*subdm)->probs[0].ds));
270dd072f5fSMatthew G. Knepley       PetscCall(PetscDSCopyBoundary(dm->probs[field].ds, 1, fnum, (*subdm)->probs[0].ds));
271dd072f5fSMatthew G. Knepley       PetscCall(PetscDSSelectEquations(dm->probs[field].ds, 1, fnum, (*subdm)->probs[0].ds));
2729e47a1b0SMatthew G. Knepley     }
2739e47a1b0SMatthew G. Knepley     if ((*subdm)->nullspaceConstructors[0] && is) {
2749e47a1b0SMatthew G. Knepley       MatNullSpace nullSpace;
2759e47a1b0SMatthew G. Knepley 
2769e47a1b0SMatthew G. Knepley       PetscCall((*(*subdm)->nullspaceConstructors[0])(*subdm, 0, 0, &nullSpace));
2779e47a1b0SMatthew G. Knepley       PetscCall(PetscObjectCompose((PetscObject)*is, "nullspace", (PetscObject)nullSpace));
2789e47a1b0SMatthew G. Knepley       PetscCall(MatNullSpaceDestroy(&nullSpace));
2799e47a1b0SMatthew G. Knepley     }
2809e47a1b0SMatthew G. Knepley     if (dm->coarseMesh) PetscCall(DMCreateSubDM(dm->coarseMesh, numFields, fields, NULL, &(*subdm)->coarseMesh));
2819e47a1b0SMatthew G. Knepley     PetscFunctionReturn(PETSC_SUCCESS);
2829e47a1b0SMatthew G. Knepley   }
2839e47a1b0SMatthew G. Knepley   PetscCall(PetscSectionCreateSubsection(section, numFields, fields, &subsection));
2849e47a1b0SMatthew G. Knepley   PetscCall(DMSetLocalSection(*subdm, subsection));
2859e47a1b0SMatthew G. Knepley   PetscCall(PetscSectionDestroy(&subsection));
2869e47a1b0SMatthew G. Knepley   for (PetscInt f = 0; f < numFields; ++f) {
2879e47a1b0SMatthew G. Knepley     (*subdm)->nullspaceConstructors[f] = dm->nullspaceConstructors[fields[f]];
2889e47a1b0SMatthew G. Knepley     if ((*subdm)->nullspaceConstructors[f]) {
2899e47a1b0SMatthew G. Knepley       haveNull = PETSC_TRUE;
2909e47a1b0SMatthew G. Knepley       nf       = f;
2919e47a1b0SMatthew G. Knepley       of       = fields[f];
2929e47a1b0SMatthew G. Knepley     }
2939e47a1b0SMatthew G. Knepley   }
2949e47a1b0SMatthew G. Knepley   if (dm->probs) {
2959e47a1b0SMatthew G. Knepley     PetscCall(DMSetNumFields(*subdm, numFields));
2969e47a1b0SMatthew G. Knepley     for (PetscInt f = 0; f < numFields; ++f) {
2979e47a1b0SMatthew G. Knepley       PetscObject disc;
2989e47a1b0SMatthew G. Knepley 
2999e47a1b0SMatthew G. Knepley       PetscCall(DMGetField(dm, fields[f], NULL, &disc));
3009e47a1b0SMatthew G. Knepley       PetscCall(DMSetField(*subdm, f, NULL, disc));
3019e47a1b0SMatthew G. Knepley     }
3029e47a1b0SMatthew G. Knepley     // TODO: if only FV, then cut down the components
3039e47a1b0SMatthew G. Knepley     PetscCall(DMCreateDS(*subdm));
3049e47a1b0SMatthew G. Knepley     if (numFields == 1 && is) {
3059e47a1b0SMatthew G. Knepley       PetscObject disc, space, pmat;
3069e47a1b0SMatthew G. Knepley 
3079e47a1b0SMatthew G. Knepley       PetscCall(DMGetField(*subdm, 0, NULL, &disc));
3089e47a1b0SMatthew G. Knepley       PetscCall(PetscObjectQuery(disc, "nullspace", &space));
3099e47a1b0SMatthew G. Knepley       if (space) PetscCall(PetscObjectCompose((PetscObject)*is, "nullspace", space));
3109e47a1b0SMatthew G. Knepley       PetscCall(PetscObjectQuery(disc, "nearnullspace", &space));
3119e47a1b0SMatthew G. Knepley       if (space) PetscCall(PetscObjectCompose((PetscObject)*is, "nearnullspace", space));
3129e47a1b0SMatthew G. Knepley       PetscCall(PetscObjectQuery(disc, "pmat", &pmat));
3139e47a1b0SMatthew G. Knepley       if (pmat) PetscCall(PetscObjectCompose((PetscObject)*is, "pmat", pmat));
3149e47a1b0SMatthew G. Knepley     }
3159e47a1b0SMatthew G. Knepley     // Check if DSes record their DM fields
316e21d936aSMatthew G. Knepley     if (dm->probs[0].fields) {
3179310035eSMatthew G. Knepley       PetscInt d, e;
3189310035eSMatthew G. Knepley 
3199310035eSMatthew G. Knepley       for (d = 0, e = 0; d < dm->Nds && e < (*subdm)->Nds; ++d) {
3209310035eSMatthew G. Knepley         const PetscInt  Nf = dm->probs[d].ds->Nf;
3219310035eSMatthew G. Knepley         const PetscInt *fld;
3229310035eSMatthew G. Knepley         PetscInt        f, g;
3239310035eSMatthew G. Knepley 
3249566063dSJacob Faibussowitsch         PetscCall(ISGetIndices(dm->probs[d].fields, &fld));
3259310035eSMatthew G. Knepley         for (f = 0; f < Nf; ++f) {
3269371c9d4SSatish Balay           for (g = 0; g < numFields; ++g)
3279371c9d4SSatish Balay             if (fld[f] == fields[g]) break;
3289310035eSMatthew G. Knepley           if (g < numFields) break;
3299310035eSMatthew G. Knepley         }
3309566063dSJacob Faibussowitsch         PetscCall(ISRestoreIndices(dm->probs[d].fields, &fld));
3319310035eSMatthew G. Knepley         if (f == Nf) continue;
3329566063dSJacob Faibussowitsch         PetscCall(PetscDSCopyConstants(dm->probs[d].ds, (*subdm)->probs[e].ds));
3339566063dSJacob Faibussowitsch         PetscCall(PetscDSCopyBoundary(dm->probs[d].ds, numFields, fields, (*subdm)->probs[e].ds));
3349e47a1b0SMatthew G. Knepley         // Translate DM fields to DS fields
3359310035eSMatthew G. Knepley         {
336e21d936aSMatthew G. Knepley           IS              infields, dsfields;
33774a61aaaSMatthew G. Knepley           const PetscInt *fld, *ofld;
33874a61aaaSMatthew G. Knepley           PetscInt       *fidx;
3399e47a1b0SMatthew G. Knepley           PetscInt        onf, nf;
340e21d936aSMatthew G. Knepley 
3419566063dSJacob Faibussowitsch           PetscCall(ISCreateGeneral(PETSC_COMM_SELF, numFields, fields, PETSC_USE_POINTER, &infields));
3429566063dSJacob Faibussowitsch           PetscCall(ISIntersect(infields, dm->probs[d].fields, &dsfields));
3439566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&infields));
3449566063dSJacob Faibussowitsch           PetscCall(ISGetLocalSize(dsfields, &nf));
3457a8be351SBarry Smith           PetscCheck(nf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "DS cannot be supported on 0 fields");
3469566063dSJacob Faibussowitsch           PetscCall(ISGetIndices(dsfields, &fld));
3479566063dSJacob Faibussowitsch           PetscCall(ISGetLocalSize(dm->probs[d].fields, &onf));
3489566063dSJacob Faibussowitsch           PetscCall(ISGetIndices(dm->probs[d].fields, &ofld));
3499566063dSJacob Faibussowitsch           PetscCall(PetscMalloc1(nf, &fidx));
3509e47a1b0SMatthew G. Knepley           for (PetscInt f = 0, g = 0; f < onf && g < nf; ++f) {
35174a61aaaSMatthew G. Knepley             if (ofld[f] == fld[g]) fidx[g++] = f;
35274a61aaaSMatthew G. Knepley           }
3539566063dSJacob Faibussowitsch           PetscCall(ISRestoreIndices(dm->probs[d].fields, &ofld));
3549566063dSJacob Faibussowitsch           PetscCall(ISRestoreIndices(dsfields, &fld));
3559566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&dsfields));
356*bb4b53efSMatthew G. Knepley           PetscCall(PetscDSSelectDiscretizations(dm->probs[0].ds, nf, fidx, PETSC_DETERMINE, PETSC_DETERMINE, (*subdm)->probs[0].ds));
3579566063dSJacob Faibussowitsch           PetscCall(PetscDSSelectEquations(dm->probs[0].ds, nf, fidx, (*subdm)->probs[0].ds));
3589566063dSJacob Faibussowitsch           PetscCall(PetscFree(fidx));
3599310035eSMatthew G. Knepley         }
3609310035eSMatthew G. Knepley         ++e;
3619310035eSMatthew G. Knepley       }
362e21d936aSMatthew G. Knepley     } else {
3639566063dSJacob Faibussowitsch       PetscCall(PetscDSCopyConstants(dm->probs[0].ds, (*subdm)->probs[0].ds));
3649566063dSJacob Faibussowitsch       PetscCall(PetscDSCopyBoundary(dm->probs[0].ds, PETSC_DETERMINE, NULL, (*subdm)->probs[0].ds));
365*bb4b53efSMatthew G. Knepley       PetscCall(PetscDSSelectDiscretizations(dm->probs[0].ds, numFields, fields, PETSC_DETERMINE, PETSC_DETERMINE, (*subdm)->probs[0].ds));
3669566063dSJacob Faibussowitsch       PetscCall(PetscDSSelectEquations(dm->probs[0].ds, numFields, fields, (*subdm)->probs[0].ds));
367d5af271aSMatthew G. Knepley     }
368e21d936aSMatthew G. Knepley   }
3698cda7954SMatthew G. Knepley   if (haveNull && is) {
3708cda7954SMatthew G. Knepley     MatNullSpace nullSpace;
3718cda7954SMatthew G. Knepley 
3729566063dSJacob Faibussowitsch     PetscCall((*(*subdm)->nullspaceConstructors[nf])(*subdm, of, nf, &nullSpace));
3739566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)*is, "nullspace", (PetscObject)nullSpace));
3749566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceDestroy(&nullSpace));
3758cda7954SMatthew G. Knepley   }
37648a46eb9SPierre Jolivet   if (dm->coarseMesh) PetscCall(DMCreateSubDM(dm->coarseMesh, numFields, fields, NULL, &(*subdm)->coarseMesh));
3779e47a1b0SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
3784d9407bcSMatthew G. Knepley }
3799e47a1b0SMatthew G. Knepley 
3805d83a8b1SBarry Smith /*@
3819e47a1b0SMatthew 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`.
3829e47a1b0SMatthew G. Knepley 
3839e47a1b0SMatthew G. Knepley   Not Collective
3849e47a1b0SMatthew G. Knepley 
3859e47a1b0SMatthew G. Knepley   Input Parameters:
3869e47a1b0SMatthew G. Knepley + dm        - The `DM` object
3879e47a1b0SMatthew G. Knepley . numFields - The number of fields to incorporate into `subdm`
388dd072f5fSMatthew G. Knepley . fields    - The field numbers of the selected fields
389dd072f5fSMatthew G. Knepley . numComps  - The number of components from each field to incorporate into `subdm`, or PETSC_DECIDE for all components
390dd072f5fSMatthew G. Knepley - comps     - The component numbers of the selected fields (omitted for PTESC_DECIDE fields)
3919e47a1b0SMatthew G. Knepley 
3929e47a1b0SMatthew G. Knepley   Output Parameters:
3939e47a1b0SMatthew G. Knepley + is    - The global indices for the subproblem or `NULL`
3949e47a1b0SMatthew G. Knepley - subdm - The `DM` for the subproblem, which must already have be cloned from `dm` or `NULL`
3959e47a1b0SMatthew G. Knepley 
3969e47a1b0SMatthew G. Knepley   Level: intermediate
3979e47a1b0SMatthew G. Knepley 
3989e47a1b0SMatthew G. Knepley   Notes:
3999e47a1b0SMatthew G. Knepley   If `is` and `subdm` are both `NULL` this does nothing
4009e47a1b0SMatthew G. Knepley 
4019e47a1b0SMatthew G. Knepley .seealso: `DMCreateSubDM()`, `DMGetLocalSection()`, `DMPlexSetMigrationSF()`, `DMView()`
4029e47a1b0SMatthew G. Knepley @*/
403dd072f5fSMatthew G. Knepley PetscErrorCode DMCreateSectionSubDM(DM dm, PetscInt numFields, const PetscInt fields[], const PetscInt numComps[], const PetscInt comps[], IS *is, DM *subdm)
4049e47a1b0SMatthew G. Knepley {
4059e47a1b0SMatthew G. Knepley   PetscSection section, sectionGlobal;
4069e47a1b0SMatthew G. Knepley   PetscInt     Nf;
4079e47a1b0SMatthew G. Knepley 
4089e47a1b0SMatthew G. Knepley   PetscFunctionBegin;
4099e47a1b0SMatthew G. Knepley   if (!numFields) PetscFunctionReturn(PETSC_SUCCESS);
4109e47a1b0SMatthew G. Knepley   PetscCall(DMGetLocalSection(dm, &section));
4119e47a1b0SMatthew G. Knepley   PetscCall(DMGetGlobalSection(dm, &sectionGlobal));
4129e47a1b0SMatthew G. Knepley   PetscCheck(section, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields");
4139e47a1b0SMatthew G. Knepley   PetscCheck(sectionGlobal, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields");
4149e47a1b0SMatthew G. Knepley   PetscCall(PetscSectionGetNumFields(section, &Nf));
4159e47a1b0SMatthew 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);
4169e47a1b0SMatthew G. Knepley 
417dd072f5fSMatthew G. Knepley   if (is) PetscCall(PetscSectionSelectFields_Private(section, sectionGlobal, numFields, fields, numComps, comps, is));
418dd072f5fSMatthew G. Knepley   if (subdm) PetscCall(DMSelectFields_Private(dm, section, numFields, fields, numComps, comps, is, subdm));
4193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4204d9407bcSMatthew G. Knepley }
4212adcc780SMatthew G. Knepley 
422792b654fSMatthew G. Knepley /*@C
4230b3275a6SBarry 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`
424792b654fSMatthew G. Knepley 
42520f4b53cSBarry Smith   Not Collective
426792b654fSMatthew G. Knepley 
427d8d19677SJose E. Roman   Input Parameters:
4280b3275a6SBarry Smith + dms - The `DM` objects, the must all have the same topology; for example obtained with `DMClone()`
4290b3275a6SBarry Smith - len - The number of `DM` in `dms`
430792b654fSMatthew G. Knepley 
431792b654fSMatthew G. Knepley   Output Parameters:
43220f4b53cSBarry Smith + is      - The global indices for the subproblem, or `NULL`
4330b3275a6SBarry Smith - superdm - The `DM` for the superproblem, which must already have be cloned and contain the same topology as the `dms`
434792b654fSMatthew G. Knepley 
435792b654fSMatthew G. Knepley   Level: intermediate
436792b654fSMatthew G. Knepley 
437a4e35b19SJacob Faibussowitsch .seealso: `DMCreateSuperDM()`, `DMGetLocalSection()`, `DMPlexSetMigrationSF()`, `DMView()`
438792b654fSMatthew G. Knepley @*/
4395d83a8b1SBarry Smith PetscErrorCode DMCreateSectionSuperDM(DM dms[], PetscInt len, IS *is[], DM *superdm)
440d71ae5a4SJacob Faibussowitsch {
4419796d432SMatthew G. Knepley   MPI_Comm     comm;
4429796d432SMatthew G. Knepley   PetscSection supersection, *sections, *sectionGlobals;
4438cda7954SMatthew G. Knepley   PetscInt    *Nfs, Nf = 0, f, supf, oldf = -1, nullf = -1, i;
4449796d432SMatthew G. Knepley   PetscBool    haveNull = PETSC_FALSE;
4452adcc780SMatthew G. Knepley 
4462adcc780SMatthew G. Knepley   PetscFunctionBegin;
4479566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dms[0], &comm));
4489796d432SMatthew G. Knepley   /* Pull out local and global sections */
4499566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(len, &Nfs, len, &sections, len, &sectionGlobals));
4502adcc780SMatthew G. Knepley   for (i = 0; i < len; ++i) {
4519566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dms[i], &sections[i]));
4529566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalSection(dms[i], &sectionGlobals[i]));
4537a8be351SBarry Smith     PetscCheck(sections[i], comm, PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields");
4547a8be351SBarry Smith     PetscCheck(sectionGlobals[i], comm, PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields");
4559566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetNumFields(sections[i], &Nfs[i]));
4562adcc780SMatthew G. Knepley     Nf += Nfs[i];
4572adcc780SMatthew G. Knepley   }
4589796d432SMatthew G. Knepley   /* Create the supersection */
4599566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreateSupersection(sections, len, &supersection));
4609566063dSJacob Faibussowitsch   PetscCall(DMSetLocalSection(*superdm, supersection));
4619796d432SMatthew G. Knepley   /* Create ISes */
4622adcc780SMatthew G. Knepley   if (is) {
4639796d432SMatthew G. Knepley     PetscSection supersectionGlobal;
4649796d432SMatthew G. Knepley     PetscInt     bs = -1, startf = 0;
4652adcc780SMatthew G. Knepley 
4669566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(len, is));
4679566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalSection(*superdm, &supersectionGlobal));
4689796d432SMatthew G. Knepley     for (i = 0; i < len; startf += Nfs[i], ++i) {
4699796d432SMatthew G. Knepley       PetscInt *subIndices;
470ec4c761aSMatthew G. Knepley       PetscInt  subSize, subOff, pStart, pEnd, p, start, end, dummy;
4712adcc780SMatthew G. Knepley 
4729566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetChart(sectionGlobals[i], &pStart, &pEnd));
4739566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetConstrainedStorageSize(sectionGlobals[i], &subSize));
4749566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(subSize, &subIndices));
4759796d432SMatthew G. Knepley       for (p = pStart, subOff = 0; p < pEnd; ++p) {
4769796d432SMatthew G. Knepley         PetscInt gdof, gcdof, gtdof, d;
4772adcc780SMatthew G. Knepley 
4789566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(sectionGlobals[i], p, &gdof));
4799566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetConstraintDof(sections[i], p, &gcdof));
4802adcc780SMatthew G. Knepley         gtdof = gdof - gcdof;
4812adcc780SMatthew G. Knepley         if (gdof > 0 && gtdof) {
4829371c9d4SSatish Balay           if (bs < 0) {
4839371c9d4SSatish Balay             bs = gtdof;
4849371c9d4SSatish Balay           } else if (bs != gtdof) {
4859371c9d4SSatish Balay             bs = 1;
4869371c9d4SSatish Balay           }
4879566063dSJacob Faibussowitsch           PetscCall(DMGetGlobalFieldOffset_Private(*superdm, p, startf, &start, &dummy));
4889566063dSJacob Faibussowitsch           PetscCall(DMGetGlobalFieldOffset_Private(*superdm, p, startf + Nfs[i] - 1, &dummy, &end));
48963a3b9bcSJacob 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);
4909796d432SMatthew G. Knepley           for (d = start; d < end; ++d, ++subOff) subIndices[subOff] = d;
4912adcc780SMatthew G. Knepley         }
4922adcc780SMatthew G. Knepley       }
4939566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(comm, subSize, subIndices, PETSC_OWN_POINTER, &(*is)[i]));
4942adcc780SMatthew G. Knepley       /* Must have same blocksize on all procs (some might have no points) */
4959796d432SMatthew G. Knepley       {
4969796d432SMatthew G. Knepley         PetscInt bs = -1, bsLocal[2], bsMinMax[2];
4979796d432SMatthew G. Knepley 
4989371c9d4SSatish Balay         bsLocal[0] = bs < 0 ? PETSC_MAX_INT : bs;
4999371c9d4SSatish Balay         bsLocal[1] = bs;
5009566063dSJacob Faibussowitsch         PetscCall(PetscGlobalMinMaxInt(comm, bsLocal, bsMinMax));
5019371c9d4SSatish Balay         if (bsMinMax[0] != bsMinMax[1]) {
5029371c9d4SSatish Balay           bs = 1;
5039371c9d4SSatish Balay         } else {
5049371c9d4SSatish Balay           bs = bsMinMax[0];
5059371c9d4SSatish Balay         }
5069566063dSJacob Faibussowitsch         PetscCall(ISSetBlockSize((*is)[i], bs));
5072adcc780SMatthew G. Knepley       }
5082adcc780SMatthew G. Knepley     }
5092adcc780SMatthew G. Knepley   }
5109796d432SMatthew G. Knepley   /* Preserve discretizations */
511e5e52638SMatthew G. Knepley   if (len && dms[0]->probs) {
5129566063dSJacob Faibussowitsch     PetscCall(DMSetNumFields(*superdm, Nf));
5139796d432SMatthew G. Knepley     for (i = 0, supf = 0; i < len; ++i) {
5149796d432SMatthew G. Knepley       for (f = 0; f < Nfs[i]; ++f, ++supf) {
5152adcc780SMatthew G. Knepley         PetscObject disc;
5162adcc780SMatthew G. Knepley 
5179566063dSJacob Faibussowitsch         PetscCall(DMGetField(dms[i], f, NULL, &disc));
5189566063dSJacob Faibussowitsch         PetscCall(DMSetField(*superdm, supf, NULL, disc));
5192adcc780SMatthew G. Knepley       }
5202adcc780SMatthew G. Knepley     }
5219566063dSJacob Faibussowitsch     PetscCall(DMCreateDS(*superdm));
5222adcc780SMatthew G. Knepley   }
5239796d432SMatthew G. Knepley   /* Preserve nullspaces */
5249796d432SMatthew G. Knepley   for (i = 0, supf = 0; i < len; ++i) {
5259796d432SMatthew G. Knepley     for (f = 0; f < Nfs[i]; ++f, ++supf) {
5269796d432SMatthew G. Knepley       (*superdm)->nullspaceConstructors[supf] = dms[i]->nullspaceConstructors[f];
5279796d432SMatthew G. Knepley       if ((*superdm)->nullspaceConstructors[supf]) {
5289796d432SMatthew G. Knepley         haveNull = PETSC_TRUE;
5299796d432SMatthew G. Knepley         nullf    = supf;
5308cda7954SMatthew G. Knepley         oldf     = f;
5312adcc780SMatthew G. Knepley       }
5329796d432SMatthew G. Knepley     }
5339796d432SMatthew G. Knepley   }
5349796d432SMatthew G. Knepley   /* Attach nullspace to IS */
5359796d432SMatthew G. Knepley   if (haveNull && is) {
5369796d432SMatthew G. Knepley     MatNullSpace nullSpace;
5379796d432SMatthew G. Knepley 
5389566063dSJacob Faibussowitsch     PetscCall((*(*superdm)->nullspaceConstructors[nullf])(*superdm, oldf, nullf, &nullSpace));
5399566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)(*is)[nullf], "nullspace", (PetscObject)nullSpace));
5409566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceDestroy(&nullSpace));
5419796d432SMatthew G. Knepley   }
5429566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&supersection));
5439566063dSJacob Faibussowitsch   PetscCall(PetscFree3(Nfs, sections, sectionGlobals));
5443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5452adcc780SMatthew G. Knepley }
546