1af0996ceSBarry Smith #include <petsc/private/dmimpl.h> /*I "petscdm.h" I*/
22764a2aaSMatthew G. Knepley #include <petscds.h>
311689aebSJed Brown
DMCreateGlobalVector_Section_Private(DM dm,Vec * vec)4d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateGlobalVector_Section_Private(DM dm, Vec *vec)
5d71ae5a4SJacob Faibussowitsch {
611689aebSJed Brown PetscSection gSection;
7cc30b04dSMatthew G Knepley PetscInt localSize, bs, blockSize = -1, pStart, pEnd, p;
85d1c0279SHong Zhang PetscInt in[2], out[2];
911689aebSJed Brown
1011689aebSJed Brown PetscFunctionBegin;
119566063dSJacob Faibussowitsch PetscCall(DMGetGlobalSection(dm, &gSection));
129566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(gSection, &pStart, &pEnd));
1311689aebSJed Brown for (p = pStart; p < pEnd; ++p) {
1411689aebSJed Brown PetscInt dof, cdof;
1511689aebSJed Brown
169566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(gSection, p, &dof));
179566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(gSection, p, &cdof));
180680ce0aSHong Zhang
195b8243e1SJed Brown if (dof - cdof > 0) {
205b8243e1SJed Brown if (blockSize < 0) {
210680ce0aSHong Zhang /* set blockSize */
220680ce0aSHong Zhang blockSize = dof - cdof;
235b8243e1SJed Brown } else {
245b8243e1SJed Brown blockSize = PetscGCD(dof - cdof, blockSize);
2511689aebSJed Brown }
2611689aebSJed Brown }
270680ce0aSHong Zhang }
285d1c0279SHong Zhang
291690c2aeSBarry Smith // You cannot negate PETSC_INT_MIN
301690c2aeSBarry Smith in[0] = blockSize < 0 ? -PETSC_INT_MAX : -blockSize;
315d1c0279SHong Zhang in[1] = blockSize;
32462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(in, out, 2, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm)));
335d1c0279SHong Zhang /* -out[0] = min(blockSize), out[1] = max(blockSize) */
345d1c0279SHong Zhang if (-out[0] == out[1]) {
355d1c0279SHong Zhang bs = out[1];
365d1c0279SHong Zhang } else bs = 1;
375d1c0279SHong Zhang
385d1c0279SHong Zhang if (bs < 0) { /* Everyone was empty */
395d1c0279SHong Zhang blockSize = 1;
405d1c0279SHong Zhang bs = 1;
415d1c0279SHong Zhang }
425d1c0279SHong Zhang
439566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstrainedStorageSize(gSection, &localSize));
447a8be351SBarry 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);
459566063dSJacob Faibussowitsch PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), vec));
469566063dSJacob Faibussowitsch PetscCall(VecSetSizes(*vec, localSize, PETSC_DETERMINE));
479566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(*vec, bs));
489566063dSJacob Faibussowitsch PetscCall(VecSetType(*vec, dm->vectype));
499566063dSJacob Faibussowitsch PetscCall(VecSetDM(*vec, dm));
509566063dSJacob Faibussowitsch /* PetscCall(VecSetLocalToGlobalMapping(*vec, dm->ltogmap)); */
513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
5211689aebSJed Brown }
5311689aebSJed Brown
DMCreateLocalVector_Section_Private(DM dm,Vec * vec)54d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateLocalVector_Section_Private(DM dm, Vec *vec)
55d71ae5a4SJacob Faibussowitsch {
56fad22124SMatthew G Knepley PetscSection section;
5711689aebSJed Brown PetscInt localSize, blockSize = -1, pStart, pEnd, p;
5811689aebSJed Brown
5911689aebSJed Brown PetscFunctionBegin;
609566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, §ion));
619566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
6211689aebSJed Brown for (p = pStart; p < pEnd; ++p) {
6311689aebSJed Brown PetscInt dof;
6411689aebSJed Brown
659566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, p, &dof));
6611689aebSJed Brown if ((blockSize < 0) && (dof > 0)) blockSize = dof;
675b8243e1SJed Brown if (dof > 0) blockSize = PetscGCD(dof, blockSize);
6811689aebSJed Brown }
699566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(section, &localSize));
709566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, vec));
719566063dSJacob Faibussowitsch PetscCall(VecSetSizes(*vec, localSize, localSize));
7258b7e2c1SStefano Zampini PetscCall(VecSetBlockSize(*vec, PetscAbs(blockSize)));
739566063dSJacob Faibussowitsch PetscCall(VecSetType(*vec, dm->vectype));
749566063dSJacob Faibussowitsch PetscCall(VecSetDM(*vec, dm));
753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
7611689aebSJed Brown }
774d9407bcSMatthew G. Knepley
PetscSectionSelectFields_Private(PetscSection s,PetscSection gs,PetscInt numFields,const PetscInt fields[],const PetscInt numComps[],const PetscInt comps[],IS * is)78dd072f5fSMatthew G. Knepley static PetscErrorCode PetscSectionSelectFields_Private(PetscSection s, PetscSection gs, PetscInt numFields, const PetscInt fields[], const PetscInt numComps[], const PetscInt comps[], IS *is)
79d71ae5a4SJacob Faibussowitsch {
8055f9ef1dSMatthew G. Knepley IS permutation;
8155f9ef1dSMatthew G. Knepley const PetscInt *perm = NULL;
824d9407bcSMatthew G. Knepley PetscInt *subIndices;
83a5e5a98fSStefano Zampini PetscInt mbs, bs = 0, bsLocal[2], bsMinMax[2];
84dd072f5fSMatthew G. Knepley PetscInt pStart, pEnd, Nc, subSize = 0, subOff = 0;
854d9407bcSMatthew G. Knepley
864d9407bcSMatthew G. Knepley PetscFunctionBegin;
879e47a1b0SMatthew G. Knepley PetscCall(PetscSectionGetChart(gs, &pStart, &pEnd));
8855f9ef1dSMatthew G. Knepley PetscCall(PetscSectionGetPermutation(s, &permutation));
8955f9ef1dSMatthew G. Knepley if (permutation) PetscCall(ISGetIndices(permutation, &perm));
90dd072f5fSMatthew G. Knepley if (numComps) {
91dd072f5fSMatthew G. Knepley for (PetscInt f = 0, off = 0; f < numFields; ++f) {
923a544194SStefano Zampini PetscInt Nc;
933a544194SStefano Zampini
94dd072f5fSMatthew G. Knepley if (numComps[f] < 0) continue;
95dd072f5fSMatthew G. Knepley PetscCall(PetscSectionGetFieldComponents(s, f, &Nc));
96dd072f5fSMatthew 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);
97dd072f5fSMatthew 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);
98dd072f5fSMatthew G. Knepley bs += numComps[f];
999e47a1b0SMatthew G. Knepley }
100dd072f5fSMatthew G. Knepley } else {
1019e47a1b0SMatthew G. Knepley for (PetscInt f = 0; f < numFields; ++f) {
1029e47a1b0SMatthew G. Knepley PetscInt Nc;
1039e47a1b0SMatthew G. Knepley
1049e47a1b0SMatthew G. Knepley PetscCall(PetscSectionGetFieldComponents(s, fields[f], &Nc));
1053a544194SStefano Zampini bs += Nc;
1063a544194SStefano Zampini }
107dd072f5fSMatthew G. Knepley }
108a5e5a98fSStefano Zampini mbs = -1; /* multiple of block size not set */
1099e47a1b0SMatthew G. Knepley for (PetscInt p = pStart; p < pEnd; ++p) {
11055f9ef1dSMatthew G. Knepley const PetscInt point = perm ? perm[p - pStart] : p;
111b9eca265Ssarens PetscInt gdof, pSubSize = 0;
1124d9407bcSMatthew G. Knepley
11355f9ef1dSMatthew G. Knepley PetscCall(PetscSectionGetDof(gs, point, &gdof));
1144d9407bcSMatthew G. Knepley if (gdof > 0) {
115dd072f5fSMatthew G. Knepley PetscInt off = 0;
1164d9407bcSMatthew G. Knepley
117dd072f5fSMatthew G. Knepley for (PetscInt f = 0; f < numFields; ++f) {
118dd072f5fSMatthew G. Knepley PetscInt fdof, fcdof, sfdof, sfcdof = 0;
119dd072f5fSMatthew G. Knepley
120dd072f5fSMatthew G. Knepley PetscCall(PetscSectionGetFieldComponents(s, f, &Nc));
12155f9ef1dSMatthew G. Knepley PetscCall(PetscSectionGetFieldDof(s, point, fields[f], &fdof));
12255f9ef1dSMatthew G. Knepley PetscCall(PetscSectionGetFieldConstraintDof(s, point, fields[f], &fcdof));
123dd072f5fSMatthew G. Knepley if (numComps && numComps[f] >= 0) {
124dd072f5fSMatthew G. Knepley const PetscInt *ind;
125dd072f5fSMatthew G. Knepley
126dd072f5fSMatthew G. Knepley // Assume sets of dofs on points are of size Nc
12755f9ef1dSMatthew 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);
128dd072f5fSMatthew G. Knepley sfdof = (fdof / Nc) * numComps[f];
12955f9ef1dSMatthew G. Knepley PetscCall(PetscSectionGetFieldConstraintIndices(s, point, fields[f], &ind));
130dd072f5fSMatthew G. Knepley for (PetscInt i = 0; i < (fdof / Nc); ++i) {
131dd072f5fSMatthew G. Knepley for (PetscInt c = 0, fcc = 0; c < Nc; ++c) {
132dd072f5fSMatthew G. Knepley if (c == comps[off + fcc]) {
133dd072f5fSMatthew G. Knepley ++fcc;
134dd072f5fSMatthew G. Knepley ++sfcdof;
135dd072f5fSMatthew G. Knepley }
136dd072f5fSMatthew G. Knepley }
137dd072f5fSMatthew G. Knepley }
138dd072f5fSMatthew G. Knepley pSubSize += sfdof - sfcdof;
139dd072f5fSMatthew G. Knepley off += numComps[f];
140dd072f5fSMatthew G. Knepley } else {
141b9eca265Ssarens pSubSize += fdof - fcdof;
142b9eca265Ssarens }
143dd072f5fSMatthew G. Knepley }
144b9eca265Ssarens subSize += pSubSize;
145a5e5a98fSStefano Zampini if (pSubSize && pSubSize % bs) {
146a5e5a98fSStefano Zampini // Layout does not admit a pointwise block size -> set mbs to 0
147a5e5a98fSStefano Zampini mbs = 0;
148a5e5a98fSStefano Zampini } else if (pSubSize) {
149a5e5a98fSStefano Zampini if (mbs == -1) mbs = pSubSize / bs;
150a5e5a98fSStefano Zampini else mbs = PetscMin(mbs, pSubSize / bs);
1514d9407bcSMatthew G. Knepley }
1524d9407bcSMatthew G. Knepley }
1534d9407bcSMatthew G. Knepley }
154a5e5a98fSStefano Zampini
1559e47a1b0SMatthew G. Knepley // Must have same blocksize on all procs (some might have no points)
156a5e5a98fSStefano Zampini bsLocal[0] = mbs < 0 ? PETSC_INT_MAX : mbs;
157a5e5a98fSStefano Zampini bsLocal[1] = mbs;
1589e47a1b0SMatthew G. Knepley PetscCall(PetscGlobalMinMaxInt(PetscObjectComm((PetscObject)gs), bsLocal, bsMinMax));
159a5e5a98fSStefano Zampini if (bsMinMax[0] != bsMinMax[1]) { /* different multiple of block size -> set bs to 1 */
1609371c9d4SSatish Balay bs = 1;
161a5e5a98fSStefano Zampini } else { /* same multiple */
162a5e5a98fSStefano Zampini mbs = bsMinMax[0];
163a5e5a98fSStefano Zampini bs *= mbs;
1649371c9d4SSatish Balay }
1659566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(subSize, &subIndices));
1669e47a1b0SMatthew G. Knepley for (PetscInt p = pStart; p < pEnd; ++p) {
16755f9ef1dSMatthew G. Knepley const PetscInt point = perm ? perm[p - pStart] : p;
1684d9407bcSMatthew G. Knepley PetscInt gdof, goff;
1694d9407bcSMatthew G. Knepley
17055f9ef1dSMatthew G. Knepley PetscCall(PetscSectionGetDof(gs, point, &gdof));
1714d9407bcSMatthew G. Knepley if (gdof > 0) {
172dd072f5fSMatthew G. Knepley PetscInt off = 0;
173dd072f5fSMatthew G. Knepley
17455f9ef1dSMatthew G. Knepley PetscCall(PetscSectionGetOffset(gs, point, &goff));
1759e47a1b0SMatthew G. Knepley for (PetscInt f = 0; f < numFields; ++f) {
1769e47a1b0SMatthew G. Knepley PetscInt fdof, fcdof, poff = 0;
1774d9407bcSMatthew G. Knepley
1784d9407bcSMatthew G. Knepley /* Can get rid of this loop by storing field information in the global section */
1799e47a1b0SMatthew G. Knepley for (PetscInt f2 = 0; f2 < fields[f]; ++f2) {
18055f9ef1dSMatthew G. Knepley PetscCall(PetscSectionGetFieldDof(s, point, f2, &fdof));
18155f9ef1dSMatthew G. Knepley PetscCall(PetscSectionGetFieldConstraintDof(s, point, f2, &fcdof));
1824d9407bcSMatthew G. Knepley poff += fdof - fcdof;
1834d9407bcSMatthew G. Knepley }
18455f9ef1dSMatthew G. Knepley PetscCall(PetscSectionGetFieldDof(s, point, fields[f], &fdof));
18555f9ef1dSMatthew G. Knepley PetscCall(PetscSectionGetFieldConstraintDof(s, point, fields[f], &fcdof));
186dd072f5fSMatthew G. Knepley
187dd072f5fSMatthew G. Knepley if (numComps && numComps[f] >= 0) {
188dd072f5fSMatthew G. Knepley const PetscInt *ind;
189dd072f5fSMatthew G. Knepley
190dd072f5fSMatthew G. Knepley // Assume sets of dofs on points are of size Nc
19155f9ef1dSMatthew G. Knepley PetscCall(PetscSectionGetFieldConstraintIndices(s, point, fields[f], &ind));
192dd072f5fSMatthew G. Knepley for (PetscInt i = 0, fcoff = 0, pfoff = 0; i < (fdof / Nc); ++i) {
193dd072f5fSMatthew G. Knepley for (PetscInt c = 0, fcc = 0; c < Nc; ++c) {
194dd072f5fSMatthew G. Knepley const PetscInt k = i * Nc + c;
195dd072f5fSMatthew G. Knepley
196dd072f5fSMatthew G. Knepley if (ind[fcoff] == k) {
197dd072f5fSMatthew G. Knepley ++fcoff;
198dd072f5fSMatthew G. Knepley continue;
199dd072f5fSMatthew G. Knepley }
200dd072f5fSMatthew G. Knepley if (c == comps[off + fcc]) {
201dd072f5fSMatthew G. Knepley ++fcc;
202dd072f5fSMatthew G. Knepley subIndices[subOff++] = goff + poff + pfoff;
203dd072f5fSMatthew G. Knepley }
204dd072f5fSMatthew G. Knepley ++pfoff;
205dd072f5fSMatthew G. Knepley }
206dd072f5fSMatthew G. Knepley }
207dd072f5fSMatthew G. Knepley off += numComps[f];
208dd072f5fSMatthew G. Knepley } else {
2099e47a1b0SMatthew G. Knepley for (PetscInt fc = 0; fc < fdof - fcdof; ++fc, ++subOff) subIndices[subOff] = goff + poff + fc;
2104d9407bcSMatthew G. Knepley }
2114d9407bcSMatthew G. Knepley }
2124d9407bcSMatthew G. Knepley }
213dd072f5fSMatthew G. Knepley }
21455f9ef1dSMatthew G. Knepley if (permutation) PetscCall(ISRestoreIndices(permutation, &perm));
2159e47a1b0SMatthew 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);
2169e47a1b0SMatthew G. Knepley PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)gs), subSize, subIndices, PETSC_OWN_POINTER, is));
217627349f5SMatthew G. Knepley if (bs > 1) {
2189e47a1b0SMatthew G. Knepley // We need to check that the block size does not come from non-contiguous fields
2199e47a1b0SMatthew G. Knepley PetscInt set = 1, rset = 1;
2209e47a1b0SMatthew G. Knepley for (PetscInt i = 0; i < subSize; i += bs) {
2219e47a1b0SMatthew G. Knepley for (PetscInt j = 0; j < bs; ++j) {
2229371c9d4SSatish Balay if (subIndices[i + j] != subIndices[i] + j) {
2239371c9d4SSatish Balay set = 0;
2249371c9d4SSatish Balay break;
2259371c9d4SSatish Balay }
226627349f5SMatthew G. Knepley }
227627349f5SMatthew G. Knepley }
228462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&set, &rset, 1, MPIU_INT, MPI_PROD, PetscObjectComm((PetscObject)gs)));
22988206212SAlexis Marboeuf if (rset) PetscCall(ISSetBlockSize(*is, bs));
230627349f5SMatthew G. Knepley }
2319e47a1b0SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
2324d9407bcSMatthew G. Knepley }
2339e47a1b0SMatthew G. Knepley
DMSelectFields_Private(DM dm,PetscSection section,PetscInt numFields,const PetscInt fields[],const PetscInt numComps[],const PetscInt comps[],IS * is,DM * subdm)234dd072f5fSMatthew 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)
2359e47a1b0SMatthew G. Knepley {
2364d9407bcSMatthew G. Knepley PetscSection subsection;
2374d9407bcSMatthew G. Knepley PetscBool haveNull = PETSC_FALSE;
2389e47a1b0SMatthew G. Knepley PetscInt nf = 0, of = 0;
2394d9407bcSMatthew G. Knepley
2409e47a1b0SMatthew G. Knepley PetscFunctionBegin;
241*4758e3ceSMatthew G. Knepley // Create nullspace constructor slots
242*4758e3ceSMatthew G. Knepley if (dm->nullspaceConstructors) {
243*4758e3ceSMatthew G. Knepley PetscCall(PetscFree2((*subdm)->nullspaceConstructors, (*subdm)->nearnullspaceConstructors));
244*4758e3ceSMatthew G. Knepley PetscCall(PetscCalloc2(numFields, &(*subdm)->nullspaceConstructors, numFields, &(*subdm)->nearnullspaceConstructors));
245*4758e3ceSMatthew G. Knepley }
246dd072f5fSMatthew G. Knepley if (numComps) {
247dd072f5fSMatthew G. Knepley const PetscInt field = fields[0];
248dd072f5fSMatthew G. Knepley
249dd072f5fSMatthew G. Knepley PetscCheck(numFields == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "We only support a single field for component selection right now");
250dd072f5fSMatthew G. Knepley PetscCall(PetscSectionCreateComponentSubsection(section, numComps[field], comps, &subsection));
2519566063dSJacob Faibussowitsch PetscCall(DMSetLocalSection(*subdm, subsection));
2529566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&subsection));
253*4758e3ceSMatthew G. Knepley if (dm->nullspaceConstructors) (*subdm)->nullspaceConstructors[field] = dm->nullspaceConstructors[field];
254e5e52638SMatthew G. Knepley if (dm->probs) {
2559e47a1b0SMatthew G. Knepley PetscFV fv, fvNew;
256dd072f5fSMatthew G. Knepley PetscInt fnum[1] = {field};
2570f21e855SMatthew G. Knepley
2589e47a1b0SMatthew G. Knepley PetscCall(DMSetNumFields(*subdm, 1));
259dd072f5fSMatthew G. Knepley PetscCall(DMGetField(dm, field, NULL, (PetscObject *)&fv));
2609e47a1b0SMatthew G. Knepley PetscCall(PetscFVClone(fv, &fvNew));
261dd072f5fSMatthew G. Knepley PetscCall(PetscFVSetNumComponents(fvNew, numComps[0]));
2629e47a1b0SMatthew G. Knepley PetscCall(DMSetField(*subdm, 0, NULL, (PetscObject)fvNew));
2639e47a1b0SMatthew G. Knepley PetscCall(PetscFVDestroy(&fvNew));
2649566063dSJacob Faibussowitsch PetscCall(DMCreateDS(*subdm));
265dd072f5fSMatthew G. Knepley if (numComps[0] == 1 && is) {
2660f21e855SMatthew G. Knepley PetscObject disc, space, pmat;
2674d9407bcSMatthew G. Knepley
268dd072f5fSMatthew G. Knepley PetscCall(DMGetField(*subdm, field, NULL, &disc));
2699566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(disc, "nullspace", &space));
2709566063dSJacob Faibussowitsch if (space) PetscCall(PetscObjectCompose((PetscObject)*is, "nullspace", space));
2719566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(disc, "nearnullspace", &space));
2729566063dSJacob Faibussowitsch if (space) PetscCall(PetscObjectCompose((PetscObject)*is, "nearnullspace", space));
2739566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(disc, "pmat", &pmat));
2749566063dSJacob Faibussowitsch if (pmat) PetscCall(PetscObjectCompose((PetscObject)*is, "pmat", pmat));
2754d9407bcSMatthew G. Knepley }
276dd072f5fSMatthew G. Knepley PetscCall(PetscDSCopyConstants(dm->probs[field].ds, (*subdm)->probs[0].ds));
277dd072f5fSMatthew G. Knepley PetscCall(PetscDSCopyBoundary(dm->probs[field].ds, 1, fnum, (*subdm)->probs[0].ds));
278dd072f5fSMatthew G. Knepley PetscCall(PetscDSSelectEquations(dm->probs[field].ds, 1, fnum, (*subdm)->probs[0].ds));
2799e47a1b0SMatthew G. Knepley }
280*4758e3ceSMatthew G. Knepley if ((*subdm)->nullspaceConstructors && (*subdm)->nullspaceConstructors[0] && is) {
2819e47a1b0SMatthew G. Knepley MatNullSpace nullSpace;
2829e47a1b0SMatthew G. Knepley
2839e47a1b0SMatthew G. Knepley PetscCall((*(*subdm)->nullspaceConstructors[0])(*subdm, 0, 0, &nullSpace));
2849e47a1b0SMatthew G. Knepley PetscCall(PetscObjectCompose((PetscObject)*is, "nullspace", (PetscObject)nullSpace));
2859e47a1b0SMatthew G. Knepley PetscCall(MatNullSpaceDestroy(&nullSpace));
2869e47a1b0SMatthew G. Knepley }
2879e47a1b0SMatthew G. Knepley if (dm->coarseMesh) PetscCall(DMCreateSubDM(dm->coarseMesh, numFields, fields, NULL, &(*subdm)->coarseMesh));
2889e47a1b0SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
2899e47a1b0SMatthew G. Knepley }
2909e47a1b0SMatthew G. Knepley PetscCall(PetscSectionCreateSubsection(section, numFields, fields, &subsection));
2919e47a1b0SMatthew G. Knepley PetscCall(DMSetLocalSection(*subdm, subsection));
2929e47a1b0SMatthew G. Knepley PetscCall(PetscSectionDestroy(&subsection));
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));
355bb4b53efSMatthew G. Knepley PetscCall(PetscDSSelectDiscretizations(dm->probs[0].ds, nf, fidx, PETSC_DETERMINE, PETSC_DETERMINE, (*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));
364bb4b53efSMatthew G. Knepley PetscCall(PetscDSSelectDiscretizations(dm->probs[0].ds, numFields, fields, PETSC_DETERMINE, PETSC_DETERMINE, (*subdm)->probs[0].ds));
3659566063dSJacob Faibussowitsch PetscCall(PetscDSSelectEquations(dm->probs[0].ds, numFields, fields, (*subdm)->probs[0].ds));
366d5af271aSMatthew G. Knepley }
367e21d936aSMatthew G. Knepley }
368*4758e3ceSMatthew G. Knepley for (PetscInt f = 0; f < numFields; ++f) {
369*4758e3ceSMatthew G. Knepley if (dm->nullspaceConstructors) {
370*4758e3ceSMatthew G. Knepley (*subdm)->nullspaceConstructors[f] = dm->nullspaceConstructors[fields[f]];
371*4758e3ceSMatthew G. Knepley if ((*subdm)->nullspaceConstructors[f]) {
372*4758e3ceSMatthew G. Knepley haveNull = PETSC_TRUE;
373*4758e3ceSMatthew G. Knepley nf = f;
374*4758e3ceSMatthew G. Knepley of = fields[f];
375*4758e3ceSMatthew G. Knepley }
376*4758e3ceSMatthew G. Knepley }
377*4758e3ceSMatthew G. Knepley }
3788cda7954SMatthew G. Knepley if (haveNull && is) {
3798cda7954SMatthew G. Knepley MatNullSpace nullSpace;
3808cda7954SMatthew G. Knepley
3819566063dSJacob Faibussowitsch PetscCall((*(*subdm)->nullspaceConstructors[nf])(*subdm, of, nf, &nullSpace));
3829566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)*is, "nullspace", (PetscObject)nullSpace));
3839566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&nullSpace));
3848cda7954SMatthew G. Knepley }
38548a46eb9SPierre Jolivet if (dm->coarseMesh) PetscCall(DMCreateSubDM(dm->coarseMesh, numFields, fields, NULL, &(*subdm)->coarseMesh));
3869e47a1b0SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
3874d9407bcSMatthew G. Knepley }
3889e47a1b0SMatthew G. Knepley
3895d83a8b1SBarry Smith /*@
3909e47a1b0SMatthew 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`.
3919e47a1b0SMatthew G. Knepley
3929e47a1b0SMatthew G. Knepley Not Collective
3939e47a1b0SMatthew G. Knepley
3949e47a1b0SMatthew G. Knepley Input Parameters:
3959e47a1b0SMatthew G. Knepley + dm - The `DM` object
3969e47a1b0SMatthew G. Knepley . numFields - The number of fields to incorporate into `subdm`
397dd072f5fSMatthew G. Knepley . fields - The field numbers of the selected fields
398dd072f5fSMatthew G. Knepley . numComps - The number of components from each field to incorporate into `subdm`, or PETSC_DECIDE for all components
399dd072f5fSMatthew G. Knepley - comps - The component numbers of the selected fields (omitted for PTESC_DECIDE fields)
4009e47a1b0SMatthew G. Knepley
4019e47a1b0SMatthew G. Knepley Output Parameters:
4029e47a1b0SMatthew G. Knepley + is - The global indices for the subproblem or `NULL`
4039e47a1b0SMatthew G. Knepley - subdm - The `DM` for the subproblem, which must already have be cloned from `dm` or `NULL`
4049e47a1b0SMatthew G. Knepley
4059e47a1b0SMatthew G. Knepley Level: intermediate
4069e47a1b0SMatthew G. Knepley
4079e47a1b0SMatthew G. Knepley Notes:
4089e47a1b0SMatthew G. Knepley If `is` and `subdm` are both `NULL` this does nothing
4099e47a1b0SMatthew G. Knepley
4109e47a1b0SMatthew G. Knepley .seealso: `DMCreateSubDM()`, `DMGetLocalSection()`, `DMPlexSetMigrationSF()`, `DMView()`
4119e47a1b0SMatthew G. Knepley @*/
DMCreateSectionSubDM(DM dm,PetscInt numFields,const PetscInt fields[],const PetscInt numComps[],const PetscInt comps[],IS * is,DM * subdm)412dd072f5fSMatthew G. Knepley PetscErrorCode DMCreateSectionSubDM(DM dm, PetscInt numFields, const PetscInt fields[], const PetscInt numComps[], const PetscInt comps[], IS *is, DM *subdm)
4139e47a1b0SMatthew G. Knepley {
4149e47a1b0SMatthew G. Knepley PetscSection section, sectionGlobal;
4159e47a1b0SMatthew G. Knepley PetscInt Nf;
4169e47a1b0SMatthew G. Knepley
4179e47a1b0SMatthew G. Knepley PetscFunctionBegin;
4189e47a1b0SMatthew G. Knepley if (!numFields) PetscFunctionReturn(PETSC_SUCCESS);
4199e47a1b0SMatthew G. Knepley PetscCall(DMGetLocalSection(dm, §ion));
4209e47a1b0SMatthew G. Knepley PetscCall(DMGetGlobalSection(dm, §ionGlobal));
4219e47a1b0SMatthew G. Knepley PetscCheck(section, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields");
4229e47a1b0SMatthew G. Knepley PetscCheck(sectionGlobal, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields");
4239e47a1b0SMatthew G. Knepley PetscCall(PetscSectionGetNumFields(section, &Nf));
4249e47a1b0SMatthew 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);
4259e47a1b0SMatthew G. Knepley
426dd072f5fSMatthew G. Knepley if (is) PetscCall(PetscSectionSelectFields_Private(section, sectionGlobal, numFields, fields, numComps, comps, is));
427dd072f5fSMatthew G. Knepley if (subdm) PetscCall(DMSelectFields_Private(dm, section, numFields, fields, numComps, comps, is, subdm));
4283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
4294d9407bcSMatthew G. Knepley }
4302adcc780SMatthew G. Knepley
431792b654fSMatthew G. Knepley /*@C
4320b3275a6SBarry 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`
433792b654fSMatthew G. Knepley
43420f4b53cSBarry Smith Not Collective
435792b654fSMatthew G. Knepley
436d8d19677SJose E. Roman Input Parameters:
4370b3275a6SBarry Smith + dms - The `DM` objects, the must all have the same topology; for example obtained with `DMClone()`
4380b3275a6SBarry Smith - len - The number of `DM` in `dms`
439792b654fSMatthew G. Knepley
440792b654fSMatthew G. Knepley Output Parameters:
44120f4b53cSBarry Smith + is - The global indices for the subproblem, or `NULL`
4420b3275a6SBarry Smith - superdm - The `DM` for the superproblem, which must already have be cloned and contain the same topology as the `dms`
443792b654fSMatthew G. Knepley
444792b654fSMatthew G. Knepley Level: intermediate
445792b654fSMatthew G. Knepley
446a4e35b19SJacob Faibussowitsch .seealso: `DMCreateSuperDM()`, `DMGetLocalSection()`, `DMPlexSetMigrationSF()`, `DMView()`
447792b654fSMatthew G. Knepley @*/
DMCreateSectionSuperDM(DM dms[],PetscInt len,IS * is[],DM * superdm)4485d83a8b1SBarry Smith PetscErrorCode DMCreateSectionSuperDM(DM dms[], PetscInt len, IS *is[], DM *superdm)
449d71ae5a4SJacob Faibussowitsch {
4509796d432SMatthew G. Knepley MPI_Comm comm;
4519796d432SMatthew G. Knepley PetscSection supersection, *sections, *sectionGlobals;
4528cda7954SMatthew G. Knepley PetscInt *Nfs, Nf = 0, f, supf, oldf = -1, nullf = -1, i;
4539796d432SMatthew G. Knepley PetscBool haveNull = PETSC_FALSE;
4542adcc780SMatthew G. Knepley
4552adcc780SMatthew G. Knepley PetscFunctionBegin;
4569566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dms[0], &comm));
4579796d432SMatthew G. Knepley /* Pull out local and global sections */
4589566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(len, &Nfs, len, §ions, len, §ionGlobals));
4592adcc780SMatthew G. Knepley for (i = 0; i < len; ++i) {
4609566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dms[i], §ions[i]));
4619566063dSJacob Faibussowitsch PetscCall(DMGetGlobalSection(dms[i], §ionGlobals[i]));
4627a8be351SBarry Smith PetscCheck(sections[i], comm, PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields");
4637a8be351SBarry Smith PetscCheck(sectionGlobals[i], comm, PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields");
4649566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(sections[i], &Nfs[i]));
4652adcc780SMatthew G. Knepley Nf += Nfs[i];
4662adcc780SMatthew G. Knepley }
4679796d432SMatthew G. Knepley /* Create the supersection */
4689566063dSJacob Faibussowitsch PetscCall(PetscSectionCreateSupersection(sections, len, &supersection));
4699566063dSJacob Faibussowitsch PetscCall(DMSetLocalSection(*superdm, supersection));
4709796d432SMatthew G. Knepley /* Create ISes */
4712adcc780SMatthew G. Knepley if (is) {
4729796d432SMatthew G. Knepley PetscSection supersectionGlobal;
4739796d432SMatthew G. Knepley PetscInt bs = -1, startf = 0;
4742adcc780SMatthew G. Knepley
4759566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(len, is));
4769566063dSJacob Faibussowitsch PetscCall(DMGetGlobalSection(*superdm, &supersectionGlobal));
4779796d432SMatthew G. Knepley for (i = 0; i < len; startf += Nfs[i], ++i) {
4789796d432SMatthew G. Knepley PetscInt *subIndices;
479ec4c761aSMatthew G. Knepley PetscInt subSize, subOff, pStart, pEnd, p, start, end, dummy;
4802adcc780SMatthew G. Knepley
4819566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(sectionGlobals[i], &pStart, &pEnd));
4829566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstrainedStorageSize(sectionGlobals[i], &subSize));
4839566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(subSize, &subIndices));
4849796d432SMatthew G. Knepley for (p = pStart, subOff = 0; p < pEnd; ++p) {
4859796d432SMatthew G. Knepley PetscInt gdof, gcdof, gtdof, d;
4862adcc780SMatthew G. Knepley
4879566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(sectionGlobals[i], p, &gdof));
4889566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(sections[i], p, &gcdof));
4892adcc780SMatthew G. Knepley gtdof = gdof - gcdof;
4902adcc780SMatthew G. Knepley if (gdof > 0 && gtdof) {
4919371c9d4SSatish Balay if (bs < 0) {
4929371c9d4SSatish Balay bs = gtdof;
4939371c9d4SSatish Balay } else if (bs != gtdof) {
4949371c9d4SSatish Balay bs = 1;
4959371c9d4SSatish Balay }
4969566063dSJacob Faibussowitsch PetscCall(DMGetGlobalFieldOffset_Private(*superdm, p, startf, &start, &dummy));
4979566063dSJacob Faibussowitsch PetscCall(DMGetGlobalFieldOffset_Private(*superdm, p, startf + Nfs[i] - 1, &dummy, &end));
49863a3b9bcSJacob 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);
4999796d432SMatthew G. Knepley for (d = start; d < end; ++d, ++subOff) subIndices[subOff] = d;
5002adcc780SMatthew G. Knepley }
5012adcc780SMatthew G. Knepley }
5029566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, subSize, subIndices, PETSC_OWN_POINTER, &(*is)[i]));
5032adcc780SMatthew G. Knepley /* Must have same blocksize on all procs (some might have no points) */
5049796d432SMatthew G. Knepley {
5059796d432SMatthew G. Knepley PetscInt bs = -1, bsLocal[2], bsMinMax[2];
5069796d432SMatthew G. Knepley
5071690c2aeSBarry Smith bsLocal[0] = bs < 0 ? PETSC_INT_MAX : bs;
5089371c9d4SSatish Balay bsLocal[1] = bs;
5099566063dSJacob Faibussowitsch PetscCall(PetscGlobalMinMaxInt(comm, bsLocal, bsMinMax));
5109371c9d4SSatish Balay if (bsMinMax[0] != bsMinMax[1]) {
5119371c9d4SSatish Balay bs = 1;
5129371c9d4SSatish Balay } else {
5139371c9d4SSatish Balay bs = bsMinMax[0];
5149371c9d4SSatish Balay }
5159566063dSJacob Faibussowitsch PetscCall(ISSetBlockSize((*is)[i], bs));
5162adcc780SMatthew G. Knepley }
5172adcc780SMatthew G. Knepley }
5182adcc780SMatthew G. Knepley }
5199796d432SMatthew G. Knepley /* Preserve discretizations */
520e5e52638SMatthew G. Knepley if (len && dms[0]->probs) {
5219566063dSJacob Faibussowitsch PetscCall(DMSetNumFields(*superdm, Nf));
5229796d432SMatthew G. Knepley for (i = 0, supf = 0; i < len; ++i) {
5239796d432SMatthew G. Knepley for (f = 0; f < Nfs[i]; ++f, ++supf) {
5242adcc780SMatthew G. Knepley PetscObject disc;
5252adcc780SMatthew G. Knepley
5269566063dSJacob Faibussowitsch PetscCall(DMGetField(dms[i], f, NULL, &disc));
5279566063dSJacob Faibussowitsch PetscCall(DMSetField(*superdm, supf, NULL, disc));
5282adcc780SMatthew G. Knepley }
5292adcc780SMatthew G. Knepley }
5309566063dSJacob Faibussowitsch PetscCall(DMCreateDS(*superdm));
5312adcc780SMatthew G. Knepley }
532*4758e3ceSMatthew G. Knepley // Create nullspace constructor slots
533*4758e3ceSMatthew G. Knepley PetscCall(PetscFree2((*superdm)->nullspaceConstructors, (*superdm)->nearnullspaceConstructors));
534*4758e3ceSMatthew G. Knepley PetscCall(PetscCalloc2(Nf, &(*superdm)->nullspaceConstructors, Nf, &(*superdm)->nearnullspaceConstructors));
5359796d432SMatthew G. Knepley /* Preserve nullspaces */
5369796d432SMatthew G. Knepley for (i = 0, supf = 0; i < len; ++i) {
5379796d432SMatthew G. Knepley for (f = 0; f < Nfs[i]; ++f, ++supf) {
538*4758e3ceSMatthew G. Knepley if (dms[i]->nullspaceConstructors) {
5399796d432SMatthew G. Knepley (*superdm)->nullspaceConstructors[supf] = dms[i]->nullspaceConstructors[f];
5409796d432SMatthew G. Knepley if ((*superdm)->nullspaceConstructors[supf]) {
5419796d432SMatthew G. Knepley haveNull = PETSC_TRUE;
5429796d432SMatthew G. Knepley nullf = supf;
5438cda7954SMatthew G. Knepley oldf = f;
5442adcc780SMatthew G. Knepley }
5459796d432SMatthew G. Knepley }
5469796d432SMatthew G. Knepley }
547*4758e3ceSMatthew G. Knepley }
5489796d432SMatthew G. Knepley /* Attach nullspace to IS */
5499796d432SMatthew G. Knepley if (haveNull && is) {
5509796d432SMatthew G. Knepley MatNullSpace nullSpace;
5519796d432SMatthew G. Knepley
5529566063dSJacob Faibussowitsch PetscCall((*(*superdm)->nullspaceConstructors[nullf])(*superdm, oldf, nullf, &nullSpace));
5539566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)(*is)[nullf], "nullspace", (PetscObject)nullSpace));
5549566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&nullSpace));
5559796d432SMatthew G. Knepley }
5569566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&supersection));
5579566063dSJacob Faibussowitsch PetscCall(PetscFree3(Nfs, sections, sectionGlobals));
5583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
5592adcc780SMatthew G. Knepley }
560