xref: /petsc/src/dm/interface/dmi.c (revision 8333ec13d514b7c7db3166bad5929e34a5079e59)
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;
1111689aebSJed Brown 
1211689aebSJed Brown   PetscFunctionBegin;
1311689aebSJed Brown   ierr = DMGetDefaultGlobalSection(dm, &gSection);CHKERRQ(ierr);
14fad22124SMatthew G Knepley   ierr = PetscSectionGetChart(gSection, &pStart, &pEnd);CHKERRQ(ierr);
1511689aebSJed Brown   for (p = pStart; p < pEnd; ++p) {
1611689aebSJed Brown     PetscInt dof, cdof;
1711689aebSJed Brown 
18fad22124SMatthew G Knepley     ierr = PetscSectionGetDof(gSection, p, &dof);CHKERRQ(ierr);
19fad22124SMatthew G Knepley     ierr = PetscSectionGetConstraintDof(gSection, p, &cdof);CHKERRQ(ierr);
201f680262SMatthew G Knepley     if ((blockSize < 0) && (dof > 0) && (dof-cdof > 0)) blockSize = dof-cdof;
2111689aebSJed Brown     if ((dof > 0) && (dof-cdof != blockSize)) {
2211689aebSJed Brown       blockSize = 1;
2311689aebSJed Brown       break;
2411689aebSJed Brown     }
2511689aebSJed Brown   }
266bf20b62SMatthew G. Knepley   if (blockSize < 0) blockSize = PETSC_MAX_INT;
2782f516ccSBarry Smith   ierr = MPI_Allreduce(&blockSize, &bs, 1, MPIU_INT, MPI_MIN, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
286bf20b62SMatthew G. Knepley   if (blockSize == PETSC_MAX_INT) blockSize = 1; /* Everyone was empty */
29fad22124SMatthew G Knepley   ierr = PetscSectionGetConstrainedStorageSize(gSection, &localSize);CHKERRQ(ierr);
3082f516ccSBarry Smith   if (localSize%blockSize) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mismatch between blocksize %d and local storage size %d", blockSize, localSize);
3182f516ccSBarry Smith   ierr = VecCreate(PetscObjectComm((PetscObject)dm), vec);CHKERRQ(ierr);
3211689aebSJed Brown   ierr = VecSetSizes(*vec, localSize, PETSC_DETERMINE);CHKERRQ(ierr);
33cc30b04dSMatthew G Knepley   ierr = VecSetBlockSize(*vec, bs);CHKERRQ(ierr);
34c0dedaeaSBarry Smith   ierr = VecSetType(*vec,dm->vectype);CHKERRQ(ierr);
3511689aebSJed Brown   ierr = VecSetDM(*vec, dm);CHKERRQ(ierr);
3611689aebSJed Brown   /* ierr = VecSetLocalToGlobalMapping(*vec, dm->ltogmap);CHKERRQ(ierr); */
3711689aebSJed Brown   PetscFunctionReturn(0);
3811689aebSJed Brown }
3911689aebSJed Brown 
4011689aebSJed Brown #undef __FUNCT__
4111689aebSJed Brown #define __FUNCT__ "DMCreateLocalVector_Section_Private"
4211689aebSJed Brown PetscErrorCode DMCreateLocalVector_Section_Private(DM dm,Vec *vec)
4311689aebSJed Brown {
44fad22124SMatthew G Knepley   PetscSection   section;
4511689aebSJed Brown   PetscInt       localSize, blockSize = -1, pStart, pEnd, p;
46fad22124SMatthew G Knepley   PetscErrorCode ierr;
4711689aebSJed Brown 
4811689aebSJed Brown   PetscFunctionBegin;
49fad22124SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
50fad22124SMatthew G Knepley   ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
5111689aebSJed Brown   for (p = pStart; p < pEnd; ++p) {
5211689aebSJed Brown     PetscInt dof;
5311689aebSJed Brown 
54fad22124SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
5511689aebSJed Brown     if ((blockSize < 0) && (dof > 0)) blockSize = dof;
5611689aebSJed Brown     if ((dof > 0) && (dof != blockSize)) {
5711689aebSJed Brown       blockSize = 1;
5811689aebSJed Brown       break;
5911689aebSJed Brown     }
6011689aebSJed Brown   }
61fad22124SMatthew G Knepley   ierr = PetscSectionGetStorageSize(section, &localSize);CHKERRQ(ierr);
6211689aebSJed Brown   ierr = VecCreate(PETSC_COMM_SELF, vec);CHKERRQ(ierr);
6311689aebSJed Brown   ierr = VecSetSizes(*vec, localSize, localSize);CHKERRQ(ierr);
6411689aebSJed Brown   ierr = VecSetBlockSize(*vec, blockSize);CHKERRQ(ierr);
65c0dedaeaSBarry Smith   ierr = VecSetType(*vec,dm->vectype);CHKERRQ(ierr);
6611689aebSJed Brown   ierr = VecSetDM(*vec, dm);CHKERRQ(ierr);
6711689aebSJed Brown   PetscFunctionReturn(0);
6811689aebSJed Brown }
694d9407bcSMatthew G. Knepley 
704d9407bcSMatthew G. Knepley #undef __FUNCT__
714d9407bcSMatthew G. Knepley #define __FUNCT__ "DMCreateSubDM_Section_Private"
724d9407bcSMatthew G. Knepley /* This assumes that the DM has been cloned prior to the call */
734d9407bcSMatthew G. Knepley PetscErrorCode DMCreateSubDM_Section_Private(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm)
744d9407bcSMatthew G. Knepley {
754d9407bcSMatthew G. Knepley   PetscSection   section, sectionGlobal;
764d9407bcSMatthew G. Knepley   PetscInt      *subIndices;
774d9407bcSMatthew G. Knepley   PetscInt       subSize = 0, subOff = 0, nF, f, pStart, pEnd, p;
784d9407bcSMatthew G. Knepley   PetscErrorCode ierr;
794d9407bcSMatthew G. Knepley 
804d9407bcSMatthew G. Knepley   PetscFunctionBegin;
814d9407bcSMatthew G. Knepley   if (!numFields) PetscFunctionReturn(0);
824d9407bcSMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
834d9407bcSMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dm, &sectionGlobal);CHKERRQ(ierr);
844d9407bcSMatthew G. Knepley   if (!section) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields");
854d9407bcSMatthew G. Knepley   if (!sectionGlobal) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields");
864d9407bcSMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &nF);CHKERRQ(ierr);
874d9407bcSMatthew 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);
884d9407bcSMatthew G. Knepley   if (is) {
894d9407bcSMatthew G. Knepley     ierr = PetscSectionGetChart(sectionGlobal, &pStart, &pEnd);CHKERRQ(ierr);
904d9407bcSMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
914d9407bcSMatthew G. Knepley       PetscInt gdof;
924d9407bcSMatthew G. Knepley 
934d9407bcSMatthew G. Knepley       ierr = PetscSectionGetDof(sectionGlobal, p, &gdof);CHKERRQ(ierr);
944d9407bcSMatthew G. Knepley       if (gdof > 0) {
954d9407bcSMatthew G. Knepley         for (f = 0; f < numFields; ++f) {
964d9407bcSMatthew G. Knepley           PetscInt fdof, fcdof;
974d9407bcSMatthew G. Knepley 
984d9407bcSMatthew G. Knepley           ierr     = PetscSectionGetFieldDof(section, p, fields[f], &fdof);CHKERRQ(ierr);
994d9407bcSMatthew G. Knepley           ierr     = PetscSectionGetFieldConstraintDof(section, p, fields[f], &fcdof);CHKERRQ(ierr);
1004d9407bcSMatthew G. Knepley           subSize += fdof-fcdof;
1014d9407bcSMatthew G. Knepley         }
1024d9407bcSMatthew G. Knepley       }
1034d9407bcSMatthew G. Knepley     }
104785e854fSJed Brown     ierr = PetscMalloc1(subSize, &subIndices);CHKERRQ(ierr);
1054d9407bcSMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
1064d9407bcSMatthew G. Knepley       PetscInt gdof, goff;
1074d9407bcSMatthew G. Knepley 
1084d9407bcSMatthew G. Knepley       ierr = PetscSectionGetDof(sectionGlobal, p, &gdof);CHKERRQ(ierr);
1094d9407bcSMatthew G. Knepley       if (gdof > 0) {
1104d9407bcSMatthew G. Knepley         ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr);
1114d9407bcSMatthew G. Knepley         for (f = 0; f < numFields; ++f) {
1124d9407bcSMatthew G. Knepley           PetscInt fdof, fcdof, fc, f2, poff = 0;
1134d9407bcSMatthew G. Knepley 
1144d9407bcSMatthew G. Knepley           /* Can get rid of this loop by storing field information in the global section */
1154d9407bcSMatthew G. Knepley           for (f2 = 0; f2 < fields[f]; ++f2) {
1164d9407bcSMatthew G. Knepley             ierr  = PetscSectionGetFieldDof(section, p, f2, &fdof);CHKERRQ(ierr);
1174d9407bcSMatthew G. Knepley             ierr  = PetscSectionGetFieldConstraintDof(section, p, f2, &fcdof);CHKERRQ(ierr);
1184d9407bcSMatthew G. Knepley             poff += fdof-fcdof;
1194d9407bcSMatthew G. Knepley           }
1204d9407bcSMatthew G. Knepley           ierr = PetscSectionGetFieldDof(section, p, fields[f], &fdof);CHKERRQ(ierr);
1214d9407bcSMatthew G. Knepley           ierr = PetscSectionGetFieldConstraintDof(section, p, fields[f], &fcdof);CHKERRQ(ierr);
1224d9407bcSMatthew G. Knepley           for (fc = 0; fc < fdof-fcdof; ++fc, ++subOff) {
1234d9407bcSMatthew G. Knepley             subIndices[subOff] = goff+poff+fc;
1244d9407bcSMatthew G. Knepley           }
1254d9407bcSMatthew G. Knepley         }
1264d9407bcSMatthew G. Knepley       }
1274d9407bcSMatthew G. Knepley     }
1284d9407bcSMatthew G. Knepley     ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm), subSize, subIndices, PETSC_OWN_POINTER, is);CHKERRQ(ierr);
1294d9407bcSMatthew G. Knepley   }
1304d9407bcSMatthew G. Knepley   if (subdm) {
1314d9407bcSMatthew G. Knepley     PetscSection subsection;
1324d9407bcSMatthew G. Knepley     PetscBool    haveNull = PETSC_FALSE;
1334d9407bcSMatthew G. Knepley     PetscInt     f, nf = 0;
1344d9407bcSMatthew G. Knepley 
1354d9407bcSMatthew G. Knepley     ierr = PetscSectionCreateSubsection(section, numFields, fields, &subsection);CHKERRQ(ierr);
1364d9407bcSMatthew G. Knepley     ierr = DMSetDefaultSection(*subdm, subsection);CHKERRQ(ierr);
1374d9407bcSMatthew G. Knepley     ierr = PetscSectionDestroy(&subsection);CHKERRQ(ierr);
1384d9407bcSMatthew G. Knepley     for (f = 0; f < numFields; ++f) {
1394d9407bcSMatthew G. Knepley       (*subdm)->nullspaceConstructors[f] = dm->nullspaceConstructors[fields[f]];
1404d9407bcSMatthew G. Knepley       if ((*subdm)->nullspaceConstructors[f]) {
1414d9407bcSMatthew G. Knepley         haveNull = PETSC_TRUE;
1424d9407bcSMatthew G. Knepley         nf       = f;
1434d9407bcSMatthew G. Knepley       }
1444d9407bcSMatthew G. Knepley     }
145f646a522SMatthew G. Knepley     if (haveNull && is) {
1464d9407bcSMatthew G. Knepley       MatNullSpace nullSpace;
1474d9407bcSMatthew G. Knepley 
1484d9407bcSMatthew G. Knepley       ierr = (*(*subdm)->nullspaceConstructors[nf])(*subdm, nf, &nullSpace);CHKERRQ(ierr);
1494d9407bcSMatthew G. Knepley       ierr = PetscObjectCompose((PetscObject) *is, "nullspace", (PetscObject) nullSpace);CHKERRQ(ierr);
1504d9407bcSMatthew G. Knepley       ierr = MatNullSpaceDestroy(&nullSpace);CHKERRQ(ierr);
1514d9407bcSMatthew G. Knepley     }
1520f21e855SMatthew G. Knepley     if (dm->prob) {
1530f21e855SMatthew G. Knepley       PetscInt Nf;
1540f21e855SMatthew G. Knepley 
1552764a2aaSMatthew G. Knepley       ierr = PetscDSGetNumFields(dm->prob, &Nf);CHKERRQ(ierr);
1560f21e855SMatthew 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);
1574d9407bcSMatthew G. Knepley       ierr = DMSetNumFields(*subdm, numFields);CHKERRQ(ierr);
1584d9407bcSMatthew G. Knepley       for (f = 0; f < numFields; ++f) {
1590f21e855SMatthew G. Knepley         PetscObject disc;
1600f21e855SMatthew G. Knepley 
1610f21e855SMatthew G. Knepley         ierr = DMGetField(dm, fields[f], &disc);CHKERRQ(ierr);
1620f21e855SMatthew G. Knepley         ierr = DMSetField(*subdm, f, disc);CHKERRQ(ierr);
1634d9407bcSMatthew G. Knepley       }
164f646a522SMatthew G. Knepley       if (numFields == 1 && is) {
1650f21e855SMatthew G. Knepley         PetscObject disc, space, pmat;
1664d9407bcSMatthew G. Knepley 
1670f21e855SMatthew G. Knepley         ierr = DMGetField(*subdm, 0, &disc);CHKERRQ(ierr);
1680f21e855SMatthew G. Knepley         ierr = PetscObjectQuery(disc, "nullspace", &space);CHKERRQ(ierr);
1690f21e855SMatthew G. Knepley         if (space) {ierr = PetscObjectCompose((PetscObject) *is, "nullspace", space);CHKERRQ(ierr);}
1700f21e855SMatthew G. Knepley         ierr = PetscObjectQuery(disc, "nearnullspace", &space);CHKERRQ(ierr);
1710f21e855SMatthew G. Knepley         if (space) {ierr = PetscObjectCompose((PetscObject) *is, "nearnullspace", space);CHKERRQ(ierr);}
1720f21e855SMatthew G. Knepley         ierr = PetscObjectQuery(disc, "pmat", &pmat);CHKERRQ(ierr);
1730f21e855SMatthew G. Knepley         if (pmat) {ierr = PetscObjectCompose((PetscObject) *is, "pmat", pmat);CHKERRQ(ierr);}
1744d9407bcSMatthew G. Knepley       }
1754d9407bcSMatthew G. Knepley     }
1764d9407bcSMatthew G. Knepley   }
177*8333ec13SMatthew G. Knepley #if 0
178*8333ec13SMatthew G. Knepley   /* We need a way to filter the original SF for given fields:
179*8333ec13SMatthew G. Knepley        - Keeping the original section around it too much I think
180*8333ec13SMatthew G. Knepley        - We could keep the distributed section, and subset it
181*8333ec13SMatthew G. Knepley    */
182*8333ec13SMatthew G. Knepley   if (dm->sfNatural) {
183*8333ec13SMatthew G. Knepley     PetscSF sfNatural;
184*8333ec13SMatthew G. Knepley 
185*8333ec13SMatthew G. Knepley     ierr = PetscSectionCreateSubsection(dm->originalSection, numFields, fields, &(*subdm)->originalSection);CHKERRQ(ierr);
186*8333ec13SMatthew G. Knepley     ierr = DMPlexCreateGlobalToNaturalPetscSF(*subdm, &sfNatural);CHKERRQ(ierr);
187*8333ec13SMatthew G. Knepley     ierr = DMPlexSetGlobalToNaturalPetscSF(*subdm, sfNatural);CHKERRQ(ierr);
188*8333ec13SMatthew G. Knepley   }
189*8333ec13SMatthew G. Knepley #endif
1904d9407bcSMatthew G. Knepley   PetscFunctionReturn(0);
1914d9407bcSMatthew G. Knepley }
192