xref: /petsc/src/dm/interface/dmi.c (revision 7a8be3513cf479fb6a672bd91de7eb6883fcfd02)
1af0996ceSBarry Smith #include <petsc/private/dmimpl.h>     /*I      "petscdm.h"     I*/
22764a2aaSMatthew G. Knepley #include <petscds.h>
311689aebSJed Brown 
45b8243e1SJed Brown // Greatest common divisor of two nonnegative integers
55b8243e1SJed Brown PetscInt PetscGCD(PetscInt a, PetscInt b) {
65b8243e1SJed Brown   while (b != 0) {
75b8243e1SJed Brown     PetscInt tmp = b;
85b8243e1SJed Brown     b = a % b;
95b8243e1SJed Brown     a = tmp;
105b8243e1SJed Brown   }
115b8243e1SJed Brown   return a;
125b8243e1SJed Brown }
135b8243e1SJed Brown 
1411689aebSJed Brown PetscErrorCode DMCreateGlobalVector_Section_Private(DM dm,Vec *vec)
1511689aebSJed Brown {
1611689aebSJed Brown   PetscSection   gSection;
17cc30b04dSMatthew G Knepley   PetscInt       localSize, bs, blockSize = -1, pStart, pEnd, p;
18fad22124SMatthew G Knepley   PetscErrorCode ierr;
195d1c0279SHong Zhang   PetscInt       in[2],out[2];
2011689aebSJed Brown 
2111689aebSJed Brown   PetscFunctionBegin;
22e87a4003SBarry Smith   ierr = DMGetGlobalSection(dm, &gSection);CHKERRQ(ierr);
23fad22124SMatthew G Knepley   ierr = PetscSectionGetChart(gSection, &pStart, &pEnd);CHKERRQ(ierr);
2411689aebSJed Brown   for (p = pStart; p < pEnd; ++p) {
2511689aebSJed Brown     PetscInt dof, cdof;
2611689aebSJed Brown 
27fad22124SMatthew G Knepley     ierr = PetscSectionGetDof(gSection, p, &dof);CHKERRQ(ierr);
28fad22124SMatthew G Knepley     ierr = PetscSectionGetConstraintDof(gSection, p, &cdof);CHKERRQ(ierr);
290680ce0aSHong Zhang 
305b8243e1SJed Brown     if (dof - cdof > 0) {
315b8243e1SJed Brown       if (blockSize < 0) {
320680ce0aSHong Zhang         /* set blockSize */
330680ce0aSHong Zhang         blockSize = dof-cdof;
345b8243e1SJed Brown       } else {
355b8243e1SJed Brown         blockSize = PetscGCD(dof - cdof, blockSize);
3611689aebSJed Brown       }
3711689aebSJed Brown     }
380680ce0aSHong Zhang   }
395d1c0279SHong Zhang 
402c8be454SMatthew G. Knepley   in[0] = blockSize < 0 ? PETSC_MIN_INT : -blockSize;
415d1c0279SHong Zhang   in[1] = blockSize;
42820f2d46SBarry Smith   ierr = MPIU_Allreduce(in,out,2,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)dm));CHKERRMPI(ierr);
435d1c0279SHong Zhang   /* -out[0] = min(blockSize), out[1] = max(blockSize) */
445d1c0279SHong Zhang   if (-out[0] == out[1]) {
455d1c0279SHong Zhang     bs = out[1];
465d1c0279SHong Zhang   } else bs = 1;
475d1c0279SHong Zhang 
485d1c0279SHong Zhang   if (bs < 0) { /* Everyone was empty */
495d1c0279SHong Zhang     blockSize = 1;
505d1c0279SHong Zhang     bs        = 1;
515d1c0279SHong Zhang   }
525d1c0279SHong Zhang 
53fad22124SMatthew G Knepley   ierr = PetscSectionGetConstrainedStorageSize(gSection, &localSize);CHKERRQ(ierr);
54*7a8be351SBarry Smith   PetscCheck(localSize%blockSize == 0,PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mismatch between blocksize %" PetscInt_FMT " and local storage size %" PetscInt_FMT, blockSize, localSize);
5582f516ccSBarry Smith   ierr = VecCreate(PetscObjectComm((PetscObject)dm), vec);CHKERRQ(ierr);
5611689aebSJed Brown   ierr = VecSetSizes(*vec, localSize, PETSC_DETERMINE);CHKERRQ(ierr);
57cc30b04dSMatthew G Knepley   ierr = VecSetBlockSize(*vec, bs);CHKERRQ(ierr);
58c0dedaeaSBarry Smith   ierr = VecSetType(*vec,dm->vectype);CHKERRQ(ierr);
5911689aebSJed Brown   ierr = VecSetDM(*vec, dm);CHKERRQ(ierr);
6011689aebSJed Brown   /* ierr = VecSetLocalToGlobalMapping(*vec, dm->ltogmap);CHKERRQ(ierr); */
6111689aebSJed Brown   PetscFunctionReturn(0);
6211689aebSJed Brown }
6311689aebSJed Brown 
6411689aebSJed Brown PetscErrorCode DMCreateLocalVector_Section_Private(DM dm,Vec *vec)
6511689aebSJed Brown {
66fad22124SMatthew G Knepley   PetscSection   section;
6711689aebSJed Brown   PetscInt       localSize, blockSize = -1, pStart, pEnd, p;
68fad22124SMatthew G Knepley   PetscErrorCode ierr;
6911689aebSJed Brown 
7011689aebSJed Brown   PetscFunctionBegin;
7192fd8e1eSJed Brown   ierr = DMGetLocalSection(dm, &section);CHKERRQ(ierr);
72fad22124SMatthew G Knepley   ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
7311689aebSJed Brown   for (p = pStart; p < pEnd; ++p) {
7411689aebSJed Brown     PetscInt dof;
7511689aebSJed Brown 
76fad22124SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
7711689aebSJed Brown     if ((blockSize < 0) && (dof > 0)) blockSize = dof;
785b8243e1SJed Brown     if (dof > 0) blockSize = PetscGCD(dof, blockSize);
7911689aebSJed Brown   }
80fad22124SMatthew G Knepley   ierr = PetscSectionGetStorageSize(section, &localSize);CHKERRQ(ierr);
8111689aebSJed Brown   ierr = VecCreate(PETSC_COMM_SELF, vec);CHKERRQ(ierr);
8211689aebSJed Brown   ierr = VecSetSizes(*vec, localSize, localSize);CHKERRQ(ierr);
8311689aebSJed Brown   ierr = VecSetBlockSize(*vec, blockSize);CHKERRQ(ierr);
84c0dedaeaSBarry Smith   ierr = VecSetType(*vec,dm->vectype);CHKERRQ(ierr);
8511689aebSJed Brown   ierr = VecSetDM(*vec, dm);CHKERRQ(ierr);
8611689aebSJed Brown   PetscFunctionReturn(0);
8711689aebSJed Brown }
884d9407bcSMatthew G. Knepley 
89792b654fSMatthew G. Knepley /*@C
90792b654fSMatthew G. Knepley   DMCreateSectionSubDM - Returns an IS and subDM+subSection encapsulating a subproblem defined by the fields in a PetscSection in the DM.
91792b654fSMatthew G. Knepley 
92792b654fSMatthew G. Knepley   Not collective
93792b654fSMatthew G. Knepley 
94792b654fSMatthew G. Knepley   Input Parameters:
95792b654fSMatthew G. Knepley + dm        - The DM object
96792b654fSMatthew G. Knepley . numFields - The number of fields in this subproblem
97792b654fSMatthew G. Knepley - fields    - The field numbers of the selected fields
98792b654fSMatthew G. Knepley 
99792b654fSMatthew G. Knepley   Output Parameters:
100792b654fSMatthew G. Knepley + is - The global indices for the subproblem
101792b654fSMatthew G. Knepley - subdm - The DM for the subproblem, which must already have be cloned from dm
102792b654fSMatthew G. Knepley 
103792b654fSMatthew 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,
104792b654fSMatthew G. Knepley   such as Plex and Forest.
105792b654fSMatthew G. Knepley 
106792b654fSMatthew G. Knepley   Level: intermediate
107792b654fSMatthew G. Knepley 
10892fd8e1eSJed Brown .seealso DMCreateSubDM(), DMGetLocalSection(), DMPlexSetMigrationSF(), DMView()
109792b654fSMatthew G. Knepley @*/
110792b654fSMatthew G. Knepley PetscErrorCode DMCreateSectionSubDM(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm)
1114d9407bcSMatthew G. Knepley {
1124d9407bcSMatthew G. Knepley   PetscSection   section, sectionGlobal;
1134d9407bcSMatthew G. Knepley   PetscInt      *subIndices;
114696188eaSMatthew G. Knepley   PetscInt       subSize = 0, subOff = 0, Nf, f, pStart, pEnd, p;
1154d9407bcSMatthew G. Knepley   PetscErrorCode ierr;
1164d9407bcSMatthew G. Knepley 
1174d9407bcSMatthew G. Knepley   PetscFunctionBegin;
1184d9407bcSMatthew G. Knepley   if (!numFields) PetscFunctionReturn(0);
11992fd8e1eSJed Brown   ierr = DMGetLocalSection(dm, &section);CHKERRQ(ierr);
120e87a4003SBarry Smith   ierr = DMGetGlobalSection(dm, &sectionGlobal);CHKERRQ(ierr);
121*7a8be351SBarry Smith   PetscCheck(section,PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields");
122*7a8be351SBarry Smith   PetscCheck(sectionGlobal,PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields");
123696188eaSMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
124*7a8be351SBarry Smith   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);
1254d9407bcSMatthew G. Knepley   if (is) {
1263a544194SStefano Zampini     PetscInt bs, bsLocal[2], bsMinMax[2];
1273a544194SStefano Zampini 
1283a544194SStefano Zampini     for (f = 0, bs = 0; f < numFields; ++f) {
1293a544194SStefano Zampini       PetscInt Nc;
1303a544194SStefano Zampini 
1313a544194SStefano Zampini       ierr = PetscSectionGetFieldComponents(section, fields[f], &Nc);CHKERRQ(ierr);
1323a544194SStefano Zampini       bs  += Nc;
1333a544194SStefano Zampini     }
1344d9407bcSMatthew G. Knepley     ierr = PetscSectionGetChart(sectionGlobal, &pStart, &pEnd);CHKERRQ(ierr);
1354d9407bcSMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
136b9eca265Ssarens       PetscInt gdof, pSubSize  = 0;
1374d9407bcSMatthew G. Knepley 
1384d9407bcSMatthew G. Knepley       ierr = PetscSectionGetDof(sectionGlobal, p, &gdof);CHKERRQ(ierr);
1394d9407bcSMatthew G. Knepley       if (gdof > 0) {
1404d9407bcSMatthew G. Knepley         for (f = 0; f < numFields; ++f) {
1414d9407bcSMatthew G. Knepley           PetscInt fdof, fcdof;
1424d9407bcSMatthew G. Knepley 
1434d9407bcSMatthew G. Knepley           ierr     = PetscSectionGetFieldDof(section, p, fields[f], &fdof);CHKERRQ(ierr);
1444d9407bcSMatthew G. Knepley           ierr     = PetscSectionGetFieldConstraintDof(section, p, fields[f], &fcdof);CHKERRQ(ierr);
145b9eca265Ssarens           pSubSize += fdof-fcdof;
146b9eca265Ssarens         }
147b9eca265Ssarens         subSize += pSubSize;
1483a544194SStefano Zampini         if (pSubSize && bs != pSubSize) {
149b9eca265Ssarens           /* Layout does not admit a pointwise block size */
150b9eca265Ssarens           bs = 1;
1514d9407bcSMatthew G. Knepley         }
1524d9407bcSMatthew G. Knepley       }
1534d9407bcSMatthew G. Knepley     }
154b9eca265Ssarens     /* Must have same blocksize on all procs (some might have no points) */
1550be3e97aSMatthew G. Knepley     bsLocal[0] = bs < 0 ? PETSC_MAX_INT : bs; bsLocal[1] = bs;
1560be3e97aSMatthew G. Knepley     ierr = PetscGlobalMinMaxInt(PetscObjectComm((PetscObject) dm), bsLocal, bsMinMax);CHKERRQ(ierr);
1570be3e97aSMatthew G. Knepley     if (bsMinMax[0] != bsMinMax[1]) {bs = 1;}
1580be3e97aSMatthew G. Knepley     else                            {bs = bsMinMax[0];}
159785e854fSJed Brown     ierr = PetscMalloc1(subSize, &subIndices);CHKERRQ(ierr);
1604d9407bcSMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
1614d9407bcSMatthew G. Knepley       PetscInt gdof, goff;
1624d9407bcSMatthew G. Knepley 
1634d9407bcSMatthew G. Knepley       ierr = PetscSectionGetDof(sectionGlobal, p, &gdof);CHKERRQ(ierr);
1644d9407bcSMatthew G. Knepley       if (gdof > 0) {
1654d9407bcSMatthew G. Knepley         ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr);
1664d9407bcSMatthew G. Knepley         for (f = 0; f < numFields; ++f) {
1674d9407bcSMatthew G. Knepley           PetscInt fdof, fcdof, fc, f2, poff = 0;
1684d9407bcSMatthew G. Knepley 
1694d9407bcSMatthew G. Knepley           /* Can get rid of this loop by storing field information in the global section */
1704d9407bcSMatthew G. Knepley           for (f2 = 0; f2 < fields[f]; ++f2) {
1714d9407bcSMatthew G. Knepley             ierr  = PetscSectionGetFieldDof(section, p, f2, &fdof);CHKERRQ(ierr);
1724d9407bcSMatthew G. Knepley             ierr  = PetscSectionGetFieldConstraintDof(section, p, f2, &fcdof);CHKERRQ(ierr);
1734d9407bcSMatthew G. Knepley             poff += fdof-fcdof;
1744d9407bcSMatthew G. Knepley           }
1754d9407bcSMatthew G. Knepley           ierr = PetscSectionGetFieldDof(section, p, fields[f], &fdof);CHKERRQ(ierr);
1764d9407bcSMatthew G. Knepley           ierr = PetscSectionGetFieldConstraintDof(section, p, fields[f], &fcdof);CHKERRQ(ierr);
1774d9407bcSMatthew G. Knepley           for (fc = 0; fc < fdof-fcdof; ++fc, ++subOff) {
1784d9407bcSMatthew G. Knepley             subIndices[subOff] = goff+poff+fc;
1794d9407bcSMatthew G. Knepley           }
1804d9407bcSMatthew G. Knepley         }
1814d9407bcSMatthew G. Knepley       }
1824d9407bcSMatthew G. Knepley     }
1834d9407bcSMatthew G. Knepley     ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm), subSize, subIndices, PETSC_OWN_POINTER, is);CHKERRQ(ierr);
184627349f5SMatthew G. Knepley     if (bs > 1) {
185627349f5SMatthew G. Knepley       /* We need to check that the block size does not come from non-contiguous fields */
186627349f5SMatthew G. Knepley       PetscInt i, j, set = 1;
187627349f5SMatthew G. Knepley       for (i = 0; i < subSize; i += bs) {
188627349f5SMatthew G. Knepley         for (j = 0; j < bs; ++j) {
189627349f5SMatthew G. Knepley           if (subIndices[i+j] != subIndices[i]+j) {set = 0; break;}
190627349f5SMatthew G. Knepley         }
191627349f5SMatthew G. Knepley       }
192627349f5SMatthew G. Knepley       if (set) {ierr = ISSetBlockSize(*is, bs);CHKERRQ(ierr);}
193627349f5SMatthew G. Knepley     }
1944d9407bcSMatthew G. Knepley   }
1954d9407bcSMatthew G. Knepley   if (subdm) {
1964d9407bcSMatthew G. Knepley     PetscSection subsection;
1974d9407bcSMatthew G. Knepley     PetscBool    haveNull = PETSC_FALSE;
1988cda7954SMatthew G. Knepley     PetscInt     f, nf = 0, of = 0;
1994d9407bcSMatthew G. Knepley 
2004d9407bcSMatthew G. Knepley     ierr = PetscSectionCreateSubsection(section, numFields, fields, &subsection);CHKERRQ(ierr);
20192fd8e1eSJed Brown     ierr = DMSetLocalSection(*subdm, subsection);CHKERRQ(ierr);
2024d9407bcSMatthew G. Knepley     ierr = PetscSectionDestroy(&subsection);CHKERRQ(ierr);
2034d9407bcSMatthew G. Knepley     for (f = 0; f < numFields; ++f) {
2044d9407bcSMatthew G. Knepley       (*subdm)->nullspaceConstructors[f] = dm->nullspaceConstructors[fields[f]];
2054d9407bcSMatthew G. Knepley       if ((*subdm)->nullspaceConstructors[f]) {
2064d9407bcSMatthew G. Knepley         haveNull = PETSC_TRUE;
2074d9407bcSMatthew G. Knepley         nf       = f;
2088cda7954SMatthew G. Knepley         of       = fields[f];
2094d9407bcSMatthew G. Knepley       }
2104d9407bcSMatthew G. Knepley     }
211e5e52638SMatthew G. Knepley     if (dm->probs) {
2124d9407bcSMatthew G. Knepley       ierr = DMSetNumFields(*subdm, numFields);CHKERRQ(ierr);
2134d9407bcSMatthew G. Knepley       for (f = 0; f < numFields; ++f) {
2140f21e855SMatthew G. Knepley         PetscObject disc;
2150f21e855SMatthew G. Knepley 
21644a7f3ddSMatthew G. Knepley         ierr = DMGetField(dm, fields[f], NULL, &disc);CHKERRQ(ierr);
21744a7f3ddSMatthew G. Knepley         ierr = DMSetField(*subdm, f, NULL, disc);CHKERRQ(ierr);
2184d9407bcSMatthew G. Knepley       }
219e5e52638SMatthew G. Knepley       ierr = DMCreateDS(*subdm);CHKERRQ(ierr);
220f646a522SMatthew G. Knepley       if (numFields == 1 && is) {
2210f21e855SMatthew G. Knepley         PetscObject disc, space, pmat;
2224d9407bcSMatthew G. Knepley 
22344a7f3ddSMatthew G. Knepley         ierr = DMGetField(*subdm, 0, NULL, &disc);CHKERRQ(ierr);
2240f21e855SMatthew G. Knepley         ierr = PetscObjectQuery(disc, "nullspace", &space);CHKERRQ(ierr);
2250f21e855SMatthew G. Knepley         if (space) {ierr = PetscObjectCompose((PetscObject) *is, "nullspace", space);CHKERRQ(ierr);}
2260f21e855SMatthew G. Knepley         ierr = PetscObjectQuery(disc, "nearnullspace", &space);CHKERRQ(ierr);
2270f21e855SMatthew G. Knepley         if (space) {ierr = PetscObjectCompose((PetscObject) *is, "nearnullspace", space);CHKERRQ(ierr);}
2280f21e855SMatthew G. Knepley         ierr = PetscObjectQuery(disc, "pmat", &pmat);CHKERRQ(ierr);
2290f21e855SMatthew G. Knepley         if (pmat) {ierr = PetscObjectCompose((PetscObject) *is, "pmat", pmat);CHKERRQ(ierr);}
2304d9407bcSMatthew G. Knepley       }
2319310035eSMatthew G. Knepley       /* Check if DSes record their DM fields */
232e21d936aSMatthew G. Knepley       if (dm->probs[0].fields) {
2339310035eSMatthew G. Knepley         PetscInt d, e;
2349310035eSMatthew G. Knepley 
2359310035eSMatthew G. Knepley         for (d = 0, e = 0; d < dm->Nds && e < (*subdm)->Nds; ++d) {
2369310035eSMatthew G. Knepley           const PetscInt  Nf = dm->probs[d].ds->Nf;
2379310035eSMatthew G. Knepley           const PetscInt *fld;
2389310035eSMatthew G. Knepley           PetscInt        f, g;
2399310035eSMatthew G. Knepley 
2409310035eSMatthew G. Knepley           ierr = ISGetIndices(dm->probs[d].fields, &fld);CHKERRQ(ierr);
2419310035eSMatthew G. Knepley           for (f = 0; f < Nf; ++f) {
2429310035eSMatthew G. Knepley             for (g = 0; g < numFields; ++g) if (fld[f] == fields[g]) break;
2439310035eSMatthew G. Knepley             if (g < numFields) break;
2449310035eSMatthew G. Knepley           }
2459310035eSMatthew G. Knepley           ierr = ISRestoreIndices(dm->probs[d].fields, &fld);CHKERRQ(ierr);
2469310035eSMatthew G. Knepley           if (f == Nf) continue;
2479310035eSMatthew G. Knepley           ierr = PetscDSCopyConstants(dm->probs[d].ds, (*subdm)->probs[e].ds);CHKERRQ(ierr);
24836951cb5SMatthew G. Knepley           ierr = PetscDSCopyBoundary(dm->probs[d].ds, numFields, fields, (*subdm)->probs[e].ds);CHKERRQ(ierr);
2499310035eSMatthew G. Knepley           /* Translate DM fields to DS fields */
2509310035eSMatthew G. Knepley           {
251e21d936aSMatthew G. Knepley             IS              infields, dsfields;
25274a61aaaSMatthew G. Knepley             const PetscInt *fld, *ofld;
25374a61aaaSMatthew G. Knepley             PetscInt       *fidx;
25474a61aaaSMatthew G. Knepley             PetscInt        onf, nf, f, g;
255e21d936aSMatthew G. Knepley 
256e21d936aSMatthew G. Knepley             ierr = ISCreateGeneral(PETSC_COMM_SELF, numFields, fields, PETSC_USE_POINTER, &infields);CHKERRQ(ierr);
2579310035eSMatthew G. Knepley             ierr = ISIntersect(infields, dm->probs[d].fields, &dsfields);CHKERRQ(ierr);
258e21d936aSMatthew G. Knepley             ierr = ISDestroy(&infields);CHKERRQ(ierr);
259e21d936aSMatthew G. Knepley             ierr = ISGetLocalSize(dsfields, &nf);CHKERRQ(ierr);
260*7a8be351SBarry Smith             PetscCheck(nf,PETSC_COMM_SELF, PETSC_ERR_PLIB, "DS cannot be supported on 0 fields");
261e21d936aSMatthew G. Knepley             ierr = ISGetIndices(dsfields, &fld);CHKERRQ(ierr);
2629310035eSMatthew G. Knepley             ierr = ISGetLocalSize(dm->probs[d].fields, &onf);CHKERRQ(ierr);
2639310035eSMatthew G. Knepley             ierr = ISGetIndices(dm->probs[d].fields, &ofld);CHKERRQ(ierr);
26474a61aaaSMatthew G. Knepley             ierr = PetscMalloc1(nf, &fidx);CHKERRQ(ierr);
2653268a0aaSMatthew G. Knepley             for (f = 0, g = 0; f < onf && g < nf; ++f) {
26674a61aaaSMatthew G. Knepley               if (ofld[f] == fld[g]) fidx[g++] = f;
26774a61aaaSMatthew G. Knepley             }
2689310035eSMatthew G. Knepley             ierr = ISRestoreIndices(dm->probs[d].fields, &ofld);CHKERRQ(ierr);
269e21d936aSMatthew G. Knepley             ierr = ISRestoreIndices(dsfields, &fld);CHKERRQ(ierr);
270e21d936aSMatthew G. Knepley             ierr = ISDestroy(&dsfields);CHKERRQ(ierr);
2716c1eb96dSMatthew G. Knepley             ierr = PetscDSSelectDiscretizations(dm->probs[0].ds, nf, fidx, (*subdm)->probs[0].ds);CHKERRQ(ierr);
27274a61aaaSMatthew G. Knepley             ierr = PetscDSSelectEquations(dm->probs[0].ds, nf, fidx, (*subdm)->probs[0].ds);CHKERRQ(ierr);
27374a61aaaSMatthew G. Knepley             ierr = PetscFree(fidx);CHKERRQ(ierr);
2749310035eSMatthew G. Knepley           }
2759310035eSMatthew G. Knepley           ++e;
2769310035eSMatthew G. Knepley         }
277e21d936aSMatthew G. Knepley       } else {
2789310035eSMatthew G. Knepley         ierr = PetscDSCopyConstants(dm->probs[0].ds, (*subdm)->probs[0].ds);CHKERRQ(ierr);
27936951cb5SMatthew G. Knepley         ierr = PetscDSCopyBoundary(dm->probs[0].ds, PETSC_DETERMINE, NULL, (*subdm)->probs[0].ds);CHKERRQ(ierr);
2806c1eb96dSMatthew G. Knepley         ierr = PetscDSSelectDiscretizations(dm->probs[0].ds, numFields, fields, (*subdm)->probs[0].ds);CHKERRQ(ierr);
281e5e52638SMatthew G. Knepley         ierr = PetscDSSelectEquations(dm->probs[0].ds, numFields, fields, (*subdm)->probs[0].ds);CHKERRQ(ierr);
282d5af271aSMatthew G. Knepley       }
283e21d936aSMatthew G. Knepley     }
2848cda7954SMatthew G. Knepley     if (haveNull && is) {
2858cda7954SMatthew G. Knepley       MatNullSpace nullSpace;
2868cda7954SMatthew G. Knepley 
2878cda7954SMatthew G. Knepley       ierr = (*(*subdm)->nullspaceConstructors[nf])(*subdm, of, nf, &nullSpace);CHKERRQ(ierr);
2888cda7954SMatthew G. Knepley       ierr = PetscObjectCompose((PetscObject) *is, "nullspace", (PetscObject) nullSpace);CHKERRQ(ierr);
2898cda7954SMatthew G. Knepley       ierr = MatNullSpaceDestroy(&nullSpace);CHKERRQ(ierr);
2908cda7954SMatthew G. Knepley     }
291d5af271aSMatthew G. Knepley     if (dm->coarseMesh) {
292d5af271aSMatthew G. Knepley       ierr = DMCreateSubDM(dm->coarseMesh, numFields, fields, NULL, &(*subdm)->coarseMesh);CHKERRQ(ierr);
2934d9407bcSMatthew G. Knepley     }
2944d9407bcSMatthew G. Knepley   }
2954d9407bcSMatthew G. Knepley   PetscFunctionReturn(0);
2964d9407bcSMatthew G. Knepley }
2972adcc780SMatthew G. Knepley 
298792b654fSMatthew G. Knepley /*@C
299792b654fSMatthew G. Knepley   DMCreateSectionSuperDM - Returns an arrays of ISes and DM+Section encapsulating a superproblem defined by the DM+Sections passed in.
300792b654fSMatthew G. Knepley 
301792b654fSMatthew G. Knepley   Not collective
302792b654fSMatthew G. Knepley 
303d8d19677SJose E. Roman   Input Parameters:
304792b654fSMatthew G. Knepley + dms - The DM objects
305792b654fSMatthew G. Knepley - len - The number of DMs
306792b654fSMatthew G. Knepley 
307792b654fSMatthew G. Knepley   Output Parameters:
308792b654fSMatthew G. Knepley + is - The global indices for the subproblem, or NULL
309792b654fSMatthew G. Knepley - superdm - The DM for the superproblem, which must already have be cloned
310792b654fSMatthew G. Knepley 
311792b654fSMatthew 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,
312792b654fSMatthew G. Knepley   such as Plex and Forest.
313792b654fSMatthew G. Knepley 
314792b654fSMatthew G. Knepley   Level: intermediate
315792b654fSMatthew G. Knepley 
31692fd8e1eSJed Brown .seealso DMCreateSuperDM(), DMGetLocalSection(), DMPlexSetMigrationSF(), DMView()
317792b654fSMatthew G. Knepley @*/
318792b654fSMatthew G. Knepley PetscErrorCode DMCreateSectionSuperDM(DM dms[], PetscInt len, IS **is, DM *superdm)
3192adcc780SMatthew G. Knepley {
3209796d432SMatthew G. Knepley   MPI_Comm       comm;
3219796d432SMatthew G. Knepley   PetscSection   supersection, *sections, *sectionGlobals;
3228cda7954SMatthew G. Knepley   PetscInt      *Nfs, Nf = 0, f, supf, oldf = -1, nullf = -1, i;
3239796d432SMatthew G. Knepley   PetscBool      haveNull = PETSC_FALSE;
3242adcc780SMatthew G. Knepley   PetscErrorCode ierr;
3252adcc780SMatthew G. Knepley 
3262adcc780SMatthew G. Knepley   PetscFunctionBegin;
3279796d432SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)dms[0], &comm);CHKERRQ(ierr);
3289796d432SMatthew G. Knepley   /* Pull out local and global sections */
3292adcc780SMatthew G. Knepley   ierr = PetscMalloc3(len, &Nfs, len, &sections, len, &sectionGlobals);CHKERRQ(ierr);
3302adcc780SMatthew G. Knepley   for (i = 0 ; i < len; ++i) {
33192fd8e1eSJed Brown     ierr = DMGetLocalSection(dms[i], &sections[i]);CHKERRQ(ierr);
332e87a4003SBarry Smith     ierr = DMGetGlobalSection(dms[i], &sectionGlobals[i]);CHKERRQ(ierr);
333*7a8be351SBarry Smith     PetscCheck(sections[i],comm, PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields");
334*7a8be351SBarry Smith     PetscCheck(sectionGlobals[i],comm, PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields");
3352adcc780SMatthew G. Knepley     ierr = PetscSectionGetNumFields(sections[i], &Nfs[i]);CHKERRQ(ierr);
3362adcc780SMatthew G. Knepley     Nf += Nfs[i];
3372adcc780SMatthew G. Knepley   }
3389796d432SMatthew G. Knepley   /* Create the supersection */
3399796d432SMatthew G. Knepley   ierr = PetscSectionCreateSupersection(sections, len, &supersection);CHKERRQ(ierr);
34092fd8e1eSJed Brown   ierr = DMSetLocalSection(*superdm, supersection);CHKERRQ(ierr);
3419796d432SMatthew G. Knepley   /* Create ISes */
3422adcc780SMatthew G. Knepley   if (is) {
3439796d432SMatthew G. Knepley     PetscSection supersectionGlobal;
3449796d432SMatthew G. Knepley     PetscInt     bs = -1, startf = 0;
3452adcc780SMatthew G. Knepley 
3462adcc780SMatthew G. Knepley     ierr = PetscMalloc1(len, is);CHKERRQ(ierr);
3476f0eb057SJed Brown     ierr = DMGetGlobalSection(*superdm, &supersectionGlobal);CHKERRQ(ierr);
3489796d432SMatthew G. Knepley     for (i = 0 ; i < len; startf += Nfs[i], ++i) {
3499796d432SMatthew G. Knepley       PetscInt *subIndices;
350ec4c761aSMatthew G. Knepley       PetscInt  subSize, subOff, pStart, pEnd, p, start, end, dummy;
3512adcc780SMatthew G. Knepley 
3522adcc780SMatthew G. Knepley       ierr = PetscSectionGetChart(sectionGlobals[i], &pStart, &pEnd);CHKERRQ(ierr);
3532adcc780SMatthew G. Knepley       ierr = PetscSectionGetConstrainedStorageSize(sectionGlobals[i], &subSize);CHKERRQ(ierr);
3542adcc780SMatthew G. Knepley       ierr = PetscMalloc1(subSize, &subIndices);CHKERRQ(ierr);
3559796d432SMatthew G. Knepley       for (p = pStart, subOff = 0; p < pEnd; ++p) {
3569796d432SMatthew G. Knepley         PetscInt gdof, gcdof, gtdof, d;
3572adcc780SMatthew G. Knepley 
3582adcc780SMatthew G. Knepley         ierr = PetscSectionGetDof(sectionGlobals[i], p, &gdof);CHKERRQ(ierr);
3592adcc780SMatthew G. Knepley         ierr = PetscSectionGetConstraintDof(sections[i], p, &gcdof);CHKERRQ(ierr);
3602adcc780SMatthew G. Knepley         gtdof = gdof - gcdof;
3612adcc780SMatthew G. Knepley         if (gdof > 0 && gtdof) {
3622adcc780SMatthew G. Knepley           if (bs < 0)           {bs = gtdof;}
3632adcc780SMatthew G. Knepley           else if (bs != gtdof) {bs = 1;}
364ec4c761aSMatthew G. Knepley           ierr = DMGetGlobalFieldOffset_Private(*superdm, p, startf, &start, &dummy);CHKERRQ(ierr);
365ec4c761aSMatthew G. Knepley           ierr = DMGetGlobalFieldOffset_Private(*superdm, p, startf+Nfs[i]-1, &dummy, &end);CHKERRQ(ierr);
366*7a8be351SBarry Smith           PetscCheck(end-start == gtdof,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);
3679796d432SMatthew G. Knepley           for (d = start; d < end; ++d, ++subOff) subIndices[subOff] = d;
3682adcc780SMatthew G. Knepley         }
3692adcc780SMatthew G. Knepley       }
3709796d432SMatthew G. Knepley       ierr = ISCreateGeneral(comm, subSize, subIndices, PETSC_OWN_POINTER, &(*is)[i]);CHKERRQ(ierr);
3712adcc780SMatthew G. Knepley       /* Must have same blocksize on all procs (some might have no points) */
3729796d432SMatthew G. Knepley       {
3739796d432SMatthew G. Knepley         PetscInt bs = -1, bsLocal[2], bsMinMax[2];
3749796d432SMatthew G. Knepley 
3752adcc780SMatthew G. Knepley         bsLocal[0] = bs < 0 ? PETSC_MAX_INT : bs; bsLocal[1] = bs;
3769796d432SMatthew G. Knepley         ierr = PetscGlobalMinMaxInt(comm, bsLocal, bsMinMax);CHKERRQ(ierr);
3772adcc780SMatthew G. Knepley         if (bsMinMax[0] != bsMinMax[1]) {bs = 1;}
3782adcc780SMatthew G. Knepley         else                            {bs = bsMinMax[0];}
3792adcc780SMatthew G. Knepley         ierr = ISSetBlockSize((*is)[i], bs);CHKERRQ(ierr);
3802adcc780SMatthew G. Knepley       }
3812adcc780SMatthew G. Knepley     }
3822adcc780SMatthew G. Knepley   }
3839796d432SMatthew G. Knepley   /* Preserve discretizations */
384e5e52638SMatthew G. Knepley   if (len && dms[0]->probs) {
3852adcc780SMatthew G. Knepley     ierr = DMSetNumFields(*superdm, Nf);CHKERRQ(ierr);
3869796d432SMatthew G. Knepley     for (i = 0, supf = 0; i < len; ++i) {
3879796d432SMatthew G. Knepley       for (f = 0; f < Nfs[i]; ++f, ++supf) {
3882adcc780SMatthew G. Knepley         PetscObject disc;
3892adcc780SMatthew G. Knepley 
39044a7f3ddSMatthew G. Knepley         ierr = DMGetField(dms[i], f, NULL, &disc);CHKERRQ(ierr);
39144a7f3ddSMatthew G. Knepley         ierr = DMSetField(*superdm, supf, NULL, disc);CHKERRQ(ierr);
3922adcc780SMatthew G. Knepley       }
3932adcc780SMatthew G. Knepley     }
394e5e52638SMatthew G. Knepley     ierr = DMCreateDS(*superdm);CHKERRQ(ierr);
3952adcc780SMatthew G. Knepley   }
3969796d432SMatthew G. Knepley   /* Preserve nullspaces */
3979796d432SMatthew G. Knepley   for (i = 0, supf = 0; i < len; ++i) {
3989796d432SMatthew G. Knepley     for (f = 0; f < Nfs[i]; ++f, ++supf) {
3999796d432SMatthew G. Knepley       (*superdm)->nullspaceConstructors[supf] = dms[i]->nullspaceConstructors[f];
4009796d432SMatthew G. Knepley       if ((*superdm)->nullspaceConstructors[supf]) {
4019796d432SMatthew G. Knepley         haveNull = PETSC_TRUE;
4029796d432SMatthew G. Knepley         nullf    = supf;
4038cda7954SMatthew G. Knepley         oldf     = f;
4042adcc780SMatthew G. Knepley       }
4059796d432SMatthew G. Knepley     }
4069796d432SMatthew G. Knepley   }
4079796d432SMatthew G. Knepley   /* Attach nullspace to IS */
4089796d432SMatthew G. Knepley   if (haveNull && is) {
4099796d432SMatthew G. Knepley     MatNullSpace nullSpace;
4109796d432SMatthew G. Knepley 
4118cda7954SMatthew G. Knepley     ierr = (*(*superdm)->nullspaceConstructors[nullf])(*superdm, oldf, nullf, &nullSpace);CHKERRQ(ierr);
4129796d432SMatthew G. Knepley     ierr = PetscObjectCompose((PetscObject) (*is)[nullf], "nullspace", (PetscObject) nullSpace);CHKERRQ(ierr);
4139796d432SMatthew G. Knepley     ierr = MatNullSpaceDestroy(&nullSpace);CHKERRQ(ierr);
4149796d432SMatthew G. Knepley   }
4159796d432SMatthew G. Knepley   ierr = PetscSectionDestroy(&supersection);CHKERRQ(ierr);
41695b1d6b7SMatthew G. Knepley   ierr = PetscFree3(Nfs, sections, sectionGlobals);CHKERRQ(ierr);
4172adcc780SMatthew G. Knepley   PetscFunctionReturn(0);
4182adcc780SMatthew G. Knepley }
419