Actual source code: dmi.c
  1: #include <petsc/private/dmimpl.h>
  2: #include <petscds.h>
  4: // Greatest common divisor of two nonnegative integers
  5: PetscInt PetscGCD(PetscInt a, PetscInt b)
  6: {
  7:   while (b != 0) {
  8:     PetscInt tmp = b;
  9:     b            = a % b;
 10:     a            = tmp;
 11:   }
 12:   return a;
 13: }
 15: PetscErrorCode DMCreateGlobalVector_Section_Private(DM dm, Vec *vec)
 16: {
 17:   PetscSection gSection;
 18:   PetscInt     localSize, bs, blockSize = -1, pStart, pEnd, p;
 19:   PetscInt     in[2], out[2];
 21:   PetscFunctionBegin;
 22:   PetscCall(DMGetGlobalSection(dm, &gSection));
 23:   PetscCall(PetscSectionGetChart(gSection, &pStart, &pEnd));
 24:   for (p = pStart; p < pEnd; ++p) {
 25:     PetscInt dof, cdof;
 27:     PetscCall(PetscSectionGetDof(gSection, p, &dof));
 28:     PetscCall(PetscSectionGetConstraintDof(gSection, p, &cdof));
 30:     if (dof - cdof > 0) {
 31:       if (blockSize < 0) {
 32:         /* set blockSize */
 33:         blockSize = dof - cdof;
 34:       } else {
 35:         blockSize = PetscGCD(dof - cdof, blockSize);
 36:       }
 37:     }
 38:   }
 40:   // You cannot negate PETSC_INT_MIN
 41:   in[0] = blockSize < 0 ? -PETSC_INT_MAX : -blockSize;
 42:   in[1] = blockSize;
 43:   PetscCallMPI(MPIU_Allreduce(in, out, 2, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm)));
 44:   /* -out[0] = min(blockSize), out[1] = max(blockSize) */
 45:   if (-out[0] == out[1]) {
 46:     bs = out[1];
 47:   } else bs = 1;
 49:   if (bs < 0) { /* Everyone was empty */
 50:     blockSize = 1;
 51:     bs        = 1;
 52:   }
 54:   PetscCall(PetscSectionGetConstrainedStorageSize(gSection, &localSize));
 55:   PetscCheck(localSize % blockSize == 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mismatch between blocksize %" PetscInt_FMT " and local storage size %" PetscInt_FMT, blockSize, localSize);
 56:   PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), vec));
 57:   PetscCall(VecSetSizes(*vec, localSize, PETSC_DETERMINE));
 58:   PetscCall(VecSetBlockSize(*vec, bs));
 59:   PetscCall(VecSetType(*vec, dm->vectype));
 60:   PetscCall(VecSetDM(*vec, dm));
 61:   /* PetscCall(VecSetLocalToGlobalMapping(*vec, dm->ltogmap)); */
 62:   PetscFunctionReturn(PETSC_SUCCESS);
 63: }
 65: PetscErrorCode DMCreateLocalVector_Section_Private(DM dm, Vec *vec)
 66: {
 67:   PetscSection section;
 68:   PetscInt     localSize, blockSize = -1, pStart, pEnd, p;
 70:   PetscFunctionBegin;
 71:   PetscCall(DMGetLocalSection(dm, §ion));
 72:   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
 73:   for (p = pStart; p < pEnd; ++p) {
 74:     PetscInt dof;
 76:     PetscCall(PetscSectionGetDof(section, p, &dof));
 77:     if ((blockSize < 0) && (dof > 0)) blockSize = dof;
 78:     if (dof > 0) blockSize = PetscGCD(dof, blockSize);
 79:   }
 80:   PetscCall(PetscSectionGetStorageSize(section, &localSize));
 81:   PetscCall(VecCreate(PETSC_COMM_SELF, vec));
 82:   PetscCall(VecSetSizes(*vec, localSize, localSize));
 83:   PetscCall(VecSetBlockSize(*vec, PetscAbs(blockSize)));
 84:   PetscCall(VecSetType(*vec, dm->vectype));
 85:   PetscCall(VecSetDM(*vec, dm));
 86:   PetscFunctionReturn(PETSC_SUCCESS);
 87: }
 89: static PetscErrorCode PetscSectionSelectFields_Private(PetscSection s, PetscSection gs, PetscInt numFields, const PetscInt fields[], const PetscInt numComps[], const PetscInt comps[], IS *is)
 90: {
 91:   PetscInt *subIndices;
 92:   PetscInt  mbs, bs = 0, bsLocal[2], bsMinMax[2];
 93:   PetscInt  pStart, pEnd, Nc, subSize = 0, subOff = 0;
 95:   PetscFunctionBegin;
 96:   PetscCall(PetscSectionGetChart(gs, &pStart, &pEnd));
 97:   if (numComps) {
 98:     for (PetscInt f = 0, off = 0; f < numFields; ++f) {
 99:       PetscInt Nc;
101:       if (numComps[f] < 0) continue;
102:       PetscCall(PetscSectionGetFieldComponents(s, f, &Nc));
103:       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);
104:       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);
105:       bs += numComps[f];
106:     }
107:   } else {
108:     for (PetscInt f = 0; f < numFields; ++f) {
109:       PetscInt Nc;
111:       PetscCall(PetscSectionGetFieldComponents(s, fields[f], &Nc));
112:       bs += Nc;
113:     }
114:   }
115:   mbs = -1; /* multiple of block size not set */
116:   for (PetscInt p = pStart; p < pEnd; ++p) {
117:     PetscInt gdof, pSubSize = 0;
119:     PetscCall(PetscSectionGetDof(gs, p, &gdof));
120:     if (gdof > 0) {
121:       PetscInt off = 0;
123:       for (PetscInt f = 0; f < numFields; ++f) {
124:         PetscInt fdof, fcdof, sfdof, sfcdof = 0;
126:         PetscCall(PetscSectionGetFieldComponents(s, f, &Nc));
127:         PetscCall(PetscSectionGetFieldDof(s, p, fields[f], &fdof));
128:         PetscCall(PetscSectionGetFieldConstraintDof(s, p, fields[f], &fcdof));
129:         if (numComps && numComps[f] >= 0) {
130:           const PetscInt *ind;
132:           // Assume sets of dofs on points are of size Nc
133:           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);
134:           sfdof = (fdof / Nc) * numComps[f];
135:           PetscCall(PetscSectionGetFieldConstraintIndices(s, p, fields[f], &ind));
136:           for (PetscInt i = 0; i < (fdof / Nc); ++i) {
137:             for (PetscInt c = 0, fcc = 0; c < Nc; ++c) {
138:               if (c == comps[off + fcc]) {
139:                 ++fcc;
140:                 ++sfcdof;
141:               }
142:             }
143:           }
144:           pSubSize += sfdof - sfcdof;
145:           off += numComps[f];
146:         } else {
147:           pSubSize += fdof - fcdof;
148:         }
149:       }
150:       subSize += pSubSize;
151:       if (pSubSize && pSubSize % bs) {
152:         // Layout does not admit a pointwise block size -> set mbs to 0
153:         mbs = 0;
154:       } else if (pSubSize) {
155:         if (mbs == -1) mbs = pSubSize / bs;
156:         else mbs = PetscMin(mbs, pSubSize / bs);
157:       }
158:     }
159:   }
161:   // Must have same blocksize on all procs (some might have no points)
162:   bsLocal[0] = mbs < 0 ? PETSC_INT_MAX : mbs;
163:   bsLocal[1] = mbs;
164:   PetscCall(PetscGlobalMinMaxInt(PetscObjectComm((PetscObject)gs), bsLocal, bsMinMax));
165:   if (bsMinMax[0] != bsMinMax[1]) { /* different multiple of block size -> set bs to 1 */
166:     bs = 1;
167:   } else { /* same multiple */
168:     mbs = bsMinMax[0];
169:     bs *= mbs;
170:   }
171:   PetscCall(PetscMalloc1(subSize, &subIndices));
172:   for (PetscInt p = pStart; p < pEnd; ++p) {
173:     PetscInt gdof, goff;
175:     PetscCall(PetscSectionGetDof(gs, p, &gdof));
176:     if (gdof > 0) {
177:       PetscInt off = 0;
179:       PetscCall(PetscSectionGetOffset(gs, p, &goff));
180:       for (PetscInt f = 0; f < numFields; ++f) {
181:         PetscInt fdof, fcdof, poff = 0;
183:         /* Can get rid of this loop by storing field information in the global section */
184:         for (PetscInt f2 = 0; f2 < fields[f]; ++f2) {
185:           PetscCall(PetscSectionGetFieldDof(s, p, f2, &fdof));
186:           PetscCall(PetscSectionGetFieldConstraintDof(s, p, f2, &fcdof));
187:           poff += fdof - fcdof;
188:         }
189:         PetscCall(PetscSectionGetFieldDof(s, p, fields[f], &fdof));
190:         PetscCall(PetscSectionGetFieldConstraintDof(s, p, fields[f], &fcdof));
192:         if (numComps && numComps[f] >= 0) {
193:           const PetscInt *ind;
195:           // Assume sets of dofs on points are of size Nc
196:           PetscCall(PetscSectionGetFieldConstraintIndices(s, p, fields[f], &ind));
197:           for (PetscInt i = 0, fcoff = 0, pfoff = 0; i < (fdof / Nc); ++i) {
198:             for (PetscInt c = 0, fcc = 0; c < Nc; ++c) {
199:               const PetscInt k = i * Nc + c;
201:               if (ind[fcoff] == k) {
202:                 ++fcoff;
203:                 continue;
204:               }
205:               if (c == comps[off + fcc]) {
206:                 ++fcc;
207:                 subIndices[subOff++] = goff + poff + pfoff;
208:               }
209:               ++pfoff;
210:             }
211:           }
212:           off += numComps[f];
213:         } else {
214:           for (PetscInt fc = 0; fc < fdof - fcdof; ++fc, ++subOff) subIndices[subOff] = goff + poff + fc;
215:         }
216:       }
217:     }
218:   }
219:   PetscCheck(subSize == subOff, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "The offset array size %" PetscInt_FMT " != %" PetscInt_FMT " the number of indices", subSize, subOff);
220:   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)gs), subSize, subIndices, PETSC_OWN_POINTER, is));
221:   if (bs > 1) {
222:     // We need to check that the block size does not come from non-contiguous fields
223:     PetscInt set = 1, rset = 1;
224:     for (PetscInt i = 0; i < subSize; i += bs) {
225:       for (PetscInt j = 0; j < bs; ++j) {
226:         if (subIndices[i + j] != subIndices[i] + j) {
227:           set = 0;
228:           break;
229:         }
230:       }
231:     }
232:     PetscCallMPI(MPIU_Allreduce(&set, &rset, 1, MPIU_INT, MPI_PROD, PetscObjectComm((PetscObject)gs)));
233:     if (rset) PetscCall(ISSetBlockSize(*is, bs));
234:   }
235:   PetscFunctionReturn(PETSC_SUCCESS);
236: }
238: static PetscErrorCode DMSelectFields_Private(DM dm, PetscSection section, PetscInt numFields, const PetscInt fields[], const PetscInt numComps[], const PetscInt comps[], IS *is, DM *subdm)
239: {
240:   PetscSection subsection;
241:   PetscBool    haveNull = PETSC_FALSE;
242:   PetscInt     nf = 0, of = 0;
244:   PetscFunctionBegin;
245:   if (numComps) {
246:     const PetscInt field = fields[0];
248:     PetscCheck(numFields == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "We only support a single field for component selection right now");
249:     PetscCall(PetscSectionCreateComponentSubsection(section, numComps[field], comps, &subsection));
250:     PetscCall(DMSetLocalSection(*subdm, subsection));
251:     PetscCall(PetscSectionDestroy(&subsection));
252:     (*subdm)->nullspaceConstructors[field] = dm->nullspaceConstructors[field];
253:     if (dm->probs) {
254:       PetscFV  fv, fvNew;
255:       PetscInt fnum[1] = {field};
257:       PetscCall(DMSetNumFields(*subdm, 1));
258:       PetscCall(DMGetField(dm, field, NULL, (PetscObject *)&fv));
259:       PetscCall(PetscFVClone(fv, &fvNew));
260:       PetscCall(PetscFVSetNumComponents(fvNew, numComps[0]));
261:       PetscCall(DMSetField(*subdm, 0, NULL, (PetscObject)fvNew));
262:       PetscCall(PetscFVDestroy(&fvNew));
263:       PetscCall(DMCreateDS(*subdm));
264:       if (numComps[0] == 1 && is) {
265:         PetscObject disc, space, pmat;
267:         PetscCall(DMGetField(*subdm, field, NULL, &disc));
268:         PetscCall(PetscObjectQuery(disc, "nullspace", &space));
269:         if (space) PetscCall(PetscObjectCompose((PetscObject)*is, "nullspace", space));
270:         PetscCall(PetscObjectQuery(disc, "nearnullspace", &space));
271:         if (space) PetscCall(PetscObjectCompose((PetscObject)*is, "nearnullspace", space));
272:         PetscCall(PetscObjectQuery(disc, "pmat", &pmat));
273:         if (pmat) PetscCall(PetscObjectCompose((PetscObject)*is, "pmat", pmat));
274:       }
275:       PetscCall(PetscDSCopyConstants(dm->probs[field].ds, (*subdm)->probs[0].ds));
276:       PetscCall(PetscDSCopyBoundary(dm->probs[field].ds, 1, fnum, (*subdm)->probs[0].ds));
277:       PetscCall(PetscDSSelectEquations(dm->probs[field].ds, 1, fnum, (*subdm)->probs[0].ds));
278:     }
279:     if ((*subdm)->nullspaceConstructors[0] && is) {
280:       MatNullSpace nullSpace;
282:       PetscCall((*(*subdm)->nullspaceConstructors[0])(*subdm, 0, 0, &nullSpace));
283:       PetscCall(PetscObjectCompose((PetscObject)*is, "nullspace", (PetscObject)nullSpace));
284:       PetscCall(MatNullSpaceDestroy(&nullSpace));
285:     }
286:     if (dm->coarseMesh) PetscCall(DMCreateSubDM(dm->coarseMesh, numFields, fields, NULL, &(*subdm)->coarseMesh));
287:     PetscFunctionReturn(PETSC_SUCCESS);
288:   }
289:   PetscCall(PetscSectionCreateSubsection(section, numFields, fields, &subsection));
290:   PetscCall(DMSetLocalSection(*subdm, subsection));
291:   PetscCall(PetscSectionDestroy(&subsection));
292:   for (PetscInt f = 0; f < numFields; ++f) {
293:     (*subdm)->nullspaceConstructors[f] = dm->nullspaceConstructors[fields[f]];
294:     if ((*subdm)->nullspaceConstructors[f]) {
295:       haveNull = PETSC_TRUE;
296:       nf       = f;
297:       of       = fields[f];
298:     }
299:   }
300:   if (dm->probs) {
301:     PetscCall(DMSetNumFields(*subdm, numFields));
302:     for (PetscInt f = 0; f < numFields; ++f) {
303:       PetscObject disc;
305:       PetscCall(DMGetField(dm, fields[f], NULL, &disc));
306:       PetscCall(DMSetField(*subdm, f, NULL, disc));
307:     }
308:     // TODO: if only FV, then cut down the components
309:     PetscCall(DMCreateDS(*subdm));
310:     if (numFields == 1 && is) {
311:       PetscObject disc, space, pmat;
313:       PetscCall(DMGetField(*subdm, 0, NULL, &disc));
314:       PetscCall(PetscObjectQuery(disc, "nullspace", &space));
315:       if (space) PetscCall(PetscObjectCompose((PetscObject)*is, "nullspace", space));
316:       PetscCall(PetscObjectQuery(disc, "nearnullspace", &space));
317:       if (space) PetscCall(PetscObjectCompose((PetscObject)*is, "nearnullspace", space));
318:       PetscCall(PetscObjectQuery(disc, "pmat", &pmat));
319:       if (pmat) PetscCall(PetscObjectCompose((PetscObject)*is, "pmat", pmat));
320:     }
321:     // Check if DSes record their DM fields
322:     if (dm->probs[0].fields) {
323:       PetscInt d, e;
325:       for (d = 0, e = 0; d < dm->Nds && e < (*subdm)->Nds; ++d) {
326:         const PetscInt  Nf = dm->probs[d].ds->Nf;
327:         const PetscInt *fld;
328:         PetscInt        f, g;
330:         PetscCall(ISGetIndices(dm->probs[d].fields, &fld));
331:         for (f = 0; f < Nf; ++f) {
332:           for (g = 0; g < numFields; ++g)
333:             if (fld[f] == fields[g]) break;
334:           if (g < numFields) break;
335:         }
336:         PetscCall(ISRestoreIndices(dm->probs[d].fields, &fld));
337:         if (f == Nf) continue;
338:         PetscCall(PetscDSCopyConstants(dm->probs[d].ds, (*subdm)->probs[e].ds));
339:         PetscCall(PetscDSCopyBoundary(dm->probs[d].ds, numFields, fields, (*subdm)->probs[e].ds));
340:         // Translate DM fields to DS fields
341:         {
342:           IS              infields, dsfields;
343:           const PetscInt *fld, *ofld;
344:           PetscInt       *fidx;
345:           PetscInt        onf, nf;
347:           PetscCall(ISCreateGeneral(PETSC_COMM_SELF, numFields, fields, PETSC_USE_POINTER, &infields));
348:           PetscCall(ISIntersect(infields, dm->probs[d].fields, &dsfields));
349:           PetscCall(ISDestroy(&infields));
350:           PetscCall(ISGetLocalSize(dsfields, &nf));
351:           PetscCheck(nf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "DS cannot be supported on 0 fields");
352:           PetscCall(ISGetIndices(dsfields, &fld));
353:           PetscCall(ISGetLocalSize(dm->probs[d].fields, &onf));
354:           PetscCall(ISGetIndices(dm->probs[d].fields, &ofld));
355:           PetscCall(PetscMalloc1(nf, &fidx));
356:           for (PetscInt f = 0, g = 0; f < onf && g < nf; ++f) {
357:             if (ofld[f] == fld[g]) fidx[g++] = f;
358:           }
359:           PetscCall(ISRestoreIndices(dm->probs[d].fields, &ofld));
360:           PetscCall(ISRestoreIndices(dsfields, &fld));
361:           PetscCall(ISDestroy(&dsfields));
362:           PetscCall(PetscDSSelectDiscretizations(dm->probs[0].ds, nf, fidx, PETSC_DETERMINE, PETSC_DETERMINE, (*subdm)->probs[0].ds));
363:           PetscCall(PetscDSSelectEquations(dm->probs[0].ds, nf, fidx, (*subdm)->probs[0].ds));
364:           PetscCall(PetscFree(fidx));
365:         }
366:         ++e;
367:       }
368:     } else {
369:       PetscCall(PetscDSCopyConstants(dm->probs[0].ds, (*subdm)->probs[0].ds));
370:       PetscCall(PetscDSCopyBoundary(dm->probs[0].ds, PETSC_DETERMINE, NULL, (*subdm)->probs[0].ds));
371:       PetscCall(PetscDSSelectDiscretizations(dm->probs[0].ds, numFields, fields, PETSC_DETERMINE, PETSC_DETERMINE, (*subdm)->probs[0].ds));
372:       PetscCall(PetscDSSelectEquations(dm->probs[0].ds, numFields, fields, (*subdm)->probs[0].ds));
373:     }
374:   }
375:   if (haveNull && is) {
376:     MatNullSpace nullSpace;
378:     PetscCall((*(*subdm)->nullspaceConstructors[nf])(*subdm, of, nf, &nullSpace));
379:     PetscCall(PetscObjectCompose((PetscObject)*is, "nullspace", (PetscObject)nullSpace));
380:     PetscCall(MatNullSpaceDestroy(&nullSpace));
381:   }
382:   if (dm->coarseMesh) PetscCall(DMCreateSubDM(dm->coarseMesh, numFields, fields, NULL, &(*subdm)->coarseMesh));
383:   PetscFunctionReturn(PETSC_SUCCESS);
384: }
386: /*@
387:   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`.
389:   Not Collective
391:   Input Parameters:
392: + dm        - The `DM` object
393: . numFields - The number of fields to incorporate into `subdm`
394: . fields    - The field numbers of the selected fields
395: . numComps  - The number of components from each field to incorporate into `subdm`, or PETSC_DECIDE for all components
396: - comps     - The component numbers of the selected fields (omitted for PTESC_DECIDE fields)
398:   Output Parameters:
399: + is    - The global indices for the subproblem or `NULL`
400: - subdm - The `DM` for the subproblem, which must already have be cloned from `dm` or `NULL`
402:   Level: intermediate
404:   Notes:
405:   If `is` and `subdm` are both `NULL` this does nothing
407: .seealso: `DMCreateSubDM()`, `DMGetLocalSection()`, `DMPlexSetMigrationSF()`, `DMView()`
408: @*/
409: PetscErrorCode DMCreateSectionSubDM(DM dm, PetscInt numFields, const PetscInt fields[], const PetscInt numComps[], const PetscInt comps[], IS *is, DM *subdm)
410: {
411:   PetscSection section, sectionGlobal;
412:   PetscInt     Nf;
414:   PetscFunctionBegin;
415:   if (!numFields) PetscFunctionReturn(PETSC_SUCCESS);
416:   PetscCall(DMGetLocalSection(dm, §ion));
417:   PetscCall(DMGetGlobalSection(dm, §ionGlobal));
418:   PetscCheck(section, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields");
419:   PetscCheck(sectionGlobal, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields");
420:   PetscCall(PetscSectionGetNumFields(section, &Nf));
421:   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);
423:   if (is) PetscCall(PetscSectionSelectFields_Private(section, sectionGlobal, numFields, fields, numComps, comps, is));
424:   if (subdm) PetscCall(DMSelectFields_Private(dm, section, numFields, fields, numComps, comps, is, subdm));
425:   PetscFunctionReturn(PETSC_SUCCESS);
426: }
428: /*@C
429:   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`
431:   Not Collective
433:   Input Parameters:
434: + dms - The `DM` objects, the must all have the same topology; for example obtained with `DMClone()`
435: - len - The number of `DM` in `dms`
437:   Output Parameters:
438: + is      - The global indices for the subproblem, or `NULL`
439: - superdm - The `DM` for the superproblem, which must already have be cloned and contain the same topology as the `dms`
441:   Level: intermediate
443: .seealso: `DMCreateSuperDM()`, `DMGetLocalSection()`, `DMPlexSetMigrationSF()`, `DMView()`
444: @*/
445: PetscErrorCode DMCreateSectionSuperDM(DM dms[], PetscInt len, IS *is[], DM *superdm)
446: {
447:   MPI_Comm     comm;
448:   PetscSection supersection, *sections, *sectionGlobals;
449:   PetscInt    *Nfs, Nf = 0, f, supf, oldf = -1, nullf = -1, i;
450:   PetscBool    haveNull = PETSC_FALSE;
452:   PetscFunctionBegin;
453:   PetscCall(PetscObjectGetComm((PetscObject)dms[0], &comm));
454:   /* Pull out local and global sections */
455:   PetscCall(PetscMalloc3(len, &Nfs, len, §ions, len, §ionGlobals));
456:   for (i = 0; i < len; ++i) {
457:     PetscCall(DMGetLocalSection(dms[i], §ions[i]));
458:     PetscCall(DMGetGlobalSection(dms[i], §ionGlobals[i]));
459:     PetscCheck(sections[i], comm, PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields");
460:     PetscCheck(sectionGlobals[i], comm, PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields");
461:     PetscCall(PetscSectionGetNumFields(sections[i], &Nfs[i]));
462:     Nf += Nfs[i];
463:   }
464:   /* Create the supersection */
465:   PetscCall(PetscSectionCreateSupersection(sections, len, &supersection));
466:   PetscCall(DMSetLocalSection(*superdm, supersection));
467:   /* Create ISes */
468:   if (is) {
469:     PetscSection supersectionGlobal;
470:     PetscInt     bs = -1, startf = 0;
472:     PetscCall(PetscMalloc1(len, is));
473:     PetscCall(DMGetGlobalSection(*superdm, &supersectionGlobal));
474:     for (i = 0; i < len; startf += Nfs[i], ++i) {
475:       PetscInt *subIndices;
476:       PetscInt  subSize, subOff, pStart, pEnd, p, start, end, dummy;
478:       PetscCall(PetscSectionGetChart(sectionGlobals[i], &pStart, &pEnd));
479:       PetscCall(PetscSectionGetConstrainedStorageSize(sectionGlobals[i], &subSize));
480:       PetscCall(PetscMalloc1(subSize, &subIndices));
481:       for (p = pStart, subOff = 0; p < pEnd; ++p) {
482:         PetscInt gdof, gcdof, gtdof, d;
484:         PetscCall(PetscSectionGetDof(sectionGlobals[i], p, &gdof));
485:         PetscCall(PetscSectionGetConstraintDof(sections[i], p, &gcdof));
486:         gtdof = gdof - gcdof;
487:         if (gdof > 0 && gtdof) {
488:           if (bs < 0) {
489:             bs = gtdof;
490:           } else if (bs != gtdof) {
491:             bs = 1;
492:           }
493:           PetscCall(DMGetGlobalFieldOffset_Private(*superdm, p, startf, &start, &dummy));
494:           PetscCall(DMGetGlobalFieldOffset_Private(*superdm, p, startf + Nfs[i] - 1, &dummy, &end));
495:           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);
496:           for (d = start; d < end; ++d, ++subOff) subIndices[subOff] = d;
497:         }
498:       }
499:       PetscCall(ISCreateGeneral(comm, subSize, subIndices, PETSC_OWN_POINTER, &(*is)[i]));
500:       /* Must have same blocksize on all procs (some might have no points) */
501:       {
502:         PetscInt bs = -1, bsLocal[2], bsMinMax[2];
504:         bsLocal[0] = bs < 0 ? PETSC_INT_MAX : bs;
505:         bsLocal[1] = bs;
506:         PetscCall(PetscGlobalMinMaxInt(comm, bsLocal, bsMinMax));
507:         if (bsMinMax[0] != bsMinMax[1]) {
508:           bs = 1;
509:         } else {
510:           bs = bsMinMax[0];
511:         }
512:         PetscCall(ISSetBlockSize((*is)[i], bs));
513:       }
514:     }
515:   }
516:   /* Preserve discretizations */
517:   if (len && dms[0]->probs) {
518:     PetscCall(DMSetNumFields(*superdm, Nf));
519:     for (i = 0, supf = 0; i < len; ++i) {
520:       for (f = 0; f < Nfs[i]; ++f, ++supf) {
521:         PetscObject disc;
523:         PetscCall(DMGetField(dms[i], f, NULL, &disc));
524:         PetscCall(DMSetField(*superdm, supf, NULL, disc));
525:       }
526:     }
527:     PetscCall(DMCreateDS(*superdm));
528:   }
529:   /* Preserve nullspaces */
530:   for (i = 0, supf = 0; i < len; ++i) {
531:     for (f = 0; f < Nfs[i]; ++f, ++supf) {
532:       (*superdm)->nullspaceConstructors[supf] = dms[i]->nullspaceConstructors[f];
533:       if ((*superdm)->nullspaceConstructors[supf]) {
534:         haveNull = PETSC_TRUE;
535:         nullf    = supf;
536:         oldf     = f;
537:       }
538:     }
539:   }
540:   /* Attach nullspace to IS */
541:   if (haveNull && is) {
542:     MatNullSpace nullSpace;
544:     PetscCall((*(*superdm)->nullspaceConstructors[nullf])(*superdm, oldf, nullf, &nullSpace));
545:     PetscCall(PetscObjectCompose((PetscObject)(*is)[nullf], "nullspace", (PetscObject)nullSpace));
546:     PetscCall(MatNullSpaceDestroy(&nullSpace));
547:   }
548:   PetscCall(PetscSectionDestroy(&supersection));
549:   PetscCall(PetscFree3(Nfs, sections, sectionGlobals));
550:   PetscFunctionReturn(PETSC_SUCCESS);
551: }