Actual source code: dmi.c
petsc-dev 2014-02-02
1: #include <petsc-private/dmimpl.h> /*I "petscdm.h" I*/
5: PetscErrorCode DMCreateGlobalVector_Section_Private(DM dm,Vec *vec)
6: {
7: PetscSection gSection;
8: PetscInt localSize, bs, blockSize = -1, pStart, pEnd, p;
12: DMGetDefaultGlobalSection(dm, &gSection);
13: PetscSectionGetChart(gSection, &pStart, &pEnd);
14: for (p = pStart; p < pEnd; ++p) {
15: PetscInt dof, cdof;
17: PetscSectionGetDof(gSection, p, &dof);
18: PetscSectionGetConstraintDof(gSection, p, &cdof);
19: if ((blockSize < 0) && (dof > 0) && (dof-cdof > 0)) blockSize = dof-cdof;
20: if ((dof > 0) && (dof-cdof != blockSize)) {
21: blockSize = 1;
22: break;
23: }
24: }
25: if (blockSize < 0) blockSize = 1;
26: MPI_Allreduce(&blockSize, &bs, 1, MPIU_INT, MPI_MIN, PetscObjectComm((PetscObject)dm));
27: PetscSectionGetConstrainedStorageSize(gSection, &localSize);
28: if (localSize%blockSize) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mismatch between blocksize %d and local storage size %d", blockSize, localSize);
29: VecCreate(PetscObjectComm((PetscObject)dm), vec);
30: VecSetSizes(*vec, localSize, PETSC_DETERMINE);
31: VecSetBlockSize(*vec, bs);
32: VecSetType(*vec,dm->vectype);
33: VecSetDM(*vec, dm);
34: /* VecSetLocalToGlobalMapping(*vec, dm->ltogmap); */
35: /* VecSetLocalToGlobalMappingBlock(*vec, dm->ltogmapb); */
36: return(0);
37: }
41: PetscErrorCode DMCreateLocalVector_Section_Private(DM dm,Vec *vec)
42: {
43: PetscSection section;
44: PetscInt localSize, blockSize = -1, pStart, pEnd, p;
48: DMGetDefaultSection(dm, §ion);
49: PetscSectionGetChart(section, &pStart, &pEnd);
50: for (p = pStart; p < pEnd; ++p) {
51: PetscInt dof;
53: PetscSectionGetDof(section, p, &dof);
54: if ((blockSize < 0) && (dof > 0)) blockSize = dof;
55: if ((dof > 0) && (dof != blockSize)) {
56: blockSize = 1;
57: break;
58: }
59: }
60: PetscSectionGetStorageSize(section, &localSize);
61: VecCreate(PETSC_COMM_SELF, vec);
62: VecSetSizes(*vec, localSize, localSize);
63: VecSetBlockSize(*vec, blockSize);
64: VecSetType(*vec,dm->vectype);
65: VecSetDM(*vec, dm);
66: return(0);
67: }
71: /* This assumes that the DM has been cloned prior to the call */
72: PetscErrorCode DMCreateSubDM_Section_Private(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm)
73: {
74: PetscSection section, sectionGlobal;
75: PetscInt *subIndices;
76: PetscInt subSize = 0, subOff = 0, nF, f, pStart, pEnd, p;
80: if (!numFields) return(0);
81: DMGetDefaultSection(dm, §ion);
82: DMGetDefaultGlobalSection(dm, §ionGlobal);
83: if (!section) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields");
84: if (!sectionGlobal) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields");
85: PetscSectionGetNumFields(section, &nF);
86: 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);
87: if (is) {
88: PetscSectionGetChart(sectionGlobal, &pStart, &pEnd);
89: for (p = pStart; p < pEnd; ++p) {
90: PetscInt gdof;
92: PetscSectionGetDof(sectionGlobal, p, &gdof);
93: if (gdof > 0) {
94: for (f = 0; f < numFields; ++f) {
95: PetscInt fdof, fcdof;
97: PetscSectionGetFieldDof(section, p, fields[f], &fdof);
98: PetscSectionGetFieldConstraintDof(section, p, fields[f], &fcdof);
99: subSize += fdof-fcdof;
100: }
101: }
102: }
103: PetscMalloc1(subSize, &subIndices);
104: for (p = pStart; p < pEnd; ++p) {
105: PetscInt gdof, goff;
107: PetscSectionGetDof(sectionGlobal, p, &gdof);
108: if (gdof > 0) {
109: PetscSectionGetOffset(sectionGlobal, p, &goff);
110: for (f = 0; f < numFields; ++f) {
111: PetscInt fdof, fcdof, fc, f2, poff = 0;
113: /* Can get rid of this loop by storing field information in the global section */
114: for (f2 = 0; f2 < fields[f]; ++f2) {
115: PetscSectionGetFieldDof(section, p, f2, &fdof);
116: PetscSectionGetFieldConstraintDof(section, p, f2, &fcdof);
117: poff += fdof-fcdof;
118: }
119: PetscSectionGetFieldDof(section, p, fields[f], &fdof);
120: PetscSectionGetFieldConstraintDof(section, p, fields[f], &fcdof);
121: for (fc = 0; fc < fdof-fcdof; ++fc, ++subOff) {
122: subIndices[subOff] = goff+poff+fc;
123: }
124: }
125: }
126: }
127: ISCreateGeneral(PetscObjectComm((PetscObject)dm), subSize, subIndices, PETSC_OWN_POINTER, is);
128: }
129: if (subdm) {
130: PetscSection subsection;
131: PetscBool haveNull = PETSC_FALSE;
132: PetscInt f, nf = 0;
134: PetscSectionCreateSubsection(section, numFields, fields, &subsection);
135: DMSetDefaultSection(*subdm, subsection);
136: PetscSectionDestroy(&subsection);
137: for (f = 0; f < numFields; ++f) {
138: (*subdm)->nullspaceConstructors[f] = dm->nullspaceConstructors[fields[f]];
139: if ((*subdm)->nullspaceConstructors[f]) {
140: haveNull = PETSC_TRUE;
141: nf = f;
142: }
143: }
144: if (haveNull) {
145: MatNullSpace nullSpace;
147: (*(*subdm)->nullspaceConstructors[nf])(*subdm, nf, &nullSpace);
148: PetscObjectCompose((PetscObject) *is, "nullspace", (PetscObject) nullSpace);
149: MatNullSpaceDestroy(&nullSpace);
150: }
151: if (dm->fields) {
152: 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);
153: DMSetNumFields(*subdm, numFields);
154: for (f = 0; f < numFields; ++f) {
155: PetscObjectListDuplicate(dm->fields[fields[f]]->olist, &(*subdm)->fields[f]->olist);
156: }
157: if (numFields == 1) {
158: MatNullSpace space;
159: Mat pmat;
161: PetscObjectQuery((*subdm)->fields[0], "nullspace", (PetscObject*) &space);
162: if (space) {PetscObjectCompose((PetscObject) *is, "nullspace", (PetscObject) space);}
163: PetscObjectQuery((*subdm)->fields[0], "nearnullspace", (PetscObject*) &space);
164: if (space) {PetscObjectCompose((PetscObject) *is, "nearnullspace", (PetscObject) space);}
165: PetscObjectQuery((*subdm)->fields[0], "pmat", (PetscObject*) &pmat);
166: if (pmat) {PetscObjectCompose((PetscObject) *is, "pmat", (PetscObject) pmat);}
167: }
168: }
169: }
170: return(0);
171: }