xref: /petsc/src/dm/interface/dmi.c (revision 5d1c02792646264def903282786fa2a600aa1432)
1af0996ceSBarry Smith #include <petsc/private/dmimpl.h>     /*I      "petscdm.h"     I*/
22764a2aaSMatthew G. Knepley #include <petscds.h>
311689aebSJed Brown 
411689aebSJed Brown #undef __FUNCT__
511689aebSJed Brown #define __FUNCT__ "DMCreateGlobalVector_Section_Private"
611689aebSJed Brown PetscErrorCode DMCreateGlobalVector_Section_Private(DM dm,Vec *vec)
711689aebSJed Brown {
811689aebSJed Brown   PetscSection   gSection;
9cc30b04dSMatthew G Knepley   PetscInt       localSize, bs, blockSize = -1, pStart, pEnd, p;
10fad22124SMatthew G Knepley   PetscErrorCode ierr;
11*5d1c0279SHong Zhang   PetscInt       in[2],out[2];
1211689aebSJed Brown 
1311689aebSJed Brown   PetscFunctionBegin;
1411689aebSJed Brown   ierr = DMGetDefaultGlobalSection(dm, &gSection);CHKERRQ(ierr);
15fad22124SMatthew G Knepley   ierr = PetscSectionGetChart(gSection, &pStart, &pEnd);CHKERRQ(ierr);
1611689aebSJed Brown   for (p = pStart; p < pEnd; ++p) {
1711689aebSJed Brown     PetscInt dof, cdof;
1811689aebSJed Brown 
19fad22124SMatthew G Knepley     ierr = PetscSectionGetDof(gSection, p, &dof);CHKERRQ(ierr);
20fad22124SMatthew G Knepley     ierr = PetscSectionGetConstraintDof(gSection, p, &cdof);CHKERRQ(ierr);
211f680262SMatthew G Knepley     if ((blockSize < 0) && (dof > 0) && (dof-cdof > 0)) blockSize = dof-cdof;
2211689aebSJed Brown     if ((dof > 0) && (dof-cdof != blockSize)) {
2311689aebSJed Brown       blockSize = 1;
2411689aebSJed Brown       break;
2511689aebSJed Brown     }
2611689aebSJed Brown   }
27*5d1c0279SHong Zhang 
28*5d1c0279SHong Zhang   in[0] = -blockSize;
29*5d1c0279SHong Zhang   in[1] = blockSize;
30*5d1c0279SHong Zhang   ierr = MPIU_Allreduce(&in,&out,2,MPIU_INT,MPIU_MAX,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
31*5d1c0279SHong Zhang   /* -out[0] = min(blockSize), out[1] = max(blockSize) */
32*5d1c0279SHong Zhang   if (-out[0] == out[1]) {
33*5d1c0279SHong Zhang     bs = out[1];
34*5d1c0279SHong Zhang   } else bs = 1;
35*5d1c0279SHong Zhang 
36*5d1c0279SHong Zhang   if (bs < 0) { /* Everyone was empty */
37*5d1c0279SHong Zhang     blockSize = 1;
38*5d1c0279SHong Zhang     bs        = 1;
39*5d1c0279SHong Zhang   }
40*5d1c0279SHong Zhang 
41fad22124SMatthew G Knepley   ierr = PetscSectionGetConstrainedStorageSize(gSection, &localSize);CHKERRQ(ierr);
4282f516ccSBarry Smith   if (localSize%blockSize) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mismatch between blocksize %d and local storage size %d", blockSize, localSize);
4382f516ccSBarry Smith   ierr = VecCreate(PetscObjectComm((PetscObject)dm), vec);CHKERRQ(ierr);
4411689aebSJed Brown   ierr = VecSetSizes(*vec, localSize, PETSC_DETERMINE);CHKERRQ(ierr);
45cc30b04dSMatthew G Knepley   ierr = VecSetBlockSize(*vec, bs);CHKERRQ(ierr);
46c0dedaeaSBarry Smith   ierr = VecSetType(*vec,dm->vectype);CHKERRQ(ierr);
4711689aebSJed Brown   ierr = VecSetDM(*vec, dm);CHKERRQ(ierr);
4811689aebSJed Brown   /* ierr = VecSetLocalToGlobalMapping(*vec, dm->ltogmap);CHKERRQ(ierr); */
4911689aebSJed Brown   PetscFunctionReturn(0);
5011689aebSJed Brown }
5111689aebSJed Brown 
5211689aebSJed Brown #undef __FUNCT__
5311689aebSJed Brown #define __FUNCT__ "DMCreateLocalVector_Section_Private"
5411689aebSJed Brown PetscErrorCode DMCreateLocalVector_Section_Private(DM dm,Vec *vec)
5511689aebSJed Brown {
56fad22124SMatthew G Knepley   PetscSection   section;
5711689aebSJed Brown   PetscInt       localSize, blockSize = -1, pStart, pEnd, p;
58fad22124SMatthew G Knepley   PetscErrorCode ierr;
5911689aebSJed Brown 
6011689aebSJed Brown   PetscFunctionBegin;
61fad22124SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
62fad22124SMatthew G Knepley   ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
6311689aebSJed Brown   for (p = pStart; p < pEnd; ++p) {
6411689aebSJed Brown     PetscInt dof;
6511689aebSJed Brown 
66fad22124SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
6711689aebSJed Brown     if ((blockSize < 0) && (dof > 0)) blockSize = dof;
6811689aebSJed Brown     if ((dof > 0) && (dof != blockSize)) {
6911689aebSJed Brown       blockSize = 1;
7011689aebSJed Brown       break;
7111689aebSJed Brown     }
7211689aebSJed Brown   }
73fad22124SMatthew G Knepley   ierr = PetscSectionGetStorageSize(section, &localSize);CHKERRQ(ierr);
7411689aebSJed Brown   ierr = VecCreate(PETSC_COMM_SELF, vec);CHKERRQ(ierr);
7511689aebSJed Brown   ierr = VecSetSizes(*vec, localSize, localSize);CHKERRQ(ierr);
7611689aebSJed Brown   ierr = VecSetBlockSize(*vec, blockSize);CHKERRQ(ierr);
77c0dedaeaSBarry Smith   ierr = VecSetType(*vec,dm->vectype);CHKERRQ(ierr);
7811689aebSJed Brown   ierr = VecSetDM(*vec, dm);CHKERRQ(ierr);
7911689aebSJed Brown   PetscFunctionReturn(0);
8011689aebSJed Brown }
814d9407bcSMatthew G. Knepley 
824d9407bcSMatthew G. Knepley #undef __FUNCT__
834d9407bcSMatthew G. Knepley #define __FUNCT__ "DMCreateSubDM_Section_Private"
844d9407bcSMatthew G. Knepley /* This assumes that the DM has been cloned prior to the call */
854d9407bcSMatthew G. Knepley PetscErrorCode DMCreateSubDM_Section_Private(DM dm, PetscInt numFields, 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);
944d9407bcSMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
954d9407bcSMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(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) {
101b9eca265Ssarens     PetscInt bs = -1, bsLocal, bsMax, bsMin;
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) */
127b9eca265Ssarens     bsLocal = bs;
128b9eca265Ssarens     ierr = MPIU_Allreduce(&bsLocal, &bsMax, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
129b9eca265Ssarens     bsLocal = bs < 0 ? bsMax : bs;
130b9eca265Ssarens     ierr = MPIU_Allreduce(&bsLocal, &bsMin, 1, MPIU_INT, MPI_MIN, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
131b9eca265Ssarens     if (bsMin != bsMax) {
132b9eca265Ssarens       bs = 1;
133b9eca265Ssarens     } else {
134b9eca265Ssarens       bs = bsMax;
135b9eca265Ssarens     }
136785e854fSJed Brown     ierr = PetscMalloc1(subSize, &subIndices);CHKERRQ(ierr);
1374d9407bcSMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
1384d9407bcSMatthew G. Knepley       PetscInt gdof, goff;
1394d9407bcSMatthew G. Knepley 
1404d9407bcSMatthew G. Knepley       ierr = PetscSectionGetDof(sectionGlobal, p, &gdof);CHKERRQ(ierr);
1414d9407bcSMatthew G. Knepley       if (gdof > 0) {
1424d9407bcSMatthew G. Knepley         ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr);
1434d9407bcSMatthew G. Knepley         for (f = 0; f < numFields; ++f) {
1444d9407bcSMatthew G. Knepley           PetscInt fdof, fcdof, fc, f2, poff = 0;
1454d9407bcSMatthew G. Knepley 
1464d9407bcSMatthew G. Knepley           /* Can get rid of this loop by storing field information in the global section */
1474d9407bcSMatthew G. Knepley           for (f2 = 0; f2 < fields[f]; ++f2) {
1484d9407bcSMatthew G. Knepley             ierr  = PetscSectionGetFieldDof(section, p, f2, &fdof);CHKERRQ(ierr);
1494d9407bcSMatthew G. Knepley             ierr  = PetscSectionGetFieldConstraintDof(section, p, f2, &fcdof);CHKERRQ(ierr);
1504d9407bcSMatthew G. Knepley             poff += fdof-fcdof;
1514d9407bcSMatthew G. Knepley           }
1524d9407bcSMatthew G. Knepley           ierr = PetscSectionGetFieldDof(section, p, fields[f], &fdof);CHKERRQ(ierr);
1534d9407bcSMatthew G. Knepley           ierr = PetscSectionGetFieldConstraintDof(section, p, fields[f], &fcdof);CHKERRQ(ierr);
1544d9407bcSMatthew G. Knepley           for (fc = 0; fc < fdof-fcdof; ++fc, ++subOff) {
1554d9407bcSMatthew G. Knepley             subIndices[subOff] = goff+poff+fc;
1564d9407bcSMatthew G. Knepley           }
1574d9407bcSMatthew G. Knepley         }
1584d9407bcSMatthew G. Knepley       }
1594d9407bcSMatthew G. Knepley     }
1604d9407bcSMatthew G. Knepley     ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm), subSize, subIndices, PETSC_OWN_POINTER, is);CHKERRQ(ierr);
161b9eca265Ssarens     ierr = ISSetBlockSize(*is, bs);CHKERRQ(ierr);
1624d9407bcSMatthew G. Knepley   }
1634d9407bcSMatthew G. Knepley   if (subdm) {
1644d9407bcSMatthew G. Knepley     PetscSection subsection;
1654d9407bcSMatthew G. Knepley     PetscBool    haveNull = PETSC_FALSE;
1664d9407bcSMatthew G. Knepley     PetscInt     f, nf = 0;
1674d9407bcSMatthew G. Knepley 
1684d9407bcSMatthew G. Knepley     ierr = PetscSectionCreateSubsection(section, numFields, fields, &subsection);CHKERRQ(ierr);
1694d9407bcSMatthew G. Knepley     ierr = DMSetDefaultSection(*subdm, subsection);CHKERRQ(ierr);
1704d9407bcSMatthew G. Knepley     ierr = PetscSectionDestroy(&subsection);CHKERRQ(ierr);
1714d9407bcSMatthew G. Knepley     for (f = 0; f < numFields; ++f) {
1724d9407bcSMatthew G. Knepley       (*subdm)->nullspaceConstructors[f] = dm->nullspaceConstructors[fields[f]];
1734d9407bcSMatthew G. Knepley       if ((*subdm)->nullspaceConstructors[f]) {
1744d9407bcSMatthew G. Knepley         haveNull = PETSC_TRUE;
1754d9407bcSMatthew G. Knepley         nf       = f;
1764d9407bcSMatthew G. Knepley       }
1774d9407bcSMatthew G. Knepley     }
178f646a522SMatthew G. Knepley     if (haveNull && is) {
1794d9407bcSMatthew G. Knepley       MatNullSpace nullSpace;
1804d9407bcSMatthew G. Knepley 
1814d9407bcSMatthew G. Knepley       ierr = (*(*subdm)->nullspaceConstructors[nf])(*subdm, nf, &nullSpace);CHKERRQ(ierr);
1824d9407bcSMatthew G. Knepley       ierr = PetscObjectCompose((PetscObject) *is, "nullspace", (PetscObject) nullSpace);CHKERRQ(ierr);
1834d9407bcSMatthew G. Knepley       ierr = MatNullSpaceDestroy(&nullSpace);CHKERRQ(ierr);
1844d9407bcSMatthew G. Knepley     }
1850f21e855SMatthew G. Knepley     if (dm->prob) {
1860f21e855SMatthew G. Knepley       PetscInt Nf;
1870f21e855SMatthew G. Knepley 
1882764a2aaSMatthew G. Knepley       ierr = PetscDSGetNumFields(dm->prob, &Nf);CHKERRQ(ierr);
1890f21e855SMatthew 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);
1904d9407bcSMatthew G. Knepley       ierr = DMSetNumFields(*subdm, numFields);CHKERRQ(ierr);
1914d9407bcSMatthew G. Knepley       for (f = 0; f < numFields; ++f) {
1920f21e855SMatthew G. Knepley         PetscObject disc;
1930f21e855SMatthew G. Knepley 
1940f21e855SMatthew G. Knepley         ierr = DMGetField(dm, fields[f], &disc);CHKERRQ(ierr);
1950f21e855SMatthew G. Knepley         ierr = DMSetField(*subdm, f, disc);CHKERRQ(ierr);
1964d9407bcSMatthew G. Knepley       }
197f646a522SMatthew G. Knepley       if (numFields == 1 && is) {
1980f21e855SMatthew G. Knepley         PetscObject disc, space, pmat;
1994d9407bcSMatthew G. Knepley 
2000f21e855SMatthew G. Knepley         ierr = DMGetField(*subdm, 0, &disc);CHKERRQ(ierr);
2010f21e855SMatthew G. Knepley         ierr = PetscObjectQuery(disc, "nullspace", &space);CHKERRQ(ierr);
2020f21e855SMatthew G. Knepley         if (space) {ierr = PetscObjectCompose((PetscObject) *is, "nullspace", space);CHKERRQ(ierr);}
2030f21e855SMatthew G. Knepley         ierr = PetscObjectQuery(disc, "nearnullspace", &space);CHKERRQ(ierr);
2040f21e855SMatthew G. Knepley         if (space) {ierr = PetscObjectCompose((PetscObject) *is, "nearnullspace", space);CHKERRQ(ierr);}
2050f21e855SMatthew G. Knepley         ierr = PetscObjectQuery(disc, "pmat", &pmat);CHKERRQ(ierr);
2060f21e855SMatthew G. Knepley         if (pmat) {ierr = PetscObjectCompose((PetscObject) *is, "pmat", pmat);CHKERRQ(ierr);}
2074d9407bcSMatthew G. Knepley       }
2084d9407bcSMatthew G. Knepley     }
2094d9407bcSMatthew G. Knepley   }
2108333ec13SMatthew G. Knepley #if 0
2118333ec13SMatthew G. Knepley   /* We need a way to filter the original SF for given fields:
2128333ec13SMatthew G. Knepley        - Keeping the original section around it too much I think
2138333ec13SMatthew G. Knepley        - We could keep the distributed section, and subset it
2148333ec13SMatthew G. Knepley    */
2158333ec13SMatthew G. Knepley   if (dm->sfNatural) {
2168333ec13SMatthew G. Knepley     PetscSF sfNatural;
2178333ec13SMatthew G. Knepley 
2188333ec13SMatthew G. Knepley     ierr = PetscSectionCreateSubsection(dm->originalSection, numFields, fields, &(*subdm)->originalSection);CHKERRQ(ierr);
2198333ec13SMatthew G. Knepley     ierr = DMPlexCreateGlobalToNaturalPetscSF(*subdm, &sfNatural);CHKERRQ(ierr);
2208333ec13SMatthew G. Knepley     ierr = DMPlexSetGlobalToNaturalPetscSF(*subdm, sfNatural);CHKERRQ(ierr);
2218333ec13SMatthew G. Knepley   }
2228333ec13SMatthew G. Knepley #endif
2234d9407bcSMatthew G. Knepley   PetscFunctionReturn(0);
2244d9407bcSMatthew G. Knepley }
225