xref: /petsc/src/dm/interface/dmi.c (revision 63a3b9bc7a1f24f247904ccba9383635fe6abade)
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;
185d1c0279SHong Zhang   PetscInt       in[2],out[2];
1911689aebSJed Brown 
2011689aebSJed Brown   PetscFunctionBegin;
219566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalSection(dm, &gSection));
229566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(gSection, &pStart, &pEnd));
2311689aebSJed Brown   for (p = pStart; p < pEnd; ++p) {
2411689aebSJed Brown     PetscInt dof, cdof;
2511689aebSJed Brown 
269566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(gSection, p, &dof));
279566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(gSection, p, &cdof));
280680ce0aSHong Zhang 
295b8243e1SJed Brown     if (dof - cdof > 0) {
305b8243e1SJed Brown       if (blockSize < 0) {
310680ce0aSHong Zhang         /* set blockSize */
320680ce0aSHong Zhang         blockSize = dof-cdof;
335b8243e1SJed Brown       } else {
345b8243e1SJed Brown         blockSize = PetscGCD(dof - cdof, blockSize);
3511689aebSJed Brown       }
3611689aebSJed Brown     }
370680ce0aSHong Zhang   }
385d1c0279SHong Zhang 
392c8be454SMatthew G. Knepley   in[0] = blockSize < 0 ? PETSC_MIN_INT : -blockSize;
405d1c0279SHong Zhang   in[1] = blockSize;
411c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(in,out,2,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)dm)));
425d1c0279SHong Zhang   /* -out[0] = min(blockSize), out[1] = max(blockSize) */
435d1c0279SHong Zhang   if (-out[0] == out[1]) {
445d1c0279SHong Zhang     bs = out[1];
455d1c0279SHong Zhang   } else bs = 1;
465d1c0279SHong Zhang 
475d1c0279SHong Zhang   if (bs < 0) { /* Everyone was empty */
485d1c0279SHong Zhang     blockSize = 1;
495d1c0279SHong Zhang     bs        = 1;
505d1c0279SHong Zhang   }
515d1c0279SHong Zhang 
529566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetConstrainedStorageSize(gSection, &localSize));
537a8be351SBarry 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);
549566063dSJacob Faibussowitsch   PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), vec));
559566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(*vec, localSize, PETSC_DETERMINE));
569566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(*vec, bs));
579566063dSJacob Faibussowitsch   PetscCall(VecSetType(*vec,dm->vectype));
589566063dSJacob Faibussowitsch   PetscCall(VecSetDM(*vec, dm));
599566063dSJacob Faibussowitsch   /* PetscCall(VecSetLocalToGlobalMapping(*vec, dm->ltogmap)); */
6011689aebSJed Brown   PetscFunctionReturn(0);
6111689aebSJed Brown }
6211689aebSJed Brown 
6311689aebSJed Brown PetscErrorCode DMCreateLocalVector_Section_Private(DM dm,Vec *vec)
6411689aebSJed Brown {
65fad22124SMatthew G Knepley   PetscSection   section;
6611689aebSJed Brown   PetscInt       localSize, blockSize = -1, pStart, pEnd, p;
6711689aebSJed Brown 
6811689aebSJed Brown   PetscFunctionBegin;
699566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
709566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
7111689aebSJed Brown   for (p = pStart; p < pEnd; ++p) {
7211689aebSJed Brown     PetscInt dof;
7311689aebSJed Brown 
749566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
7511689aebSJed Brown     if ((blockSize < 0) && (dof > 0)) blockSize = dof;
765b8243e1SJed Brown     if (dof > 0) blockSize = PetscGCD(dof, blockSize);
7711689aebSJed Brown   }
789566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(section, &localSize));
799566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF, vec));
809566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(*vec, localSize, localSize));
819566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(*vec, blockSize));
829566063dSJacob Faibussowitsch   PetscCall(VecSetType(*vec,dm->vectype));
839566063dSJacob Faibussowitsch   PetscCall(VecSetDM(*vec, dm));
8411689aebSJed Brown   PetscFunctionReturn(0);
8511689aebSJed Brown }
864d9407bcSMatthew G. Knepley 
87792b654fSMatthew G. Knepley /*@C
88792b654fSMatthew G. Knepley   DMCreateSectionSubDM - Returns an IS and subDM+subSection encapsulating a subproblem defined by the fields in a PetscSection in the DM.
89792b654fSMatthew G. Knepley 
90792b654fSMatthew G. Knepley   Not collective
91792b654fSMatthew G. Knepley 
92792b654fSMatthew G. Knepley   Input Parameters:
93792b654fSMatthew G. Knepley + dm        - The DM object
94792b654fSMatthew G. Knepley . numFields - The number of fields in this subproblem
95792b654fSMatthew G. Knepley - fields    - The field numbers of the selected fields
96792b654fSMatthew G. Knepley 
97792b654fSMatthew G. Knepley   Output Parameters:
98792b654fSMatthew G. Knepley + is - The global indices for the subproblem
99792b654fSMatthew G. Knepley - subdm - The DM for the subproblem, which must already have be cloned from dm
100792b654fSMatthew G. Knepley 
101792b654fSMatthew 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,
102792b654fSMatthew G. Knepley   such as Plex and Forest.
103792b654fSMatthew G. Knepley 
104792b654fSMatthew G. Knepley   Level: intermediate
105792b654fSMatthew G. Knepley 
10692fd8e1eSJed Brown .seealso DMCreateSubDM(), DMGetLocalSection(), DMPlexSetMigrationSF(), DMView()
107792b654fSMatthew G. Knepley @*/
108792b654fSMatthew G. Knepley PetscErrorCode DMCreateSectionSubDM(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm)
1094d9407bcSMatthew G. Knepley {
1104d9407bcSMatthew G. Knepley   PetscSection   section, sectionGlobal;
1114d9407bcSMatthew G. Knepley   PetscInt      *subIndices;
112696188eaSMatthew G. Knepley   PetscInt       subSize = 0, subOff = 0, Nf, f, pStart, pEnd, p;
1134d9407bcSMatthew G. Knepley 
1144d9407bcSMatthew G. Knepley   PetscFunctionBegin;
1154d9407bcSMatthew G. Knepley   if (!numFields) PetscFunctionReturn(0);
1169566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
1179566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalSection(dm, &sectionGlobal));
1187a8be351SBarry Smith   PetscCheck(section,PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields");
1197a8be351SBarry Smith   PetscCheck(sectionGlobal,PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields");
1209566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetNumFields(section, &Nf));
1217a8be351SBarry 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);
1224d9407bcSMatthew G. Knepley   if (is) {
1233a544194SStefano Zampini     PetscInt bs, bsLocal[2], bsMinMax[2];
1243a544194SStefano Zampini 
1253a544194SStefano Zampini     for (f = 0, bs = 0; f < numFields; ++f) {
1263a544194SStefano Zampini       PetscInt Nc;
1273a544194SStefano Zampini 
1289566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldComponents(section, fields[f], &Nc));
1293a544194SStefano Zampini       bs  += Nc;
1303a544194SStefano Zampini     }
1319566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetChart(sectionGlobal, &pStart, &pEnd));
1324d9407bcSMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
133b9eca265Ssarens       PetscInt gdof, pSubSize  = 0;
1344d9407bcSMatthew G. Knepley 
1359566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(sectionGlobal, p, &gdof));
1364d9407bcSMatthew G. Knepley       if (gdof > 0) {
1374d9407bcSMatthew G. Knepley         for (f = 0; f < numFields; ++f) {
1384d9407bcSMatthew G. Knepley           PetscInt fdof, fcdof;
1394d9407bcSMatthew G. Knepley 
1409566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetFieldDof(section, p, fields[f], &fdof));
1419566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetFieldConstraintDof(section, p, fields[f], &fcdof));
142b9eca265Ssarens           pSubSize += fdof-fcdof;
143b9eca265Ssarens         }
144b9eca265Ssarens         subSize += pSubSize;
1453a544194SStefano Zampini         if (pSubSize && bs != pSubSize) {
146b9eca265Ssarens           /* Layout does not admit a pointwise block size */
147b9eca265Ssarens           bs = 1;
1484d9407bcSMatthew G. Knepley         }
1494d9407bcSMatthew G. Knepley       }
1504d9407bcSMatthew G. Knepley     }
151b9eca265Ssarens     /* Must have same blocksize on all procs (some might have no points) */
1520be3e97aSMatthew G. Knepley     bsLocal[0] = bs < 0 ? PETSC_MAX_INT : bs; bsLocal[1] = bs;
1539566063dSJacob Faibussowitsch     PetscCall(PetscGlobalMinMaxInt(PetscObjectComm((PetscObject) dm), bsLocal, bsMinMax));
1540be3e97aSMatthew G. Knepley     if (bsMinMax[0] != bsMinMax[1]) {bs = 1;}
1550be3e97aSMatthew G. Knepley     else                            {bs = bsMinMax[0];}
1569566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(subSize, &subIndices));
1574d9407bcSMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
1584d9407bcSMatthew G. Knepley       PetscInt gdof, goff;
1594d9407bcSMatthew G. Knepley 
1609566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(sectionGlobal, p, &gdof));
1614d9407bcSMatthew G. Knepley       if (gdof > 0) {
1629566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetOffset(sectionGlobal, p, &goff));
1634d9407bcSMatthew G. Knepley         for (f = 0; f < numFields; ++f) {
1644d9407bcSMatthew G. Knepley           PetscInt fdof, fcdof, fc, f2, poff = 0;
1654d9407bcSMatthew G. Knepley 
1664d9407bcSMatthew G. Knepley           /* Can get rid of this loop by storing field information in the global section */
1674d9407bcSMatthew G. Knepley           for (f2 = 0; f2 < fields[f]; ++f2) {
1689566063dSJacob Faibussowitsch             PetscCall(PetscSectionGetFieldDof(section, p, f2, &fdof));
1699566063dSJacob Faibussowitsch             PetscCall(PetscSectionGetFieldConstraintDof(section, p, f2, &fcdof));
1704d9407bcSMatthew G. Knepley             poff += fdof-fcdof;
1714d9407bcSMatthew G. Knepley           }
1729566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetFieldDof(section, p, fields[f], &fdof));
1739566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetFieldConstraintDof(section, p, fields[f], &fcdof));
1744d9407bcSMatthew G. Knepley           for (fc = 0; fc < fdof-fcdof; ++fc, ++subOff) {
1754d9407bcSMatthew G. Knepley             subIndices[subOff] = goff+poff+fc;
1764d9407bcSMatthew G. Knepley           }
1774d9407bcSMatthew G. Knepley         }
1784d9407bcSMatthew G. Knepley       }
1794d9407bcSMatthew G. Knepley     }
1809566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), subSize, subIndices, PETSC_OWN_POINTER, is));
181627349f5SMatthew G. Knepley     if (bs > 1) {
182627349f5SMatthew G. Knepley       /* We need to check that the block size does not come from non-contiguous fields */
183627349f5SMatthew G. Knepley       PetscInt i, j, set = 1;
184627349f5SMatthew G. Knepley       for (i = 0; i < subSize; i += bs) {
185627349f5SMatthew G. Knepley         for (j = 0; j < bs; ++j) {
186627349f5SMatthew G. Knepley           if (subIndices[i+j] != subIndices[i]+j) {set = 0; break;}
187627349f5SMatthew G. Knepley         }
188627349f5SMatthew G. Knepley       }
1899566063dSJacob Faibussowitsch       if (set) PetscCall(ISSetBlockSize(*is, bs));
190627349f5SMatthew G. Knepley     }
1914d9407bcSMatthew G. Knepley   }
1924d9407bcSMatthew G. Knepley   if (subdm) {
1934d9407bcSMatthew G. Knepley     PetscSection subsection;
1944d9407bcSMatthew G. Knepley     PetscBool    haveNull = PETSC_FALSE;
1958cda7954SMatthew G. Knepley     PetscInt     f, nf = 0, of = 0;
1964d9407bcSMatthew G. Knepley 
1979566063dSJacob Faibussowitsch     PetscCall(PetscSectionCreateSubsection(section, numFields, fields, &subsection));
1989566063dSJacob Faibussowitsch     PetscCall(DMSetLocalSection(*subdm, subsection));
1999566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&subsection));
2004d9407bcSMatthew G. Knepley     for (f = 0; f < numFields; ++f) {
2014d9407bcSMatthew G. Knepley       (*subdm)->nullspaceConstructors[f] = dm->nullspaceConstructors[fields[f]];
2024d9407bcSMatthew G. Knepley       if ((*subdm)->nullspaceConstructors[f]) {
2034d9407bcSMatthew G. Knepley         haveNull = PETSC_TRUE;
2044d9407bcSMatthew G. Knepley         nf       = f;
2058cda7954SMatthew G. Knepley         of       = fields[f];
2064d9407bcSMatthew G. Knepley       }
2074d9407bcSMatthew G. Knepley     }
208e5e52638SMatthew G. Knepley     if (dm->probs) {
2099566063dSJacob Faibussowitsch       PetscCall(DMSetNumFields(*subdm, numFields));
2104d9407bcSMatthew G. Knepley       for (f = 0; f < numFields; ++f) {
2110f21e855SMatthew G. Knepley         PetscObject disc;
2120f21e855SMatthew G. Knepley 
2139566063dSJacob Faibussowitsch         PetscCall(DMGetField(dm, fields[f], NULL, &disc));
2149566063dSJacob Faibussowitsch         PetscCall(DMSetField(*subdm, f, NULL, disc));
2154d9407bcSMatthew G. Knepley       }
2169566063dSJacob Faibussowitsch       PetscCall(DMCreateDS(*subdm));
217f646a522SMatthew G. Knepley       if (numFields == 1 && is) {
2180f21e855SMatthew G. Knepley         PetscObject disc, space, pmat;
2194d9407bcSMatthew G. Knepley 
2209566063dSJacob Faibussowitsch         PetscCall(DMGetField(*subdm, 0, NULL, &disc));
2219566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "nullspace", &space));
2229566063dSJacob Faibussowitsch         if (space) PetscCall(PetscObjectCompose((PetscObject) *is, "nullspace", space));
2239566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "nearnullspace", &space));
2249566063dSJacob Faibussowitsch         if (space) PetscCall(PetscObjectCompose((PetscObject) *is, "nearnullspace", space));
2259566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "pmat", &pmat));
2269566063dSJacob Faibussowitsch         if (pmat) PetscCall(PetscObjectCompose((PetscObject) *is, "pmat", pmat));
2274d9407bcSMatthew G. Knepley       }
2289310035eSMatthew G. Knepley       /* Check if DSes record their DM fields */
229e21d936aSMatthew G. Knepley       if (dm->probs[0].fields) {
2309310035eSMatthew G. Knepley         PetscInt d, e;
2319310035eSMatthew G. Knepley 
2329310035eSMatthew G. Knepley         for (d = 0, e = 0; d < dm->Nds && e < (*subdm)->Nds; ++d) {
2339310035eSMatthew G. Knepley           const PetscInt  Nf = dm->probs[d].ds->Nf;
2349310035eSMatthew G. Knepley           const PetscInt *fld;
2359310035eSMatthew G. Knepley           PetscInt        f, g;
2369310035eSMatthew G. Knepley 
2379566063dSJacob Faibussowitsch           PetscCall(ISGetIndices(dm->probs[d].fields, &fld));
2389310035eSMatthew G. Knepley           for (f = 0; f < Nf; ++f) {
2399310035eSMatthew G. Knepley             for (g = 0; g < numFields; ++g) if (fld[f] == fields[g]) break;
2409310035eSMatthew G. Knepley             if (g < numFields) break;
2419310035eSMatthew G. Knepley           }
2429566063dSJacob Faibussowitsch           PetscCall(ISRestoreIndices(dm->probs[d].fields, &fld));
2439310035eSMatthew G. Knepley           if (f == Nf) continue;
2449566063dSJacob Faibussowitsch           PetscCall(PetscDSCopyConstants(dm->probs[d].ds, (*subdm)->probs[e].ds));
2459566063dSJacob Faibussowitsch           PetscCall(PetscDSCopyBoundary(dm->probs[d].ds, numFields, fields, (*subdm)->probs[e].ds));
2469310035eSMatthew G. Knepley           /* Translate DM fields to DS fields */
2479310035eSMatthew G. Knepley           {
248e21d936aSMatthew G. Knepley             IS              infields, dsfields;
24974a61aaaSMatthew G. Knepley             const PetscInt *fld, *ofld;
25074a61aaaSMatthew G. Knepley             PetscInt       *fidx;
25174a61aaaSMatthew G. Knepley             PetscInt        onf, nf, f, g;
252e21d936aSMatthew G. Knepley 
2539566063dSJacob Faibussowitsch             PetscCall(ISCreateGeneral(PETSC_COMM_SELF, numFields, fields, PETSC_USE_POINTER, &infields));
2549566063dSJacob Faibussowitsch             PetscCall(ISIntersect(infields, dm->probs[d].fields, &dsfields));
2559566063dSJacob Faibussowitsch             PetscCall(ISDestroy(&infields));
2569566063dSJacob Faibussowitsch             PetscCall(ISGetLocalSize(dsfields, &nf));
2577a8be351SBarry Smith             PetscCheck(nf,PETSC_COMM_SELF, PETSC_ERR_PLIB, "DS cannot be supported on 0 fields");
2589566063dSJacob Faibussowitsch             PetscCall(ISGetIndices(dsfields, &fld));
2599566063dSJacob Faibussowitsch             PetscCall(ISGetLocalSize(dm->probs[d].fields, &onf));
2609566063dSJacob Faibussowitsch             PetscCall(ISGetIndices(dm->probs[d].fields, &ofld));
2619566063dSJacob Faibussowitsch             PetscCall(PetscMalloc1(nf, &fidx));
2623268a0aaSMatthew G. Knepley             for (f = 0, g = 0; f < onf && g < nf; ++f) {
26374a61aaaSMatthew G. Knepley               if (ofld[f] == fld[g]) fidx[g++] = f;
26474a61aaaSMatthew G. Knepley             }
2659566063dSJacob Faibussowitsch             PetscCall(ISRestoreIndices(dm->probs[d].fields, &ofld));
2669566063dSJacob Faibussowitsch             PetscCall(ISRestoreIndices(dsfields, &fld));
2679566063dSJacob Faibussowitsch             PetscCall(ISDestroy(&dsfields));
2689566063dSJacob Faibussowitsch             PetscCall(PetscDSSelectDiscretizations(dm->probs[0].ds, nf, fidx, (*subdm)->probs[0].ds));
2699566063dSJacob Faibussowitsch             PetscCall(PetscDSSelectEquations(dm->probs[0].ds, nf, fidx, (*subdm)->probs[0].ds));
2709566063dSJacob Faibussowitsch             PetscCall(PetscFree(fidx));
2719310035eSMatthew G. Knepley           }
2729310035eSMatthew G. Knepley           ++e;
2739310035eSMatthew G. Knepley         }
274e21d936aSMatthew G. Knepley       } else {
2759566063dSJacob Faibussowitsch         PetscCall(PetscDSCopyConstants(dm->probs[0].ds, (*subdm)->probs[0].ds));
2769566063dSJacob Faibussowitsch         PetscCall(PetscDSCopyBoundary(dm->probs[0].ds, PETSC_DETERMINE, NULL, (*subdm)->probs[0].ds));
2779566063dSJacob Faibussowitsch         PetscCall(PetscDSSelectDiscretizations(dm->probs[0].ds, numFields, fields, (*subdm)->probs[0].ds));
2789566063dSJacob Faibussowitsch         PetscCall(PetscDSSelectEquations(dm->probs[0].ds, numFields, fields, (*subdm)->probs[0].ds));
279d5af271aSMatthew G. Knepley       }
280e21d936aSMatthew G. Knepley     }
2818cda7954SMatthew G. Knepley     if (haveNull && is) {
2828cda7954SMatthew G. Knepley       MatNullSpace nullSpace;
2838cda7954SMatthew G. Knepley 
2849566063dSJacob Faibussowitsch       PetscCall((*(*subdm)->nullspaceConstructors[nf])(*subdm, of, nf, &nullSpace));
2859566063dSJacob Faibussowitsch       PetscCall(PetscObjectCompose((PetscObject) *is, "nullspace", (PetscObject) nullSpace));
2869566063dSJacob Faibussowitsch       PetscCall(MatNullSpaceDestroy(&nullSpace));
2878cda7954SMatthew G. Knepley     }
288d5af271aSMatthew G. Knepley     if (dm->coarseMesh) {
2899566063dSJacob Faibussowitsch       PetscCall(DMCreateSubDM(dm->coarseMesh, numFields, fields, NULL, &(*subdm)->coarseMesh));
2904d9407bcSMatthew G. Knepley     }
2914d9407bcSMatthew G. Knepley   }
2924d9407bcSMatthew G. Knepley   PetscFunctionReturn(0);
2934d9407bcSMatthew G. Knepley }
2942adcc780SMatthew G. Knepley 
295792b654fSMatthew G. Knepley /*@C
296792b654fSMatthew G. Knepley   DMCreateSectionSuperDM - Returns an arrays of ISes and DM+Section encapsulating a superproblem defined by the DM+Sections passed in.
297792b654fSMatthew G. Knepley 
298792b654fSMatthew G. Knepley   Not collective
299792b654fSMatthew G. Knepley 
300d8d19677SJose E. Roman   Input Parameters:
301792b654fSMatthew G. Knepley + dms - The DM objects
302792b654fSMatthew G. Knepley - len - The number of DMs
303792b654fSMatthew G. Knepley 
304792b654fSMatthew G. Knepley   Output Parameters:
305792b654fSMatthew G. Knepley + is - The global indices for the subproblem, or NULL
306792b654fSMatthew G. Knepley - superdm - The DM for the superproblem, which must already have be cloned
307792b654fSMatthew G. Knepley 
308792b654fSMatthew 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,
309792b654fSMatthew G. Knepley   such as Plex and Forest.
310792b654fSMatthew G. Knepley 
311792b654fSMatthew G. Knepley   Level: intermediate
312792b654fSMatthew G. Knepley 
31392fd8e1eSJed Brown .seealso DMCreateSuperDM(), DMGetLocalSection(), DMPlexSetMigrationSF(), DMView()
314792b654fSMatthew G. Knepley @*/
315792b654fSMatthew G. Knepley PetscErrorCode DMCreateSectionSuperDM(DM dms[], PetscInt len, IS **is, DM *superdm)
3162adcc780SMatthew G. Knepley {
3179796d432SMatthew G. Knepley   MPI_Comm       comm;
3189796d432SMatthew G. Knepley   PetscSection   supersection, *sections, *sectionGlobals;
3198cda7954SMatthew G. Knepley   PetscInt      *Nfs, Nf = 0, f, supf, oldf = -1, nullf = -1, i;
3209796d432SMatthew G. Knepley   PetscBool      haveNull = PETSC_FALSE;
3212adcc780SMatthew G. Knepley 
3222adcc780SMatthew G. Knepley   PetscFunctionBegin;
3239566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dms[0], &comm));
3249796d432SMatthew G. Knepley   /* Pull out local and global sections */
3259566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(len, &Nfs, len, &sections, len, &sectionGlobals));
3262adcc780SMatthew G. Knepley   for (i = 0 ; i < len; ++i) {
3279566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dms[i], &sections[i]));
3289566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalSection(dms[i], &sectionGlobals[i]));
3297a8be351SBarry Smith     PetscCheck(sections[i],comm, PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields");
3307a8be351SBarry Smith     PetscCheck(sectionGlobals[i],comm, PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields");
3319566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetNumFields(sections[i], &Nfs[i]));
3322adcc780SMatthew G. Knepley     Nf += Nfs[i];
3332adcc780SMatthew G. Knepley   }
3349796d432SMatthew G. Knepley   /* Create the supersection */
3359566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreateSupersection(sections, len, &supersection));
3369566063dSJacob Faibussowitsch   PetscCall(DMSetLocalSection(*superdm, supersection));
3379796d432SMatthew G. Knepley   /* Create ISes */
3382adcc780SMatthew G. Knepley   if (is) {
3399796d432SMatthew G. Knepley     PetscSection supersectionGlobal;
3409796d432SMatthew G. Knepley     PetscInt     bs = -1, startf = 0;
3412adcc780SMatthew G. Knepley 
3429566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(len, is));
3439566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalSection(*superdm, &supersectionGlobal));
3449796d432SMatthew G. Knepley     for (i = 0 ; i < len; startf += Nfs[i], ++i) {
3459796d432SMatthew G. Knepley       PetscInt *subIndices;
346ec4c761aSMatthew G. Knepley       PetscInt  subSize, subOff, pStart, pEnd, p, start, end, dummy;
3472adcc780SMatthew G. Knepley 
3489566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetChart(sectionGlobals[i], &pStart, &pEnd));
3499566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetConstrainedStorageSize(sectionGlobals[i], &subSize));
3509566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(subSize, &subIndices));
3519796d432SMatthew G. Knepley       for (p = pStart, subOff = 0; p < pEnd; ++p) {
3529796d432SMatthew G. Knepley         PetscInt gdof, gcdof, gtdof, d;
3532adcc780SMatthew G. Knepley 
3549566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(sectionGlobals[i], p, &gdof));
3559566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetConstraintDof(sections[i], p, &gcdof));
3562adcc780SMatthew G. Knepley         gtdof = gdof - gcdof;
3572adcc780SMatthew G. Knepley         if (gdof > 0 && gtdof) {
3582adcc780SMatthew G. Knepley           if (bs < 0)           {bs = gtdof;}
3592adcc780SMatthew G. Knepley           else if (bs != gtdof) {bs = 1;}
3609566063dSJacob Faibussowitsch           PetscCall(DMGetGlobalFieldOffset_Private(*superdm, p, startf, &start, &dummy));
3619566063dSJacob Faibussowitsch           PetscCall(DMGetGlobalFieldOffset_Private(*superdm, p, startf+Nfs[i]-1, &dummy, &end));
362*63a3b9bcSJacob Faibussowitsch           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);
3639796d432SMatthew G. Knepley           for (d = start; d < end; ++d, ++subOff) subIndices[subOff] = d;
3642adcc780SMatthew G. Knepley         }
3652adcc780SMatthew G. Knepley       }
3669566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(comm, subSize, subIndices, PETSC_OWN_POINTER, &(*is)[i]));
3672adcc780SMatthew G. Knepley       /* Must have same blocksize on all procs (some might have no points) */
3689796d432SMatthew G. Knepley       {
3699796d432SMatthew G. Knepley         PetscInt bs = -1, bsLocal[2], bsMinMax[2];
3709796d432SMatthew G. Knepley 
3712adcc780SMatthew G. Knepley         bsLocal[0] = bs < 0 ? PETSC_MAX_INT : bs; bsLocal[1] = bs;
3729566063dSJacob Faibussowitsch         PetscCall(PetscGlobalMinMaxInt(comm, bsLocal, bsMinMax));
3732adcc780SMatthew G. Knepley         if (bsMinMax[0] != bsMinMax[1]) {bs = 1;}
3742adcc780SMatthew G. Knepley         else                            {bs = bsMinMax[0];}
3759566063dSJacob Faibussowitsch         PetscCall(ISSetBlockSize((*is)[i], bs));
3762adcc780SMatthew G. Knepley       }
3772adcc780SMatthew G. Knepley     }
3782adcc780SMatthew G. Knepley   }
3799796d432SMatthew G. Knepley   /* Preserve discretizations */
380e5e52638SMatthew G. Knepley   if (len && dms[0]->probs) {
3819566063dSJacob Faibussowitsch     PetscCall(DMSetNumFields(*superdm, Nf));
3829796d432SMatthew G. Knepley     for (i = 0, supf = 0; i < len; ++i) {
3839796d432SMatthew G. Knepley       for (f = 0; f < Nfs[i]; ++f, ++supf) {
3842adcc780SMatthew G. Knepley         PetscObject disc;
3852adcc780SMatthew G. Knepley 
3869566063dSJacob Faibussowitsch         PetscCall(DMGetField(dms[i], f, NULL, &disc));
3879566063dSJacob Faibussowitsch         PetscCall(DMSetField(*superdm, supf, NULL, disc));
3882adcc780SMatthew G. Knepley       }
3892adcc780SMatthew G. Knepley     }
3909566063dSJacob Faibussowitsch     PetscCall(DMCreateDS(*superdm));
3912adcc780SMatthew G. Knepley   }
3929796d432SMatthew G. Knepley   /* Preserve nullspaces */
3939796d432SMatthew G. Knepley   for (i = 0, supf = 0; i < len; ++i) {
3949796d432SMatthew G. Knepley     for (f = 0; f < Nfs[i]; ++f, ++supf) {
3959796d432SMatthew G. Knepley       (*superdm)->nullspaceConstructors[supf] = dms[i]->nullspaceConstructors[f];
3969796d432SMatthew G. Knepley       if ((*superdm)->nullspaceConstructors[supf]) {
3979796d432SMatthew G. Knepley         haveNull = PETSC_TRUE;
3989796d432SMatthew G. Knepley         nullf    = supf;
3998cda7954SMatthew G. Knepley         oldf     = f;
4002adcc780SMatthew G. Knepley       }
4019796d432SMatthew G. Knepley     }
4029796d432SMatthew G. Knepley   }
4039796d432SMatthew G. Knepley   /* Attach nullspace to IS */
4049796d432SMatthew G. Knepley   if (haveNull && is) {
4059796d432SMatthew G. Knepley     MatNullSpace nullSpace;
4069796d432SMatthew G. Knepley 
4079566063dSJacob Faibussowitsch     PetscCall((*(*superdm)->nullspaceConstructors[nullf])(*superdm, oldf, nullf, &nullSpace));
4089566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject) (*is)[nullf], "nullspace", (PetscObject) nullSpace));
4099566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceDestroy(&nullSpace));
4109796d432SMatthew G. Knepley   }
4119566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&supersection));
4129566063dSJacob Faibussowitsch   PetscCall(PetscFree3(Nfs, sections, sectionGlobals));
4132adcc780SMatthew G. Knepley   PetscFunctionReturn(0);
4142adcc780SMatthew G. Knepley }
415