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