Actual source code: dmi.c
petsc-3.8.4 2018-03-24
1: #include <petsc/private/dmimpl.h>
2: #include <petscds.h>
4: PetscErrorCode DMCreateGlobalVector_Section_Private(DM dm,Vec *vec)
5: {
6: PetscSection gSection;
7: PetscInt localSize, bs, blockSize = -1, pStart, pEnd, p;
9: PetscInt in[2],out[2];
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);
20: if (dof > 0) {
21: if (blockSize < 0 && dof-cdof > 0) {
22: /* set blockSize */
23: blockSize = dof-cdof;
24: } else if (dof-cdof != blockSize) {
25: /* non-identical blockSize, set it as 1 */
26: blockSize = 1;
27: break;
28: }
29: }
30: }
32: in[0] = blockSize < 0 ? PETSC_MIN_INT : -blockSize;
33: in[1] = blockSize;
34: MPIU_Allreduce(in,out,2,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)dm));
35: /* -out[0] = min(blockSize), out[1] = max(blockSize) */
36: if (-out[0] == out[1]) {
37: bs = out[1];
38: } else bs = 1;
40: if (bs < 0) { /* Everyone was empty */
41: blockSize = 1;
42: bs = 1;
43: }
45: PetscSectionGetConstrainedStorageSize(gSection, &localSize);
46: if (localSize%blockSize) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mismatch between blocksize %d and local storage size %d", blockSize, localSize);
47: VecCreate(PetscObjectComm((PetscObject)dm), vec);
48: VecSetSizes(*vec, localSize, PETSC_DETERMINE);
49: VecSetBlockSize(*vec, bs);
50: VecSetType(*vec,dm->vectype);
51: VecSetDM(*vec, dm);
52: /* VecSetLocalToGlobalMapping(*vec, dm->ltogmap); */
53: return(0);
54: }
56: PetscErrorCode DMCreateLocalVector_Section_Private(DM dm,Vec *vec)
57: {
58: PetscSection section;
59: PetscInt localSize, blockSize = -1, pStart, pEnd, p;
63: DMGetDefaultSection(dm, §ion);
64: PetscSectionGetChart(section, &pStart, &pEnd);
65: for (p = pStart; p < pEnd; ++p) {
66: PetscInt dof;
68: PetscSectionGetDof(section, p, &dof);
69: if ((blockSize < 0) && (dof > 0)) blockSize = dof;
70: if ((dof > 0) && (dof != blockSize)) {
71: blockSize = 1;
72: break;
73: }
74: }
75: PetscSectionGetStorageSize(section, &localSize);
76: VecCreate(PETSC_COMM_SELF, vec);
77: VecSetSizes(*vec, localSize, localSize);
78: VecSetBlockSize(*vec, blockSize);
79: VecSetType(*vec,dm->vectype);
80: VecSetDM(*vec, dm);
81: return(0);
82: }
84: /* This assumes that the DM has been cloned prior to the call */
85: PetscErrorCode DMCreateSubDM_Section_Private(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm)
86: {
87: PetscSection section, sectionGlobal;
88: PetscInt *subIndices;
89: PetscInt subSize = 0, subOff = 0, nF, f, pStart, pEnd, p;
93: if (!numFields) return(0);
94: DMGetDefaultSection(dm, §ion);
95: DMGetDefaultGlobalSection(dm, §ionGlobal);
96: if (!section) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields");
97: if (!sectionGlobal) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields");
98: PetscSectionGetNumFields(section, &nF);
99: 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);
100: if (is) {
101: PetscInt bs = -1, bsLocal, bsMax, bsMin;
102: PetscSectionGetChart(sectionGlobal, &pStart, &pEnd);
103: for (p = pStart; p < pEnd; ++p) {
104: PetscInt gdof, pSubSize = 0;
106: PetscSectionGetDof(sectionGlobal, p, &gdof);
107: if (gdof > 0) {
108: for (f = 0; f < numFields; ++f) {
109: PetscInt fdof, fcdof;
111: PetscSectionGetFieldDof(section, p, fields[f], &fdof);
112: PetscSectionGetFieldConstraintDof(section, p, fields[f], &fcdof);
113: pSubSize += fdof-fcdof;
114: }
115: subSize += pSubSize;
116: if (pSubSize) {
117: if (bs < 0) {
118: bs = pSubSize;
119: } else if (bs != pSubSize) {
120: /* Layout does not admit a pointwise block size */
121: bs = 1;
122: }
123: }
124: }
125: }
126: /* Must have same blocksize on all procs (some might have no points) */
127: bsLocal = bs;
128: MPIU_Allreduce(&bsLocal, &bsMax, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm));
129: bsLocal = bs < 0 ? bsMax : bs;
130: MPIU_Allreduce(&bsLocal, &bsMin, 1, MPIU_INT, MPI_MIN, PetscObjectComm((PetscObject)dm));
131: if (bsMin != bsMax) {
132: bs = 1;
133: } else {
134: bs = bsMax;
135: }
136: PetscMalloc1(subSize, &subIndices);
137: for (p = pStart; p < pEnd; ++p) {
138: PetscInt gdof, goff;
140: PetscSectionGetDof(sectionGlobal, p, &gdof);
141: if (gdof > 0) {
142: PetscSectionGetOffset(sectionGlobal, p, &goff);
143: for (f = 0; f < numFields; ++f) {
144: PetscInt fdof, fcdof, fc, f2, poff = 0;
146: /* Can get rid of this loop by storing field information in the global section */
147: for (f2 = 0; f2 < fields[f]; ++f2) {
148: PetscSectionGetFieldDof(section, p, f2, &fdof);
149: PetscSectionGetFieldConstraintDof(section, p, f2, &fcdof);
150: poff += fdof-fcdof;
151: }
152: PetscSectionGetFieldDof(section, p, fields[f], &fdof);
153: PetscSectionGetFieldConstraintDof(section, p, fields[f], &fcdof);
154: for (fc = 0; fc < fdof-fcdof; ++fc, ++subOff) {
155: subIndices[subOff] = goff+poff+fc;
156: }
157: }
158: }
159: }
160: ISCreateGeneral(PetscObjectComm((PetscObject)dm), subSize, subIndices, PETSC_OWN_POINTER, is);
161: ISSetBlockSize(*is, bs);
162: }
163: if (subdm) {
164: PetscSection subsection;
165: PetscBool haveNull = PETSC_FALSE;
166: PetscInt f, nf = 0;
168: PetscSectionCreateSubsection(section, numFields, fields, &subsection);
169: DMSetDefaultSection(*subdm, subsection);
170: PetscSectionDestroy(&subsection);
171: for (f = 0; f < numFields; ++f) {
172: (*subdm)->nullspaceConstructors[f] = dm->nullspaceConstructors[fields[f]];
173: if ((*subdm)->nullspaceConstructors[f]) {
174: haveNull = PETSC_TRUE;
175: nf = f;
176: }
177: }
178: if (haveNull && is) {
179: MatNullSpace nullSpace;
181: (*(*subdm)->nullspaceConstructors[nf])(*subdm, nf, &nullSpace);
182: PetscObjectCompose((PetscObject) *is, "nullspace", (PetscObject) nullSpace);
183: MatNullSpaceDestroy(&nullSpace);
184: }
185: if (dm->prob) {
186: PetscInt Nf;
188: PetscDSGetNumFields(dm->prob, &Nf);
189: 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);
190: DMSetNumFields(*subdm, numFields);
191: for (f = 0; f < numFields; ++f) {
192: PetscObject disc;
194: DMGetField(dm, fields[f], &disc);
195: DMSetField(*subdm, f, disc);
196: }
197: if (numFields == 1 && is) {
198: PetscObject disc, space, pmat;
200: DMGetField(*subdm, 0, &disc);
201: PetscObjectQuery(disc, "nullspace", &space);
202: if (space) {PetscObjectCompose((PetscObject) *is, "nullspace", space);}
203: PetscObjectQuery(disc, "nearnullspace", &space);
204: if (space) {PetscObjectCompose((PetscObject) *is, "nearnullspace", space);}
205: PetscObjectQuery(disc, "pmat", &pmat);
206: if (pmat) {PetscObjectCompose((PetscObject) *is, "pmat", pmat);}
207: }
208: }
209: }
210: #if 0
211: /* We need a way to filter the original SF for given fields:
212: - Keeping the original section around it too much I think
213: - We could keep the distributed section, and subset it
214: */
215: if (dm->sfNatural) {
216: PetscSF sfNatural;
218: PetscSectionCreateSubsection(dm->originalSection, numFields, fields, &(*subdm)->originalSection);
219: DMPlexCreateGlobalToNaturalPetscSF(*subdm, &sfNatural);
220: DMPlexSetGlobalToNaturalPetscSF(*subdm, sfNatural);
221: }
222: #endif
223: return(0);
224: }