xref: /petsc/src/dm/interface/dmi.c (revision d8d19677bbccf95218448bee62e6b87f4513e133)
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;
12e87a4003SBarry 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;
34820f2d46SBarry Smith   ierr = MPIU_Allreduce(in,out,2,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)dm));CHKERRMPI(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;
6392fd8e1eSJed Brown   ierr = DMGetLocalSection(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 
84792b654fSMatthew G. Knepley /*@C
85792b654fSMatthew G. Knepley   DMCreateSectionSubDM - Returns an IS and subDM+subSection encapsulating a subproblem defined by the fields in a PetscSection in the DM.
86792b654fSMatthew G. Knepley 
87792b654fSMatthew G. Knepley   Not collective
88792b654fSMatthew G. Knepley 
89792b654fSMatthew G. Knepley   Input Parameters:
90792b654fSMatthew G. Knepley + dm        - The DM object
91792b654fSMatthew G. Knepley . numFields - The number of fields in this subproblem
92792b654fSMatthew G. Knepley - fields    - The field numbers of the selected fields
93792b654fSMatthew G. Knepley 
94792b654fSMatthew G. Knepley   Output Parameters:
95792b654fSMatthew G. Knepley + is - The global indices for the subproblem
96792b654fSMatthew G. Knepley - subdm - The DM for the subproblem, which must already have be cloned from dm
97792b654fSMatthew G. Knepley 
98792b654fSMatthew G. Knepley   Note: This handles all information in the DM class and the PetscSection. This is used as the basis for creating subDMs in specialized classes,
99792b654fSMatthew G. Knepley   such as Plex and Forest.
100792b654fSMatthew G. Knepley 
101792b654fSMatthew G. Knepley   Level: intermediate
102792b654fSMatthew G. Knepley 
10392fd8e1eSJed Brown .seealso DMCreateSubDM(), DMGetLocalSection(), DMPlexSetMigrationSF(), DMView()
104792b654fSMatthew G. Knepley @*/
105792b654fSMatthew G. Knepley PetscErrorCode DMCreateSectionSubDM(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm)
1064d9407bcSMatthew G. Knepley {
1074d9407bcSMatthew G. Knepley   PetscSection   section, sectionGlobal;
1084d9407bcSMatthew G. Knepley   PetscInt      *subIndices;
109696188eaSMatthew G. Knepley   PetscInt       subSize = 0, subOff = 0, Nf, f, pStart, pEnd, p;
1104d9407bcSMatthew G. Knepley   PetscErrorCode ierr;
1114d9407bcSMatthew G. Knepley 
1124d9407bcSMatthew G. Knepley   PetscFunctionBegin;
1134d9407bcSMatthew G. Knepley   if (!numFields) PetscFunctionReturn(0);
11492fd8e1eSJed Brown   ierr = DMGetLocalSection(dm, &section);CHKERRQ(ierr);
115e87a4003SBarry Smith   ierr = DMGetGlobalSection(dm, &sectionGlobal);CHKERRQ(ierr);
1164d9407bcSMatthew G. Knepley   if (!section) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields");
1174d9407bcSMatthew G. Knepley   if (!sectionGlobal) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields");
118696188eaSMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
119696188eaSMatthew 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);
1204d9407bcSMatthew G. Knepley   if (is) {
1213a544194SStefano Zampini     PetscInt bs, bsLocal[2], bsMinMax[2];
1223a544194SStefano Zampini 
1233a544194SStefano Zampini     for (f = 0, bs = 0; f < numFields; ++f) {
1243a544194SStefano Zampini       PetscInt Nc;
1253a544194SStefano Zampini 
1263a544194SStefano Zampini       ierr = PetscSectionGetFieldComponents(section, fields[f], &Nc);CHKERRQ(ierr);
1273a544194SStefano Zampini       bs  += Nc;
1283a544194SStefano Zampini     }
1294d9407bcSMatthew G. Knepley     ierr = PetscSectionGetChart(sectionGlobal, &pStart, &pEnd);CHKERRQ(ierr);
1304d9407bcSMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
131b9eca265Ssarens       PetscInt gdof, pSubSize  = 0;
1324d9407bcSMatthew G. Knepley 
1334d9407bcSMatthew G. Knepley       ierr = PetscSectionGetDof(sectionGlobal, p, &gdof);CHKERRQ(ierr);
1344d9407bcSMatthew G. Knepley       if (gdof > 0) {
1354d9407bcSMatthew G. Knepley         for (f = 0; f < numFields; ++f) {
1364d9407bcSMatthew G. Knepley           PetscInt fdof, fcdof;
1374d9407bcSMatthew G. Knepley 
1384d9407bcSMatthew G. Knepley           ierr     = PetscSectionGetFieldDof(section, p, fields[f], &fdof);CHKERRQ(ierr);
1394d9407bcSMatthew G. Knepley           ierr     = PetscSectionGetFieldConstraintDof(section, p, fields[f], &fcdof);CHKERRQ(ierr);
140b9eca265Ssarens           pSubSize += fdof-fcdof;
141b9eca265Ssarens         }
142b9eca265Ssarens         subSize += pSubSize;
1433a544194SStefano Zampini         if (pSubSize && bs != pSubSize) {
144b9eca265Ssarens           /* Layout does not admit a pointwise block size */
145b9eca265Ssarens           bs = 1;
1464d9407bcSMatthew G. Knepley         }
1474d9407bcSMatthew G. Knepley       }
1484d9407bcSMatthew G. Knepley     }
149b9eca265Ssarens     /* Must have same blocksize on all procs (some might have no points) */
1500be3e97aSMatthew G. Knepley     bsLocal[0] = bs < 0 ? PETSC_MAX_INT : bs; bsLocal[1] = bs;
1510be3e97aSMatthew G. Knepley     ierr = PetscGlobalMinMaxInt(PetscObjectComm((PetscObject) dm), bsLocal, bsMinMax);CHKERRQ(ierr);
1520be3e97aSMatthew G. Knepley     if (bsMinMax[0] != bsMinMax[1]) {bs = 1;}
1530be3e97aSMatthew G. Knepley     else                            {bs = bsMinMax[0];}
154785e854fSJed Brown     ierr = PetscMalloc1(subSize, &subIndices);CHKERRQ(ierr);
1554d9407bcSMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
1564d9407bcSMatthew G. Knepley       PetscInt gdof, goff;
1574d9407bcSMatthew G. Knepley 
1584d9407bcSMatthew G. Knepley       ierr = PetscSectionGetDof(sectionGlobal, p, &gdof);CHKERRQ(ierr);
1594d9407bcSMatthew G. Knepley       if (gdof > 0) {
1604d9407bcSMatthew G. Knepley         ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr);
1614d9407bcSMatthew G. Knepley         for (f = 0; f < numFields; ++f) {
1624d9407bcSMatthew G. Knepley           PetscInt fdof, fcdof, fc, f2, poff = 0;
1634d9407bcSMatthew G. Knepley 
1644d9407bcSMatthew G. Knepley           /* Can get rid of this loop by storing field information in the global section */
1654d9407bcSMatthew G. Knepley           for (f2 = 0; f2 < fields[f]; ++f2) {
1664d9407bcSMatthew G. Knepley             ierr  = PetscSectionGetFieldDof(section, p, f2, &fdof);CHKERRQ(ierr);
1674d9407bcSMatthew G. Knepley             ierr  = PetscSectionGetFieldConstraintDof(section, p, f2, &fcdof);CHKERRQ(ierr);
1684d9407bcSMatthew G. Knepley             poff += fdof-fcdof;
1694d9407bcSMatthew G. Knepley           }
1704d9407bcSMatthew G. Knepley           ierr = PetscSectionGetFieldDof(section, p, fields[f], &fdof);CHKERRQ(ierr);
1714d9407bcSMatthew G. Knepley           ierr = PetscSectionGetFieldConstraintDof(section, p, fields[f], &fcdof);CHKERRQ(ierr);
1724d9407bcSMatthew G. Knepley           for (fc = 0; fc < fdof-fcdof; ++fc, ++subOff) {
1734d9407bcSMatthew G. Knepley             subIndices[subOff] = goff+poff+fc;
1744d9407bcSMatthew G. Knepley           }
1754d9407bcSMatthew G. Knepley         }
1764d9407bcSMatthew G. Knepley       }
1774d9407bcSMatthew G. Knepley     }
1784d9407bcSMatthew G. Knepley     ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm), subSize, subIndices, PETSC_OWN_POINTER, is);CHKERRQ(ierr);
179627349f5SMatthew G. Knepley     if (bs > 1) {
180627349f5SMatthew G. Knepley       /* We need to check that the block size does not come from non-contiguous fields */
181627349f5SMatthew G. Knepley       PetscInt i, j, set = 1;
182627349f5SMatthew G. Knepley       for (i = 0; i < subSize; i += bs) {
183627349f5SMatthew G. Knepley         for (j = 0; j < bs; ++j) {
184627349f5SMatthew G. Knepley           if (subIndices[i+j] != subIndices[i]+j) {set = 0; break;}
185627349f5SMatthew G. Knepley         }
186627349f5SMatthew G. Knepley       }
187627349f5SMatthew G. Knepley       if (set) {ierr = ISSetBlockSize(*is, bs);CHKERRQ(ierr);}
188627349f5SMatthew G. Knepley     }
1894d9407bcSMatthew G. Knepley   }
1904d9407bcSMatthew G. Knepley   if (subdm) {
1914d9407bcSMatthew G. Knepley     PetscSection subsection;
1924d9407bcSMatthew G. Knepley     PetscBool    haveNull = PETSC_FALSE;
1938cda7954SMatthew G. Knepley     PetscInt     f, nf = 0, of = 0;
1944d9407bcSMatthew G. Knepley 
1954d9407bcSMatthew G. Knepley     ierr = PetscSectionCreateSubsection(section, numFields, fields, &subsection);CHKERRQ(ierr);
19692fd8e1eSJed Brown     ierr = DMSetLocalSection(*subdm, subsection);CHKERRQ(ierr);
1974d9407bcSMatthew G. Knepley     ierr = PetscSectionDestroy(&subsection);CHKERRQ(ierr);
1984d9407bcSMatthew G. Knepley     for (f = 0; f < numFields; ++f) {
1994d9407bcSMatthew G. Knepley       (*subdm)->nullspaceConstructors[f] = dm->nullspaceConstructors[fields[f]];
2004d9407bcSMatthew G. Knepley       if ((*subdm)->nullspaceConstructors[f]) {
2014d9407bcSMatthew G. Knepley         haveNull = PETSC_TRUE;
2024d9407bcSMatthew G. Knepley         nf       = f;
2038cda7954SMatthew G. Knepley         of       = fields[f];
2044d9407bcSMatthew G. Knepley       }
2054d9407bcSMatthew G. Knepley     }
206e5e52638SMatthew G. Knepley     if (dm->probs) {
2074d9407bcSMatthew G. Knepley       ierr = DMSetNumFields(*subdm, numFields);CHKERRQ(ierr);
2084d9407bcSMatthew G. Knepley       for (f = 0; f < numFields; ++f) {
2090f21e855SMatthew G. Knepley         PetscObject disc;
2100f21e855SMatthew G. Knepley 
21144a7f3ddSMatthew G. Knepley         ierr = DMGetField(dm, fields[f], NULL, &disc);CHKERRQ(ierr);
21244a7f3ddSMatthew G. Knepley         ierr = DMSetField(*subdm, f, NULL, disc);CHKERRQ(ierr);
2134d9407bcSMatthew G. Knepley       }
214e5e52638SMatthew G. Knepley       ierr = DMCreateDS(*subdm);CHKERRQ(ierr);
215f646a522SMatthew G. Knepley       if (numFields == 1 && is) {
2160f21e855SMatthew G. Knepley         PetscObject disc, space, pmat;
2174d9407bcSMatthew G. Knepley 
21844a7f3ddSMatthew G. Knepley         ierr = DMGetField(*subdm, 0, NULL, &disc);CHKERRQ(ierr);
2190f21e855SMatthew G. Knepley         ierr = PetscObjectQuery(disc, "nullspace", &space);CHKERRQ(ierr);
2200f21e855SMatthew G. Knepley         if (space) {ierr = PetscObjectCompose((PetscObject) *is, "nullspace", space);CHKERRQ(ierr);}
2210f21e855SMatthew G. Knepley         ierr = PetscObjectQuery(disc, "nearnullspace", &space);CHKERRQ(ierr);
2220f21e855SMatthew G. Knepley         if (space) {ierr = PetscObjectCompose((PetscObject) *is, "nearnullspace", space);CHKERRQ(ierr);}
2230f21e855SMatthew G. Knepley         ierr = PetscObjectQuery(disc, "pmat", &pmat);CHKERRQ(ierr);
2240f21e855SMatthew G. Knepley         if (pmat) {ierr = PetscObjectCompose((PetscObject) *is, "pmat", pmat);CHKERRQ(ierr);}
2254d9407bcSMatthew G. Knepley       }
2269310035eSMatthew G. Knepley       /* Check if DSes record their DM fields */
227e21d936aSMatthew G. Knepley       if (dm->probs[0].fields) {
2289310035eSMatthew G. Knepley         PetscInt d, e;
2299310035eSMatthew G. Knepley 
2309310035eSMatthew G. Knepley         for (d = 0, e = 0; d < dm->Nds && e < (*subdm)->Nds; ++d) {
2319310035eSMatthew G. Knepley           const PetscInt  Nf = dm->probs[d].ds->Nf;
2329310035eSMatthew G. Knepley           const PetscInt *fld;
2339310035eSMatthew G. Knepley           PetscInt        f, g;
2349310035eSMatthew G. Knepley 
2359310035eSMatthew G. Knepley           ierr = ISGetIndices(dm->probs[d].fields, &fld);CHKERRQ(ierr);
2369310035eSMatthew G. Knepley           for (f = 0; f < Nf; ++f) {
2379310035eSMatthew G. Knepley             for (g = 0; g < numFields; ++g) if (fld[f] == fields[g]) break;
2389310035eSMatthew G. Knepley             if (g < numFields) break;
2399310035eSMatthew G. Knepley           }
2409310035eSMatthew G. Knepley           ierr = ISRestoreIndices(dm->probs[d].fields, &fld);CHKERRQ(ierr);
2419310035eSMatthew G. Knepley           if (f == Nf) continue;
2429310035eSMatthew G. Knepley           ierr = PetscDSCopyConstants(dm->probs[d].ds, (*subdm)->probs[e].ds);CHKERRQ(ierr);
24336951cb5SMatthew G. Knepley           ierr = PetscDSCopyBoundary(dm->probs[d].ds, numFields, fields, (*subdm)->probs[e].ds);CHKERRQ(ierr);
2449310035eSMatthew G. Knepley           /* Translate DM fields to DS fields */
2459310035eSMatthew G. Knepley           {
246e21d936aSMatthew G. Knepley             IS              infields, dsfields;
24774a61aaaSMatthew G. Knepley             const PetscInt *fld, *ofld;
24874a61aaaSMatthew G. Knepley             PetscInt       *fidx;
24974a61aaaSMatthew G. Knepley             PetscInt        onf, nf, f, g;
250e21d936aSMatthew G. Knepley 
251e21d936aSMatthew G. Knepley             ierr = ISCreateGeneral(PETSC_COMM_SELF, numFields, fields, PETSC_USE_POINTER, &infields);CHKERRQ(ierr);
2529310035eSMatthew G. Knepley             ierr = ISIntersect(infields, dm->probs[d].fields, &dsfields);CHKERRQ(ierr);
253e21d936aSMatthew G. Knepley             ierr = ISDestroy(&infields);CHKERRQ(ierr);
254e21d936aSMatthew G. Knepley             ierr = ISGetLocalSize(dsfields, &nf);CHKERRQ(ierr);
255e21d936aSMatthew G. Knepley             if (!nf) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "DS cannot be supported on 0 fields");
256e21d936aSMatthew G. Knepley             ierr = ISGetIndices(dsfields, &fld);CHKERRQ(ierr);
2579310035eSMatthew G. Knepley             ierr = ISGetLocalSize(dm->probs[d].fields, &onf);CHKERRQ(ierr);
2589310035eSMatthew G. Knepley             ierr = ISGetIndices(dm->probs[d].fields, &ofld);CHKERRQ(ierr);
25974a61aaaSMatthew G. Knepley             ierr = PetscMalloc1(nf, &fidx);CHKERRQ(ierr);
2603268a0aaSMatthew G. Knepley             for (f = 0, g = 0; f < onf && g < nf; ++f) {
26174a61aaaSMatthew G. Knepley               if (ofld[f] == fld[g]) fidx[g++] = f;
26274a61aaaSMatthew G. Knepley             }
2639310035eSMatthew G. Knepley             ierr = ISRestoreIndices(dm->probs[d].fields, &ofld);CHKERRQ(ierr);
264e21d936aSMatthew G. Knepley             ierr = ISRestoreIndices(dsfields, &fld);CHKERRQ(ierr);
265e21d936aSMatthew G. Knepley             ierr = ISDestroy(&dsfields);CHKERRQ(ierr);
2666c1eb96dSMatthew G. Knepley             ierr = PetscDSSelectDiscretizations(dm->probs[0].ds, nf, fidx, (*subdm)->probs[0].ds);CHKERRQ(ierr);
26774a61aaaSMatthew G. Knepley             ierr = PetscDSSelectEquations(dm->probs[0].ds, nf, fidx, (*subdm)->probs[0].ds);CHKERRQ(ierr);
26874a61aaaSMatthew G. Knepley             ierr = PetscFree(fidx);CHKERRQ(ierr);
2699310035eSMatthew G. Knepley           }
2709310035eSMatthew G. Knepley           ++e;
2719310035eSMatthew G. Knepley         }
272e21d936aSMatthew G. Knepley       } else {
2739310035eSMatthew G. Knepley         ierr = PetscDSCopyConstants(dm->probs[0].ds, (*subdm)->probs[0].ds);CHKERRQ(ierr);
27436951cb5SMatthew G. Knepley         ierr = PetscDSCopyBoundary(dm->probs[0].ds, PETSC_DETERMINE, NULL, (*subdm)->probs[0].ds);CHKERRQ(ierr);
2756c1eb96dSMatthew G. Knepley         ierr = PetscDSSelectDiscretizations(dm->probs[0].ds, numFields, fields, (*subdm)->probs[0].ds);CHKERRQ(ierr);
276e5e52638SMatthew G. Knepley         ierr = PetscDSSelectEquations(dm->probs[0].ds, numFields, fields, (*subdm)->probs[0].ds);CHKERRQ(ierr);
277d5af271aSMatthew G. Knepley       }
278e21d936aSMatthew G. Knepley     }
2798cda7954SMatthew G. Knepley     if (haveNull && is) {
2808cda7954SMatthew G. Knepley       MatNullSpace nullSpace;
2818cda7954SMatthew G. Knepley 
2828cda7954SMatthew G. Knepley       ierr = (*(*subdm)->nullspaceConstructors[nf])(*subdm, of, nf, &nullSpace);CHKERRQ(ierr);
2838cda7954SMatthew G. Knepley       ierr = PetscObjectCompose((PetscObject) *is, "nullspace", (PetscObject) nullSpace);CHKERRQ(ierr);
2848cda7954SMatthew G. Knepley       ierr = MatNullSpaceDestroy(&nullSpace);CHKERRQ(ierr);
2858cda7954SMatthew G. Knepley     }
286d5af271aSMatthew G. Knepley     if (dm->coarseMesh) {
287d5af271aSMatthew G. Knepley       ierr = DMCreateSubDM(dm->coarseMesh, numFields, fields, NULL, &(*subdm)->coarseMesh);CHKERRQ(ierr);
2884d9407bcSMatthew G. Knepley     }
2894d9407bcSMatthew G. Knepley   }
2904d9407bcSMatthew G. Knepley   PetscFunctionReturn(0);
2914d9407bcSMatthew G. Knepley }
2922adcc780SMatthew G. Knepley 
293792b654fSMatthew G. Knepley /*@C
294792b654fSMatthew G. Knepley   DMCreateSectionSuperDM - Returns an arrays of ISes and DM+Section encapsulating a superproblem defined by the DM+Sections passed in.
295792b654fSMatthew G. Knepley 
296792b654fSMatthew G. Knepley   Not collective
297792b654fSMatthew G. Knepley 
298*d8d19677SJose E. Roman   Input Parameters:
299792b654fSMatthew G. Knepley + dms - The DM objects
300792b654fSMatthew G. Knepley - len - The number of DMs
301792b654fSMatthew G. Knepley 
302792b654fSMatthew G. Knepley   Output Parameters:
303792b654fSMatthew G. Knepley + is - The global indices for the subproblem, or NULL
304792b654fSMatthew G. Knepley - superdm - The DM for the superproblem, which must already have be cloned
305792b654fSMatthew G. Knepley 
306792b654fSMatthew G. Knepley   Note: This handles all information in the DM class and the PetscSection. This is used as the basis for creating subDMs in specialized classes,
307792b654fSMatthew G. Knepley   such as Plex and Forest.
308792b654fSMatthew G. Knepley 
309792b654fSMatthew G. Knepley   Level: intermediate
310792b654fSMatthew G. Knepley 
31192fd8e1eSJed Brown .seealso DMCreateSuperDM(), DMGetLocalSection(), DMPlexSetMigrationSF(), DMView()
312792b654fSMatthew G. Knepley @*/
313792b654fSMatthew G. Knepley PetscErrorCode DMCreateSectionSuperDM(DM dms[], PetscInt len, IS **is, DM *superdm)
3142adcc780SMatthew G. Knepley {
3159796d432SMatthew G. Knepley   MPI_Comm       comm;
3169796d432SMatthew G. Knepley   PetscSection   supersection, *sections, *sectionGlobals;
3178cda7954SMatthew G. Knepley   PetscInt      *Nfs, Nf = 0, f, supf, oldf = -1, nullf = -1, i;
3189796d432SMatthew G. Knepley   PetscBool      haveNull = PETSC_FALSE;
3192adcc780SMatthew G. Knepley   PetscErrorCode ierr;
3202adcc780SMatthew G. Knepley 
3212adcc780SMatthew G. Knepley   PetscFunctionBegin;
3229796d432SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)dms[0], &comm);CHKERRQ(ierr);
3239796d432SMatthew G. Knepley   /* Pull out local and global sections */
3242adcc780SMatthew G. Knepley   ierr = PetscMalloc3(len, &Nfs, len, &sections, len, &sectionGlobals);CHKERRQ(ierr);
3252adcc780SMatthew G. Knepley   for (i = 0 ; i < len; ++i) {
32692fd8e1eSJed Brown     ierr = DMGetLocalSection(dms[i], &sections[i]);CHKERRQ(ierr);
327e87a4003SBarry Smith     ierr = DMGetGlobalSection(dms[i], &sectionGlobals[i]);CHKERRQ(ierr);
3289796d432SMatthew G. Knepley     if (!sections[i]) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields");
3299796d432SMatthew G. Knepley     if (!sectionGlobals[i]) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields");
3302adcc780SMatthew G. Knepley     ierr = PetscSectionGetNumFields(sections[i], &Nfs[i]);CHKERRQ(ierr);
3312adcc780SMatthew G. Knepley     Nf += Nfs[i];
3322adcc780SMatthew G. Knepley   }
3339796d432SMatthew G. Knepley   /* Create the supersection */
3349796d432SMatthew G. Knepley   ierr = PetscSectionCreateSupersection(sections, len, &supersection);CHKERRQ(ierr);
33592fd8e1eSJed Brown   ierr = DMSetLocalSection(*superdm, supersection);CHKERRQ(ierr);
3369796d432SMatthew G. Knepley   /* Create ISes */
3372adcc780SMatthew G. Knepley   if (is) {
3389796d432SMatthew G. Knepley     PetscSection supersectionGlobal;
3399796d432SMatthew G. Knepley     PetscInt     bs = -1, startf = 0;
3402adcc780SMatthew G. Knepley 
3412adcc780SMatthew G. Knepley     ierr = PetscMalloc1(len, is);CHKERRQ(ierr);
3426f0eb057SJed Brown     ierr = DMGetGlobalSection(*superdm, &supersectionGlobal);CHKERRQ(ierr);
3439796d432SMatthew G. Knepley     for (i = 0 ; i < len; startf += Nfs[i], ++i) {
3449796d432SMatthew G. Knepley       PetscInt *subIndices;
345ec4c761aSMatthew G. Knepley       PetscInt  subSize, subOff, pStart, pEnd, p, start, end, dummy;
3462adcc780SMatthew G. Knepley 
3472adcc780SMatthew G. Knepley       ierr = PetscSectionGetChart(sectionGlobals[i], &pStart, &pEnd);CHKERRQ(ierr);
3482adcc780SMatthew G. Knepley       ierr = PetscSectionGetConstrainedStorageSize(sectionGlobals[i], &subSize);CHKERRQ(ierr);
3492adcc780SMatthew G. Knepley       ierr = PetscMalloc1(subSize, &subIndices);CHKERRQ(ierr);
3509796d432SMatthew G. Knepley       for (p = pStart, subOff = 0; p < pEnd; ++p) {
3519796d432SMatthew G. Knepley         PetscInt gdof, gcdof, gtdof, d;
3522adcc780SMatthew G. Knepley 
3532adcc780SMatthew G. Knepley         ierr = PetscSectionGetDof(sectionGlobals[i], p, &gdof);CHKERRQ(ierr);
3542adcc780SMatthew G. Knepley         ierr = PetscSectionGetConstraintDof(sections[i], p, &gcdof);CHKERRQ(ierr);
3552adcc780SMatthew G. Knepley         gtdof = gdof - gcdof;
3562adcc780SMatthew G. Knepley         if (gdof > 0 && gtdof) {
3572adcc780SMatthew G. Knepley           if (bs < 0)           {bs = gtdof;}
3582adcc780SMatthew G. Knepley           else if (bs != gtdof) {bs = 1;}
359ec4c761aSMatthew G. Knepley           ierr = DMGetGlobalFieldOffset_Private(*superdm, p, startf, &start, &dummy);CHKERRQ(ierr);
360ec4c761aSMatthew G. Knepley           ierr = DMGetGlobalFieldOffset_Private(*superdm, p, startf+Nfs[i]-1, &dummy, &end);CHKERRQ(ierr);
3619796d432SMatthew G. Knepley           if (end-start != gtdof) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid number of global dofs %D != %D for dm %D on point %D", end-start, gtdof, i, p);
3629796d432SMatthew G. Knepley           for (d = start; d < end; ++d, ++subOff) subIndices[subOff] = d;
3632adcc780SMatthew G. Knepley         }
3642adcc780SMatthew G. Knepley       }
3659796d432SMatthew G. Knepley       ierr = ISCreateGeneral(comm, subSize, subIndices, PETSC_OWN_POINTER, &(*is)[i]);CHKERRQ(ierr);
3662adcc780SMatthew G. Knepley       /* Must have same blocksize on all procs (some might have no points) */
3679796d432SMatthew G. Knepley       {
3689796d432SMatthew G. Knepley         PetscInt bs = -1, bsLocal[2], bsMinMax[2];
3699796d432SMatthew G. Knepley 
3702adcc780SMatthew G. Knepley         bsLocal[0] = bs < 0 ? PETSC_MAX_INT : bs; bsLocal[1] = bs;
3719796d432SMatthew G. Knepley         ierr = PetscGlobalMinMaxInt(comm, bsLocal, bsMinMax);CHKERRQ(ierr);
3722adcc780SMatthew G. Knepley         if (bsMinMax[0] != bsMinMax[1]) {bs = 1;}
3732adcc780SMatthew G. Knepley         else                            {bs = bsMinMax[0];}
3742adcc780SMatthew G. Knepley         ierr = ISSetBlockSize((*is)[i], bs);CHKERRQ(ierr);
3752adcc780SMatthew G. Knepley       }
3762adcc780SMatthew G. Knepley     }
3772adcc780SMatthew G. Knepley   }
3789796d432SMatthew G. Knepley   /* Preserve discretizations */
379e5e52638SMatthew G. Knepley   if (len && dms[0]->probs) {
3802adcc780SMatthew G. Knepley     ierr = DMSetNumFields(*superdm, Nf);CHKERRQ(ierr);
3819796d432SMatthew G. Knepley     for (i = 0, supf = 0; i < len; ++i) {
3829796d432SMatthew G. Knepley       for (f = 0; f < Nfs[i]; ++f, ++supf) {
3832adcc780SMatthew G. Knepley         PetscObject disc;
3842adcc780SMatthew G. Knepley 
38544a7f3ddSMatthew G. Knepley         ierr = DMGetField(dms[i], f, NULL, &disc);CHKERRQ(ierr);
38644a7f3ddSMatthew G. Knepley         ierr = DMSetField(*superdm, supf, NULL, disc);CHKERRQ(ierr);
3872adcc780SMatthew G. Knepley       }
3882adcc780SMatthew G. Knepley     }
389e5e52638SMatthew G. Knepley     ierr = DMCreateDS(*superdm);CHKERRQ(ierr);
3902adcc780SMatthew G. Knepley   }
3919796d432SMatthew G. Knepley   /* Preserve nullspaces */
3929796d432SMatthew G. Knepley   for (i = 0, supf = 0; i < len; ++i) {
3939796d432SMatthew G. Knepley     for (f = 0; f < Nfs[i]; ++f, ++supf) {
3949796d432SMatthew G. Knepley       (*superdm)->nullspaceConstructors[supf] = dms[i]->nullspaceConstructors[f];
3959796d432SMatthew G. Knepley       if ((*superdm)->nullspaceConstructors[supf]) {
3969796d432SMatthew G. Knepley         haveNull = PETSC_TRUE;
3979796d432SMatthew G. Knepley         nullf    = supf;
3988cda7954SMatthew G. Knepley         oldf     = f;
3992adcc780SMatthew G. Knepley       }
4009796d432SMatthew G. Knepley     }
4019796d432SMatthew G. Knepley   }
4029796d432SMatthew G. Knepley   /* Attach nullspace to IS */
4039796d432SMatthew G. Knepley   if (haveNull && is) {
4049796d432SMatthew G. Knepley     MatNullSpace nullSpace;
4059796d432SMatthew G. Knepley 
4068cda7954SMatthew G. Knepley     ierr = (*(*superdm)->nullspaceConstructors[nullf])(*superdm, oldf, nullf, &nullSpace);CHKERRQ(ierr);
4079796d432SMatthew G. Knepley     ierr = PetscObjectCompose((PetscObject) (*is)[nullf], "nullspace", (PetscObject) nullSpace);CHKERRQ(ierr);
4089796d432SMatthew G. Knepley     ierr = MatNullSpaceDestroy(&nullSpace);CHKERRQ(ierr);
4099796d432SMatthew G. Knepley   }
4109796d432SMatthew G. Knepley   ierr = PetscSectionDestroy(&supersection);CHKERRQ(ierr);
41195b1d6b7SMatthew G. Knepley   ierr = PetscFree3(Nfs, sections, sectionGlobals);CHKERRQ(ierr);
4122adcc780SMatthew G. Knepley   PetscFunctionReturn(0);
4132adcc780SMatthew G. Knepley }
414