111689aebSJed Brown #include <petsc-private/dmimpl.h> /*I "petscdm.h" I*/ 211689aebSJed Brown 311689aebSJed Brown #undef __FUNCT__ 411689aebSJed Brown #define __FUNCT__ "DMCreateGlobalVector_Section_Private" 511689aebSJed Brown PetscErrorCode DMCreateGlobalVector_Section_Private(DM dm,Vec *vec) 611689aebSJed Brown { 711689aebSJed Brown PetscSection gSection; 8cc30b04dSMatthew G Knepley PetscInt localSize, bs, blockSize = -1, pStart, pEnd, p; 9fad22124SMatthew G Knepley PetscErrorCode ierr; 1011689aebSJed Brown 1111689aebSJed Brown PetscFunctionBegin; 1211689aebSJed Brown ierr = DMGetDefaultGlobalSection(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); 191f680262SMatthew G Knepley if ((blockSize < 0) && (dof > 0) && (dof-cdof > 0)) blockSize = dof-cdof; 2011689aebSJed Brown if ((dof > 0) && (dof-cdof != blockSize)) { 2111689aebSJed Brown blockSize = 1; 2211689aebSJed Brown break; 2311689aebSJed Brown } 2411689aebSJed Brown } 25cc30b04dSMatthew G Knepley if (blockSize < 0) blockSize = 1; 2682f516ccSBarry Smith ierr = MPI_Allreduce(&blockSize, &bs, 1, MPIU_INT, MPI_MIN, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 27fad22124SMatthew G Knepley ierr = PetscSectionGetConstrainedStorageSize(gSection, &localSize);CHKERRQ(ierr); 2882f516ccSBarry Smith if (localSize%blockSize) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mismatch between blocksize %d and local storage size %d", blockSize, localSize); 2982f516ccSBarry Smith ierr = VecCreate(PetscObjectComm((PetscObject)dm), vec);CHKERRQ(ierr); 3011689aebSJed Brown ierr = VecSetSizes(*vec, localSize, PETSC_DETERMINE);CHKERRQ(ierr); 31cc30b04dSMatthew G Knepley ierr = VecSetBlockSize(*vec, bs);CHKERRQ(ierr); 32c0dedaeaSBarry Smith ierr = VecSetType(*vec,dm->vectype);CHKERRQ(ierr); 3311689aebSJed Brown ierr = VecSetDM(*vec, dm);CHKERRQ(ierr); 3411689aebSJed Brown /* ierr = VecSetLocalToGlobalMapping(*vec, dm->ltogmap);CHKERRQ(ierr); */ 3511689aebSJed Brown /* ierr = VecSetLocalToGlobalMappingBlock(*vec, dm->ltogmapb);CHKERRQ(ierr); */ 3611689aebSJed Brown PetscFunctionReturn(0); 3711689aebSJed Brown } 3811689aebSJed Brown 3911689aebSJed Brown #undef __FUNCT__ 4011689aebSJed Brown #define __FUNCT__ "DMCreateLocalVector_Section_Private" 4111689aebSJed Brown PetscErrorCode DMCreateLocalVector_Section_Private(DM dm,Vec *vec) 4211689aebSJed Brown { 43fad22124SMatthew G Knepley PetscSection section; 4411689aebSJed Brown PetscInt localSize, blockSize = -1, pStart, pEnd, p; 45fad22124SMatthew G Knepley PetscErrorCode ierr; 4611689aebSJed Brown 4711689aebSJed Brown PetscFunctionBegin; 48fad22124SMatthew G Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 49fad22124SMatthew G Knepley ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr); 5011689aebSJed Brown for (p = pStart; p < pEnd; ++p) { 5111689aebSJed Brown PetscInt dof; 5211689aebSJed Brown 53fad22124SMatthew G Knepley ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 5411689aebSJed Brown if ((blockSize < 0) && (dof > 0)) blockSize = dof; 5511689aebSJed Brown if ((dof > 0) && (dof != blockSize)) { 5611689aebSJed Brown blockSize = 1; 5711689aebSJed Brown break; 5811689aebSJed Brown } 5911689aebSJed Brown } 60fad22124SMatthew G Knepley ierr = PetscSectionGetStorageSize(section, &localSize);CHKERRQ(ierr); 6111689aebSJed Brown ierr = VecCreate(PETSC_COMM_SELF, vec);CHKERRQ(ierr); 6211689aebSJed Brown ierr = VecSetSizes(*vec, localSize, localSize);CHKERRQ(ierr); 6311689aebSJed Brown ierr = VecSetBlockSize(*vec, blockSize);CHKERRQ(ierr); 64c0dedaeaSBarry Smith ierr = VecSetType(*vec,dm->vectype);CHKERRQ(ierr); 6511689aebSJed Brown ierr = VecSetDM(*vec, dm);CHKERRQ(ierr); 6611689aebSJed Brown PetscFunctionReturn(0); 6711689aebSJed Brown } 684d9407bcSMatthew G. Knepley 694d9407bcSMatthew G. Knepley #undef __FUNCT__ 704d9407bcSMatthew G. Knepley #define __FUNCT__ "DMCreateSubDM_Section_Private" 714d9407bcSMatthew G. Knepley /* This assumes that the DM has been cloned prior to the call */ 724d9407bcSMatthew G. Knepley PetscErrorCode DMCreateSubDM_Section_Private(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm) 734d9407bcSMatthew G. Knepley { 744d9407bcSMatthew G. Knepley PetscSection section, sectionGlobal; 754d9407bcSMatthew G. Knepley PetscInt *subIndices; 764d9407bcSMatthew G. Knepley PetscInt subSize = 0, subOff = 0, nF, f, pStart, pEnd, p; 774d9407bcSMatthew G. Knepley PetscErrorCode ierr; 784d9407bcSMatthew G. Knepley 794d9407bcSMatthew G. Knepley PetscFunctionBegin; 804d9407bcSMatthew G. Knepley if (!numFields) PetscFunctionReturn(0); 814d9407bcSMatthew G. Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 824d9407bcSMatthew G. Knepley ierr = DMGetDefaultGlobalSection(dm, §ionGlobal);CHKERRQ(ierr); 834d9407bcSMatthew G. Knepley if (!section) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields"); 844d9407bcSMatthew G. Knepley if (!sectionGlobal) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields"); 854d9407bcSMatthew G. Knepley ierr = PetscSectionGetNumFields(section, &nF);CHKERRQ(ierr); 864d9407bcSMatthew 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); 874d9407bcSMatthew G. Knepley if (is) { 884d9407bcSMatthew G. Knepley ierr = PetscSectionGetChart(sectionGlobal, &pStart, &pEnd);CHKERRQ(ierr); 894d9407bcSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 904d9407bcSMatthew G. Knepley PetscInt gdof; 914d9407bcSMatthew G. Knepley 924d9407bcSMatthew G. Knepley ierr = PetscSectionGetDof(sectionGlobal, p, &gdof);CHKERRQ(ierr); 934d9407bcSMatthew G. Knepley if (gdof > 0) { 944d9407bcSMatthew G. Knepley for (f = 0; f < numFields; ++f) { 954d9407bcSMatthew G. Knepley PetscInt fdof, fcdof; 964d9407bcSMatthew G. Knepley 974d9407bcSMatthew G. Knepley ierr = PetscSectionGetFieldDof(section, p, fields[f], &fdof);CHKERRQ(ierr); 984d9407bcSMatthew G. Knepley ierr = PetscSectionGetFieldConstraintDof(section, p, fields[f], &fcdof);CHKERRQ(ierr); 994d9407bcSMatthew G. Knepley subSize += fdof-fcdof; 1004d9407bcSMatthew G. Knepley } 1014d9407bcSMatthew G. Knepley } 1024d9407bcSMatthew G. Knepley } 103*785e854fSJed Brown ierr = PetscMalloc1(subSize, &subIndices);CHKERRQ(ierr); 1044d9407bcSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 1054d9407bcSMatthew G. Knepley PetscInt gdof, goff; 1064d9407bcSMatthew G. Knepley 1074d9407bcSMatthew G. Knepley ierr = PetscSectionGetDof(sectionGlobal, p, &gdof);CHKERRQ(ierr); 1084d9407bcSMatthew G. Knepley if (gdof > 0) { 1094d9407bcSMatthew G. Knepley ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr); 1104d9407bcSMatthew G. Knepley for (f = 0; f < numFields; ++f) { 1114d9407bcSMatthew G. Knepley PetscInt fdof, fcdof, fc, f2, poff = 0; 1124d9407bcSMatthew G. Knepley 1134d9407bcSMatthew G. Knepley /* Can get rid of this loop by storing field information in the global section */ 1144d9407bcSMatthew G. Knepley for (f2 = 0; f2 < fields[f]; ++f2) { 1154d9407bcSMatthew G. Knepley ierr = PetscSectionGetFieldDof(section, p, f2, &fdof);CHKERRQ(ierr); 1164d9407bcSMatthew G. Knepley ierr = PetscSectionGetFieldConstraintDof(section, p, f2, &fcdof);CHKERRQ(ierr); 1174d9407bcSMatthew G. Knepley poff += fdof-fcdof; 1184d9407bcSMatthew G. Knepley } 1194d9407bcSMatthew G. Knepley ierr = PetscSectionGetFieldDof(section, p, fields[f], &fdof);CHKERRQ(ierr); 1204d9407bcSMatthew G. Knepley ierr = PetscSectionGetFieldConstraintDof(section, p, fields[f], &fcdof);CHKERRQ(ierr); 1214d9407bcSMatthew G. Knepley for (fc = 0; fc < fdof-fcdof; ++fc, ++subOff) { 1224d9407bcSMatthew G. Knepley subIndices[subOff] = goff+poff+fc; 1234d9407bcSMatthew G. Knepley } 1244d9407bcSMatthew G. Knepley } 1254d9407bcSMatthew G. Knepley } 1264d9407bcSMatthew G. Knepley } 1274d9407bcSMatthew G. Knepley ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm), subSize, subIndices, PETSC_OWN_POINTER, is);CHKERRQ(ierr); 1284d9407bcSMatthew G. Knepley } 1294d9407bcSMatthew G. Knepley if (subdm) { 1304d9407bcSMatthew G. Knepley PetscSection subsection; 1314d9407bcSMatthew G. Knepley PetscBool haveNull = PETSC_FALSE; 1324d9407bcSMatthew G. Knepley PetscInt f, nf = 0; 1334d9407bcSMatthew G. Knepley 1344d9407bcSMatthew G. Knepley ierr = PetscSectionCreateSubsection(section, numFields, fields, &subsection);CHKERRQ(ierr); 1354d9407bcSMatthew G. Knepley ierr = DMSetDefaultSection(*subdm, subsection);CHKERRQ(ierr); 1364d9407bcSMatthew G. Knepley ierr = PetscSectionDestroy(&subsection);CHKERRQ(ierr); 1374d9407bcSMatthew G. Knepley for (f = 0; f < numFields; ++f) { 1384d9407bcSMatthew G. Knepley (*subdm)->nullspaceConstructors[f] = dm->nullspaceConstructors[fields[f]]; 1394d9407bcSMatthew G. Knepley if ((*subdm)->nullspaceConstructors[f]) { 1404d9407bcSMatthew G. Knepley haveNull = PETSC_TRUE; 1414d9407bcSMatthew G. Knepley nf = f; 1424d9407bcSMatthew G. Knepley } 1434d9407bcSMatthew G. Knepley } 1444d9407bcSMatthew G. Knepley if (haveNull) { 1454d9407bcSMatthew G. Knepley MatNullSpace nullSpace; 1464d9407bcSMatthew G. Knepley 1474d9407bcSMatthew G. Knepley ierr = (*(*subdm)->nullspaceConstructors[nf])(*subdm, nf, &nullSpace);CHKERRQ(ierr); 1484d9407bcSMatthew G. Knepley ierr = PetscObjectCompose((PetscObject) *is, "nullspace", (PetscObject) nullSpace);CHKERRQ(ierr); 1494d9407bcSMatthew G. Knepley ierr = MatNullSpaceDestroy(&nullSpace);CHKERRQ(ierr); 1504d9407bcSMatthew G. Knepley } 1514d9407bcSMatthew G. Knepley if (dm->fields) { 1524d9407bcSMatthew G. Knepley if (nF != dm->numFields) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "The number of DM fields %d does not match the number of Section fields %d", dm->numFields, nF); 1534d9407bcSMatthew G. Knepley ierr = DMSetNumFields(*subdm, numFields);CHKERRQ(ierr); 1544d9407bcSMatthew G. Knepley for (f = 0; f < numFields; ++f) { 1554d9407bcSMatthew G. Knepley ierr = PetscObjectListDuplicate(dm->fields[fields[f]]->olist, &(*subdm)->fields[f]->olist);CHKERRQ(ierr); 1564d9407bcSMatthew G. Knepley } 1574d9407bcSMatthew G. Knepley if (numFields == 1) { 1584d9407bcSMatthew G. Knepley MatNullSpace space; 1594d9407bcSMatthew G. Knepley Mat pmat; 1604d9407bcSMatthew G. Knepley 1614d9407bcSMatthew G. Knepley ierr = PetscObjectQuery((*subdm)->fields[0], "nullspace", (PetscObject*) &space);CHKERRQ(ierr); 1624d9407bcSMatthew G. Knepley if (space) {ierr = PetscObjectCompose((PetscObject) *is, "nullspace", (PetscObject) space);CHKERRQ(ierr);} 1634d9407bcSMatthew G. Knepley ierr = PetscObjectQuery((*subdm)->fields[0], "nearnullspace", (PetscObject*) &space);CHKERRQ(ierr); 1644d9407bcSMatthew G. Knepley if (space) {ierr = PetscObjectCompose((PetscObject) *is, "nearnullspace", (PetscObject) space);CHKERRQ(ierr);} 1654d9407bcSMatthew G. Knepley ierr = PetscObjectQuery((*subdm)->fields[0], "pmat", (PetscObject*) &pmat);CHKERRQ(ierr); 1664d9407bcSMatthew G. Knepley if (pmat) {ierr = PetscObjectCompose((PetscObject) *is, "pmat", (PetscObject) pmat);CHKERRQ(ierr);} 1674d9407bcSMatthew G. Knepley } 1684d9407bcSMatthew G. Knepley } 1694d9407bcSMatthew G. Knepley } 1704d9407bcSMatthew G. Knepley PetscFunctionReturn(0); 1714d9407bcSMatthew G. Knepley } 172