xref: /petsc/src/dm/interface/dmi.c (revision e87a4003e3e772bdbe5118ac0e5c086db060906d)
1af0996ceSBarry Smith #include <petsc/private/dmimpl.h>     /*I      "petscdm.h"     I*/
22764a2aaSMatthew G. Knepley #include <petscds.h>
311689aebSJed Brown 
411689aebSJed Brown PetscErrorCode DMCreateGlobalVector_Section_Private(DM dm,Vec *vec)
511689aebSJed Brown {
611689aebSJed Brown   PetscSection   gSection;
7cc30b04dSMatthew G Knepley   PetscInt       localSize, bs, blockSize = -1, pStart, pEnd, p;
8fad22124SMatthew G Knepley   PetscErrorCode ierr;
95d1c0279SHong Zhang   PetscInt       in[2],out[2];
1011689aebSJed Brown 
1111689aebSJed Brown   PetscFunctionBegin;
12*e87a4003SBarry Smith   ierr = DMGetGlobalSection(dm, &gSection);CHKERRQ(ierr);
13fad22124SMatthew G Knepley   ierr = PetscSectionGetChart(gSection, &pStart, &pEnd);CHKERRQ(ierr);
1411689aebSJed Brown   for (p = pStart; p < pEnd; ++p) {
1511689aebSJed Brown     PetscInt dof, cdof;
1611689aebSJed Brown 
17fad22124SMatthew G Knepley     ierr = PetscSectionGetDof(gSection, p, &dof);CHKERRQ(ierr);
18fad22124SMatthew G Knepley     ierr = PetscSectionGetConstraintDof(gSection, p, &cdof);CHKERRQ(ierr);
190680ce0aSHong Zhang 
200680ce0aSHong Zhang     if (dof > 0) {
210680ce0aSHong Zhang       if (blockSize < 0 && dof-cdof > 0) {
220680ce0aSHong Zhang         /* set blockSize */
230680ce0aSHong Zhang         blockSize = dof-cdof;
240680ce0aSHong Zhang       } else if (dof-cdof != blockSize) {
250680ce0aSHong Zhang         /* non-identical blockSize, set it as 1 */
2611689aebSJed Brown         blockSize = 1;
2711689aebSJed Brown         break;
2811689aebSJed Brown       }
2911689aebSJed Brown     }
300680ce0aSHong Zhang   }
315d1c0279SHong Zhang 
322c8be454SMatthew G. Knepley   in[0] = blockSize < 0 ? PETSC_MIN_INT : -blockSize;
335d1c0279SHong Zhang   in[1] = blockSize;
34bdfec8b8SHong Zhang   ierr = MPIU_Allreduce(in,out,2,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
355d1c0279SHong Zhang   /* -out[0] = min(blockSize), out[1] = max(blockSize) */
365d1c0279SHong Zhang   if (-out[0] == out[1]) {
375d1c0279SHong Zhang     bs = out[1];
385d1c0279SHong Zhang   } else bs = 1;
395d1c0279SHong Zhang 
405d1c0279SHong Zhang   if (bs < 0) { /* Everyone was empty */
415d1c0279SHong Zhang     blockSize = 1;
425d1c0279SHong Zhang     bs        = 1;
435d1c0279SHong Zhang   }
445d1c0279SHong Zhang 
45fad22124SMatthew G Knepley   ierr = PetscSectionGetConstrainedStorageSize(gSection, &localSize);CHKERRQ(ierr);
4682f516ccSBarry Smith   if (localSize%blockSize) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mismatch between blocksize %d and local storage size %d", blockSize, localSize);
4782f516ccSBarry Smith   ierr = VecCreate(PetscObjectComm((PetscObject)dm), vec);CHKERRQ(ierr);
4811689aebSJed Brown   ierr = VecSetSizes(*vec, localSize, PETSC_DETERMINE);CHKERRQ(ierr);
49cc30b04dSMatthew G Knepley   ierr = VecSetBlockSize(*vec, bs);CHKERRQ(ierr);
50c0dedaeaSBarry Smith   ierr = VecSetType(*vec,dm->vectype);CHKERRQ(ierr);
5111689aebSJed Brown   ierr = VecSetDM(*vec, dm);CHKERRQ(ierr);
5211689aebSJed Brown   /* ierr = VecSetLocalToGlobalMapping(*vec, dm->ltogmap);CHKERRQ(ierr); */
5311689aebSJed Brown   PetscFunctionReturn(0);
5411689aebSJed Brown }
5511689aebSJed Brown 
5611689aebSJed Brown PetscErrorCode DMCreateLocalVector_Section_Private(DM dm,Vec *vec)
5711689aebSJed Brown {
58fad22124SMatthew G Knepley   PetscSection   section;
5911689aebSJed Brown   PetscInt       localSize, blockSize = -1, pStart, pEnd, p;
60fad22124SMatthew G Knepley   PetscErrorCode ierr;
6111689aebSJed Brown 
6211689aebSJed Brown   PetscFunctionBegin;
63*e87a4003SBarry Smith   ierr = DMGetSection(dm, &section);CHKERRQ(ierr);
64fad22124SMatthew G Knepley   ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
6511689aebSJed Brown   for (p = pStart; p < pEnd; ++p) {
6611689aebSJed Brown     PetscInt dof;
6711689aebSJed Brown 
68fad22124SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
6911689aebSJed Brown     if ((blockSize < 0) && (dof > 0)) blockSize = dof;
7011689aebSJed Brown     if ((dof > 0) && (dof != blockSize)) {
7111689aebSJed Brown       blockSize = 1;
7211689aebSJed Brown       break;
7311689aebSJed Brown     }
7411689aebSJed Brown   }
75fad22124SMatthew G Knepley   ierr = PetscSectionGetStorageSize(section, &localSize);CHKERRQ(ierr);
7611689aebSJed Brown   ierr = VecCreate(PETSC_COMM_SELF, vec);CHKERRQ(ierr);
7711689aebSJed Brown   ierr = VecSetSizes(*vec, localSize, localSize);CHKERRQ(ierr);
7811689aebSJed Brown   ierr = VecSetBlockSize(*vec, blockSize);CHKERRQ(ierr);
79c0dedaeaSBarry Smith   ierr = VecSetType(*vec,dm->vectype);CHKERRQ(ierr);
8011689aebSJed Brown   ierr = VecSetDM(*vec, dm);CHKERRQ(ierr);
8111689aebSJed Brown   PetscFunctionReturn(0);
8211689aebSJed Brown }
834d9407bcSMatthew G. Knepley 
844d9407bcSMatthew G. Knepley /* This assumes that the DM has been cloned prior to the call */
85276c5506SMatthew G. Knepley PetscErrorCode DMCreateSubDM_Section_Private(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm)
864d9407bcSMatthew G. Knepley {
874d9407bcSMatthew G. Knepley   PetscSection   section, sectionGlobal;
884d9407bcSMatthew G. Knepley   PetscInt      *subIndices;
894d9407bcSMatthew G. Knepley   PetscInt       subSize = 0, subOff = 0, nF, f, pStart, pEnd, p;
904d9407bcSMatthew G. Knepley   PetscErrorCode ierr;
914d9407bcSMatthew G. Knepley 
924d9407bcSMatthew G. Knepley   PetscFunctionBegin;
934d9407bcSMatthew G. Knepley   if (!numFields) PetscFunctionReturn(0);
94*e87a4003SBarry Smith   ierr = DMGetSection(dm, &section);CHKERRQ(ierr);
95*e87a4003SBarry Smith   ierr = DMGetGlobalSection(dm, &sectionGlobal);CHKERRQ(ierr);
964d9407bcSMatthew G. Knepley   if (!section) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields");
974d9407bcSMatthew G. Knepley   if (!sectionGlobal) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields");
984d9407bcSMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &nF);CHKERRQ(ierr);
994d9407bcSMatthew G. Knepley   if (numFields > nF) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Number of requested fields %d greater than number of DM fields %d", numFields, nF);
1004d9407bcSMatthew G. Knepley   if (is) {
1010be3e97aSMatthew G. Knepley     PetscInt bs = -1, bsLocal[2], bsMinMax[2];
1024d9407bcSMatthew G. Knepley     ierr = PetscSectionGetChart(sectionGlobal, &pStart, &pEnd);CHKERRQ(ierr);
1034d9407bcSMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
104b9eca265Ssarens       PetscInt gdof, pSubSize  = 0;
1054d9407bcSMatthew G. Knepley 
1064d9407bcSMatthew G. Knepley       ierr = PetscSectionGetDof(sectionGlobal, p, &gdof);CHKERRQ(ierr);
1074d9407bcSMatthew G. Knepley       if (gdof > 0) {
1084d9407bcSMatthew G. Knepley         for (f = 0; f < numFields; ++f) {
1094d9407bcSMatthew G. Knepley           PetscInt fdof, fcdof;
1104d9407bcSMatthew G. Knepley 
1114d9407bcSMatthew G. Knepley           ierr     = PetscSectionGetFieldDof(section, p, fields[f], &fdof);CHKERRQ(ierr);
1124d9407bcSMatthew G. Knepley           ierr     = PetscSectionGetFieldConstraintDof(section, p, fields[f], &fcdof);CHKERRQ(ierr);
113b9eca265Ssarens           pSubSize += fdof-fcdof;
114b9eca265Ssarens         }
115b9eca265Ssarens         subSize += pSubSize;
116b9eca265Ssarens         if (pSubSize) {
117b9eca265Ssarens           if (bs < 0) {
118b9eca265Ssarens             bs = pSubSize;
119b9eca265Ssarens           } else if (bs != pSubSize) {
120b9eca265Ssarens             /* Layout does not admit a pointwise block size */
121b9eca265Ssarens             bs = 1;
1224d9407bcSMatthew G. Knepley           }
1234d9407bcSMatthew G. Knepley         }
1244d9407bcSMatthew G. Knepley       }
125b9eca265Ssarens     }
126b9eca265Ssarens     /* Must have same blocksize on all procs (some might have no points) */
1270be3e97aSMatthew G. Knepley     bsLocal[0] = bs < 0 ? PETSC_MAX_INT : bs; bsLocal[1] = bs;
1280be3e97aSMatthew G. Knepley     ierr = PetscGlobalMinMaxInt(PetscObjectComm((PetscObject) dm), bsLocal, bsMinMax);CHKERRQ(ierr);
1290be3e97aSMatthew G. Knepley     if (bsMinMax[0] != bsMinMax[1]) {bs = 1;}
1300be3e97aSMatthew G. Knepley     else                            {bs = bsMinMax[0];}
131785e854fSJed Brown     ierr = PetscMalloc1(subSize, &subIndices);CHKERRQ(ierr);
1324d9407bcSMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
1334d9407bcSMatthew G. Knepley       PetscInt gdof, goff;
1344d9407bcSMatthew G. Knepley 
1354d9407bcSMatthew G. Knepley       ierr = PetscSectionGetDof(sectionGlobal, p, &gdof);CHKERRQ(ierr);
1364d9407bcSMatthew G. Knepley       if (gdof > 0) {
1374d9407bcSMatthew G. Knepley         ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr);
1384d9407bcSMatthew G. Knepley         for (f = 0; f < numFields; ++f) {
1394d9407bcSMatthew G. Knepley           PetscInt fdof, fcdof, fc, f2, poff = 0;
1404d9407bcSMatthew G. Knepley 
1414d9407bcSMatthew G. Knepley           /* Can get rid of this loop by storing field information in the global section */
1424d9407bcSMatthew G. Knepley           for (f2 = 0; f2 < fields[f]; ++f2) {
1434d9407bcSMatthew G. Knepley             ierr  = PetscSectionGetFieldDof(section, p, f2, &fdof);CHKERRQ(ierr);
1444d9407bcSMatthew G. Knepley             ierr  = PetscSectionGetFieldConstraintDof(section, p, f2, &fcdof);CHKERRQ(ierr);
1454d9407bcSMatthew G. Knepley             poff += fdof-fcdof;
1464d9407bcSMatthew G. Knepley           }
1474d9407bcSMatthew G. Knepley           ierr = PetscSectionGetFieldDof(section, p, fields[f], &fdof);CHKERRQ(ierr);
1484d9407bcSMatthew G. Knepley           ierr = PetscSectionGetFieldConstraintDof(section, p, fields[f], &fcdof);CHKERRQ(ierr);
1494d9407bcSMatthew G. Knepley           for (fc = 0; fc < fdof-fcdof; ++fc, ++subOff) {
1504d9407bcSMatthew G. Knepley             subIndices[subOff] = goff+poff+fc;
1514d9407bcSMatthew G. Knepley           }
1524d9407bcSMatthew G. Knepley         }
1534d9407bcSMatthew G. Knepley       }
1544d9407bcSMatthew G. Knepley     }
1554d9407bcSMatthew G. Knepley     ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm), subSize, subIndices, PETSC_OWN_POINTER, is);CHKERRQ(ierr);
156627349f5SMatthew G. Knepley     if (bs > 1) {
157627349f5SMatthew G. Knepley       /* We need to check that the block size does not come from non-contiguous fields */
158627349f5SMatthew G. Knepley       PetscInt i, j, set = 1;
159627349f5SMatthew G. Knepley       for (i = 0; i < subSize; i += bs) {
160627349f5SMatthew G. Knepley         for (j = 0; j < bs; ++j) {
161627349f5SMatthew G. Knepley           if (subIndices[i+j] != subIndices[i]+j) {set = 0; break;}
162627349f5SMatthew G. Knepley         }
163627349f5SMatthew G. Knepley       }
164627349f5SMatthew G. Knepley       if (set) {ierr = ISSetBlockSize(*is, bs);CHKERRQ(ierr);}
165627349f5SMatthew G. Knepley     }
1664d9407bcSMatthew G. Knepley   }
1674d9407bcSMatthew G. Knepley   if (subdm) {
1684d9407bcSMatthew G. Knepley     PetscSection subsection;
1694d9407bcSMatthew G. Knepley     PetscBool    haveNull = PETSC_FALSE;
1704d9407bcSMatthew G. Knepley     PetscInt     f, nf = 0;
1714d9407bcSMatthew G. Knepley 
1724d9407bcSMatthew G. Knepley     ierr = PetscSectionCreateSubsection(section, numFields, fields, &subsection);CHKERRQ(ierr);
173*e87a4003SBarry Smith     ierr = DMSetSection(*subdm, subsection);CHKERRQ(ierr);
1744d9407bcSMatthew G. Knepley     ierr = PetscSectionDestroy(&subsection);CHKERRQ(ierr);
1754d9407bcSMatthew G. Knepley     for (f = 0; f < numFields; ++f) {
1764d9407bcSMatthew G. Knepley       (*subdm)->nullspaceConstructors[f] = dm->nullspaceConstructors[fields[f]];
1774d9407bcSMatthew G. Knepley       if ((*subdm)->nullspaceConstructors[f]) {
1784d9407bcSMatthew G. Knepley         haveNull = PETSC_TRUE;
1794d9407bcSMatthew G. Knepley         nf       = f;
1804d9407bcSMatthew G. Knepley       }
1814d9407bcSMatthew G. Knepley     }
182f646a522SMatthew G. Knepley     if (haveNull && is) {
1834d9407bcSMatthew G. Knepley       MatNullSpace nullSpace;
1844d9407bcSMatthew G. Knepley 
1854d9407bcSMatthew G. Knepley       ierr = (*(*subdm)->nullspaceConstructors[nf])(*subdm, nf, &nullSpace);CHKERRQ(ierr);
1864d9407bcSMatthew G. Knepley       ierr = PetscObjectCompose((PetscObject) *is, "nullspace", (PetscObject) nullSpace);CHKERRQ(ierr);
1874d9407bcSMatthew G. Knepley       ierr = MatNullSpaceDestroy(&nullSpace);CHKERRQ(ierr);
1884d9407bcSMatthew G. Knepley     }
1890f21e855SMatthew G. Knepley     if (dm->prob) {
1900f21e855SMatthew G. Knepley       PetscInt Nf;
1910f21e855SMatthew G. Knepley 
1922764a2aaSMatthew G. Knepley       ierr = PetscDSGetNumFields(dm->prob, &Nf);CHKERRQ(ierr);
1930f21e855SMatthew G. Knepley       if (nF != Nf) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "The number of DM fields %d does not match the number of Section fields %d", Nf, nF);
1944d9407bcSMatthew G. Knepley       ierr = DMSetNumFields(*subdm, numFields);CHKERRQ(ierr);
1954d9407bcSMatthew G. Knepley       for (f = 0; f < numFields; ++f) {
1960f21e855SMatthew G. Knepley         PetscObject disc;
1970f21e855SMatthew G. Knepley 
1980f21e855SMatthew G. Knepley         ierr = DMGetField(dm, fields[f], &disc);CHKERRQ(ierr);
1990f21e855SMatthew G. Knepley         ierr = DMSetField(*subdm, f, disc);CHKERRQ(ierr);
2004d9407bcSMatthew G. Knepley       }
201f646a522SMatthew G. Knepley       if (numFields == 1 && is) {
2020f21e855SMatthew G. Knepley         PetscObject disc, space, pmat;
2034d9407bcSMatthew G. Knepley 
2040f21e855SMatthew G. Knepley         ierr = DMGetField(*subdm, 0, &disc);CHKERRQ(ierr);
2050f21e855SMatthew G. Knepley         ierr = PetscObjectQuery(disc, "nullspace", &space);CHKERRQ(ierr);
2060f21e855SMatthew G. Knepley         if (space) {ierr = PetscObjectCompose((PetscObject) *is, "nullspace", space);CHKERRQ(ierr);}
2070f21e855SMatthew G. Knepley         ierr = PetscObjectQuery(disc, "nearnullspace", &space);CHKERRQ(ierr);
2080f21e855SMatthew G. Knepley         if (space) {ierr = PetscObjectCompose((PetscObject) *is, "nearnullspace", space);CHKERRQ(ierr);}
2090f21e855SMatthew G. Knepley         ierr = PetscObjectQuery(disc, "pmat", &pmat);CHKERRQ(ierr);
2100f21e855SMatthew G. Knepley         if (pmat) {ierr = PetscObjectCompose((PetscObject) *is, "pmat", pmat);CHKERRQ(ierr);}
2114d9407bcSMatthew G. Knepley       }
2124d9407bcSMatthew G. Knepley     }
2134d9407bcSMatthew G. Knepley   }
2144d9407bcSMatthew G. Knepley   PetscFunctionReturn(0);
2154d9407bcSMatthew G. Knepley }
2162adcc780SMatthew G. Knepley 
2172adcc780SMatthew G. Knepley /* This assumes that the DM has been cloned prior to the call */
2182adcc780SMatthew G. Knepley PetscErrorCode DMCreateSuperDM_Section_Private(DM dms[], PetscInt len, IS **is, DM *superdm)
2192adcc780SMatthew G. Knepley {
2202adcc780SMatthew G. Knepley   PetscSection  *sections, *sectionGlobals;
2212adcc780SMatthew G. Knepley   PetscInt      *Nfs, Nf = 0, *subIndices, i;
2222adcc780SMatthew G. Knepley   PetscErrorCode ierr;
2232adcc780SMatthew G. Knepley 
2242adcc780SMatthew G. Knepley   PetscFunctionBegin;
2252adcc780SMatthew G. Knepley   ierr = PetscMalloc3(len, &Nfs, len, &sections, len, &sectionGlobals);CHKERRQ(ierr);
2262adcc780SMatthew G. Knepley   for (i = 0 ; i < len; ++i) {
227*e87a4003SBarry Smith     ierr = DMGetSection(dms[i], &sections[i]);CHKERRQ(ierr);
228*e87a4003SBarry Smith     ierr = DMGetGlobalSection(dms[i], &sectionGlobals[i]);CHKERRQ(ierr);
2292adcc780SMatthew G. Knepley     if (!sections[i]) SETERRQ(PetscObjectComm((PetscObject)dms[0]), PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields");
2302adcc780SMatthew G. Knepley     if (!sectionGlobals[i]) SETERRQ(PetscObjectComm((PetscObject)dms[0]), PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields");
2312adcc780SMatthew G. Knepley     ierr = PetscSectionGetNumFields(sections[i], &Nfs[i]);CHKERRQ(ierr);
2322adcc780SMatthew G. Knepley     Nf += Nfs[i];
2332adcc780SMatthew G. Knepley   }
2342adcc780SMatthew G. Knepley   if (is) {
2352adcc780SMatthew G. Knepley     PetscInt *offs, *globalOffs, iOff = 0;
2362adcc780SMatthew G. Knepley 
2372adcc780SMatthew G. Knepley     ierr = PetscMalloc1(len, is);CHKERRQ(ierr);
2382adcc780SMatthew G. Knepley     ierr = PetscCalloc2(len+1, &offs, len+1, &globalOffs);CHKERRQ(ierr);
2392adcc780SMatthew G. Knepley     for (i = 0 ; i < len; ++i) {
2402adcc780SMatthew G. Knepley       ierr = PetscSectionGetConstrainedStorageSize(sectionGlobals[i], &offs[i]);CHKERRQ(ierr);
2412adcc780SMatthew G. Knepley       offs[len] += offs[i];
2422adcc780SMatthew G. Knepley     }
2432adcc780SMatthew G. Knepley     ierr = MPI_Scan(offs, globalOffs, len+1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) dms[0]));CHKERRQ(ierr);
2442adcc780SMatthew G. Knepley     for (i = 0 ; i <= len; ++i) globalOffs[i] -= offs[i];
2452adcc780SMatthew G. Knepley     for (i = 0 ; i < len; ++i, iOff += offs[i-1]) {
2462adcc780SMatthew G. Knepley       PetscInt bs = -1, bsLocal[2], bsMinMax[2];
2472adcc780SMatthew G. Knepley       PetscInt subSize = 0, subOff = 0, gtoff = globalOffs[len] - globalOffs[i], pStart, pEnd, p;
2482adcc780SMatthew G. Knepley 
2492adcc780SMatthew G. Knepley       ierr = PetscSectionGetChart(sectionGlobals[i], &pStart, &pEnd);CHKERRQ(ierr);
2502adcc780SMatthew G. Knepley       ierr = PetscSectionGetConstrainedStorageSize(sectionGlobals[i], &subSize);CHKERRQ(ierr);
2512adcc780SMatthew G. Knepley       ierr = PetscMalloc1(subSize, &subIndices);CHKERRQ(ierr);
2522adcc780SMatthew G. Knepley       for (p = pStart; p < pEnd; ++p) {
2532adcc780SMatthew G. Knepley         PetscInt gdof, gcdof, gtdof, goff, d;
2542adcc780SMatthew G. Knepley 
2552adcc780SMatthew G. Knepley         ierr = PetscSectionGetDof(sectionGlobals[i], p, &gdof);CHKERRQ(ierr);
2562adcc780SMatthew G. Knepley         ierr = PetscSectionGetConstraintDof(sections[i], p, &gcdof);CHKERRQ(ierr);
2572adcc780SMatthew G. Knepley         gtdof = gdof-gcdof;
2582adcc780SMatthew G. Knepley         if (gdof > 0 && gtdof) {
2592adcc780SMatthew G. Knepley           if (bs < 0)           {bs = gtdof;}
2602adcc780SMatthew G. Knepley           else if (bs != gtdof) {bs = 1;}
2612adcc780SMatthew G. Knepley           ierr = PetscSectionGetOffset(sectionGlobals[i], p, &goff);CHKERRQ(ierr);
2622adcc780SMatthew G. Knepley           for (d = 0; d < gtdof; ++d, ++subOff) {
2632adcc780SMatthew G. Knepley             subIndices[subOff] = goff+gtoff+d+iOff;
2642adcc780SMatthew G. Knepley           }
2652adcc780SMatthew G. Knepley         }
2662adcc780SMatthew G. Knepley       }
2672adcc780SMatthew G. Knepley       ierr = ISCreateGeneral(PetscObjectComm((PetscObject) dms[0]), subSize, subIndices, PETSC_OWN_POINTER, &(*is)[i]);CHKERRQ(ierr);
2682adcc780SMatthew G. Knepley       /* Must have same blocksize on all procs (some might have no points) */
2692adcc780SMatthew G. Knepley       bsLocal[0] = bs < 0 ? PETSC_MAX_INT : bs; bsLocal[1] = bs;
2702adcc780SMatthew G. Knepley       ierr = PetscGlobalMinMaxInt(PetscObjectComm((PetscObject) dms[0]), bsLocal, bsMinMax);CHKERRQ(ierr);
2712adcc780SMatthew G. Knepley       if (bsMinMax[0] != bsMinMax[1]) {bs = 1;}
2722adcc780SMatthew G. Knepley       else                            {bs = bsMinMax[0];}
2732adcc780SMatthew G. Knepley       ierr = ISSetBlockSize((*is)[i], bs);CHKERRQ(ierr);
2742adcc780SMatthew G. Knepley     }
27595b1d6b7SMatthew G. Knepley     ierr = PetscFree2(offs, globalOffs);CHKERRQ(ierr);
2762adcc780SMatthew G. Knepley   }
2772adcc780SMatthew G. Knepley   if (superdm) {
2782adcc780SMatthew G. Knepley     PetscSection supersection;
2792adcc780SMatthew G. Knepley     PetscBool    haveNull = PETSC_FALSE;
2802adcc780SMatthew G. Knepley     PetscInt     field, f, nf = 0;
2812adcc780SMatthew G. Knepley 
2822adcc780SMatthew G. Knepley     ierr = PetscSectionCreateSupersection(sections, len, &supersection);CHKERRQ(ierr);
283*e87a4003SBarry Smith     ierr = DMSetSection(*superdm, supersection);CHKERRQ(ierr);
2842adcc780SMatthew G. Knepley     ierr = PetscSectionDestroy(&supersection);CHKERRQ(ierr);
2852adcc780SMatthew G. Knepley     for (i = 0, field = 0; i < len; ++i) {
2862adcc780SMatthew G. Knepley       for (f = 0; f < Nfs[i]; ++f, ++field) {
2872adcc780SMatthew G. Knepley         (*superdm)->nullspaceConstructors[field] = dms[i]->nullspaceConstructors[f];
2882adcc780SMatthew G. Knepley         if ((*superdm)->nullspaceConstructors[field]) {
2892adcc780SMatthew G. Knepley           haveNull = PETSC_TRUE;
2902adcc780SMatthew G. Knepley           nf       = field;
2912adcc780SMatthew G. Knepley         }
2922adcc780SMatthew G. Knepley       }
2932adcc780SMatthew G. Knepley     }
2942adcc780SMatthew G. Knepley     if (haveNull && is) {
2952adcc780SMatthew G. Knepley       MatNullSpace nullSpace;
2962adcc780SMatthew G. Knepley 
2972adcc780SMatthew G. Knepley       ierr = (*(*superdm)->nullspaceConstructors[nf])(*superdm, nf, &nullSpace);CHKERRQ(ierr);
2982adcc780SMatthew G. Knepley       ierr = PetscObjectCompose((PetscObject) (*is)[nf], "nullspace", (PetscObject) nullSpace);CHKERRQ(ierr);
2992adcc780SMatthew G. Knepley       ierr = MatNullSpaceDestroy(&nullSpace);CHKERRQ(ierr);
3002adcc780SMatthew G. Knepley     }
3012adcc780SMatthew G. Knepley     if (len && dms[0]->prob) {
3022adcc780SMatthew G. Knepley       ierr = DMSetNumFields(*superdm, Nf);CHKERRQ(ierr);
3032adcc780SMatthew G. Knepley       for (i = 0, field = 0; i < len; ++i) {
3042adcc780SMatthew G. Knepley         for (f = 0; f < Nfs[i]; ++f, ++field) {
3052adcc780SMatthew G. Knepley           PetscObject disc;
3062adcc780SMatthew G. Knepley 
3072adcc780SMatthew G. Knepley           ierr = DMGetField(dms[i], f, &disc);CHKERRQ(ierr);
3082adcc780SMatthew G. Knepley           ierr = DMSetField(*superdm, field, disc);CHKERRQ(ierr);
3092adcc780SMatthew G. Knepley         }
3102adcc780SMatthew G. Knepley       }
3112adcc780SMatthew G. Knepley     }
3122adcc780SMatthew G. Knepley   }
31395b1d6b7SMatthew G. Knepley   ierr = PetscFree3(Nfs, sections, sectionGlobals);CHKERRQ(ierr);
3142adcc780SMatthew G. Knepley   PetscFunctionReturn(0);
3152adcc780SMatthew G. Knepley }
316