Actual source code: dmi.c
petsc-3.12.5 2020-03-29
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: DMGetGlobalSection(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: DMGetLocalSection(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: /*@C
85: DMCreateSectionSubDM - Returns an IS and subDM+subSection encapsulating a subproblem defined by the fields in a PetscSection in the DM.
87: Not collective
89: Input Parameters:
90: + dm - The DM object
91: . numFields - The number of fields in this subproblem
92: - fields - The field numbers of the selected fields
94: Output Parameters:
95: + is - The global indices for the subproblem
96: - subdm - The DM for the subproblem, which must already have be cloned from dm
98: Note: This handles all information in the DM class and the PetscSection. This is used as the basis for creating subDMs in specialized classes,
99: such as Plex and Forest.
101: Level: intermediate
103: .seealso DMCreateSubDM(), DMGetLocalSection(), DMPlexSetMigrationSF(), DMView()
104: @*/
105: PetscErrorCode DMCreateSectionSubDM(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm)
106: {
107: PetscSection section, sectionGlobal;
108: PetscInt *subIndices;
109: PetscInt subSize = 0, subOff = 0, Nf, f, pStart, pEnd, p;
113: if (!numFields) return(0);
114: DMGetLocalSection(dm, §ion);
115: DMGetGlobalSection(dm, §ionGlobal);
116: if (!section) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields");
117: if (!sectionGlobal) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields");
118: PetscSectionGetNumFields(section, &Nf);
119: 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);
120: if (is) {
121: PetscInt bs, bsLocal[2], bsMinMax[2];
123: for (f = 0, bs = 0; f < numFields; ++f) {
124: PetscInt Nc;
126: PetscSectionGetFieldComponents(section, fields[f], &Nc);
127: bs += Nc;
128: }
129: PetscSectionGetChart(sectionGlobal, &pStart, &pEnd);
130: for (p = pStart; p < pEnd; ++p) {
131: PetscInt gdof, pSubSize = 0;
133: PetscSectionGetDof(sectionGlobal, p, &gdof);
134: if (gdof > 0) {
135: for (f = 0; f < numFields; ++f) {
136: PetscInt fdof, fcdof;
138: PetscSectionGetFieldDof(section, p, fields[f], &fdof);
139: PetscSectionGetFieldConstraintDof(section, p, fields[f], &fcdof);
140: pSubSize += fdof-fcdof;
141: }
142: subSize += pSubSize;
143: if (pSubSize && bs != pSubSize) {
144: /* Layout does not admit a pointwise block size */
145: bs = 1;
146: }
147: }
148: }
149: /* Must have same blocksize on all procs (some might have no points) */
150: bsLocal[0] = bs < 0 ? PETSC_MAX_INT : bs; bsLocal[1] = bs;
151: PetscGlobalMinMaxInt(PetscObjectComm((PetscObject) dm), bsLocal, bsMinMax);
152: if (bsMinMax[0] != bsMinMax[1]) {bs = 1;}
153: else {bs = bsMinMax[0];}
154: PetscMalloc1(subSize, &subIndices);
155: for (p = pStart; p < pEnd; ++p) {
156: PetscInt gdof, goff;
158: PetscSectionGetDof(sectionGlobal, p, &gdof);
159: if (gdof > 0) {
160: PetscSectionGetOffset(sectionGlobal, p, &goff);
161: for (f = 0; f < numFields; ++f) {
162: PetscInt fdof, fcdof, fc, f2, poff = 0;
164: /* Can get rid of this loop by storing field information in the global section */
165: for (f2 = 0; f2 < fields[f]; ++f2) {
166: PetscSectionGetFieldDof(section, p, f2, &fdof);
167: PetscSectionGetFieldConstraintDof(section, p, f2, &fcdof);
168: poff += fdof-fcdof;
169: }
170: PetscSectionGetFieldDof(section, p, fields[f], &fdof);
171: PetscSectionGetFieldConstraintDof(section, p, fields[f], &fcdof);
172: for (fc = 0; fc < fdof-fcdof; ++fc, ++subOff) {
173: subIndices[subOff] = goff+poff+fc;
174: }
175: }
176: }
177: }
178: ISCreateGeneral(PetscObjectComm((PetscObject)dm), subSize, subIndices, PETSC_OWN_POINTER, is);
179: if (bs > 1) {
180: /* We need to check that the block size does not come from non-contiguous fields */
181: PetscInt i, j, set = 1;
182: for (i = 0; i < subSize; i += bs) {
183: for (j = 0; j < bs; ++j) {
184: if (subIndices[i+j] != subIndices[i]+j) {set = 0; break;}
185: }
186: }
187: if (set) {ISSetBlockSize(*is, bs);}
188: }
189: }
190: if (subdm) {
191: PetscSection subsection;
192: PetscBool haveNull = PETSC_FALSE;
193: PetscInt f, nf = 0;
195: PetscSectionCreateSubsection(section, numFields, fields, &subsection);
196: DMSetLocalSection(*subdm, subsection);
197: PetscSectionDestroy(&subsection);
198: for (f = 0; f < numFields; ++f) {
199: (*subdm)->nullspaceConstructors[f] = dm->nullspaceConstructors[fields[f]];
200: if ((*subdm)->nullspaceConstructors[f]) {
201: haveNull = PETSC_TRUE;
202: nf = f;
203: }
204: }
205: if (haveNull && is) {
206: MatNullSpace nullSpace;
208: (*(*subdm)->nullspaceConstructors[nf])(*subdm, nf, &nullSpace);
209: PetscObjectCompose((PetscObject) *is, "nullspace", (PetscObject) nullSpace);
210: MatNullSpaceDestroy(&nullSpace);
211: }
212: if (dm->probs) {
213: DMSetNumFields(*subdm, numFields);
214: for (f = 0; f < numFields; ++f) {
215: PetscObject disc;
217: DMGetField(dm, fields[f], NULL, &disc);
218: DMSetField(*subdm, f, NULL, disc);
219: }
220: DMCreateDS(*subdm);
221: if (numFields == 1 && is) {
222: PetscObject disc, space, pmat;
224: DMGetField(*subdm, 0, NULL, &disc);
225: PetscObjectQuery(disc, "nullspace", &space);
226: if (space) {PetscObjectCompose((PetscObject) *is, "nullspace", space);}
227: PetscObjectQuery(disc, "nearnullspace", &space);
228: if (space) {PetscObjectCompose((PetscObject) *is, "nearnullspace", space);}
229: PetscObjectQuery(disc, "pmat", &pmat);
230: if (pmat) {PetscObjectCompose((PetscObject) *is, "pmat", pmat);}
231: }
232: PetscDSCopyConstants(dm->probs[0].ds, (*subdm)->probs[0].ds);
233: PetscDSCopyBoundary(dm->probs[0].ds, (*subdm)->probs[0].ds);
234: /* Translate DM fields to DS fields */
235: if (dm->probs[0].fields) {
236: IS infields, dsfields;
237: const PetscInt *fld, *ofld;
238: PetscInt *fidx;
239: PetscInt onf, nf, f, g;
241: ISCreateGeneral(PETSC_COMM_SELF, numFields, fields, PETSC_USE_POINTER, &infields);
242: ISIntersect(infields, dm->probs[0].fields, &dsfields);
243: ISDestroy(&infields);
244: ISGetLocalSize(dsfields, &nf);
245: if (!nf) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "DS cannot be supported on 0 fields");
246: ISGetIndices(dsfields, &fld);
247: ISGetLocalSize(dm->probs[0].fields, &onf);
248: ISGetIndices(dm->probs[0].fields, &ofld);
249: PetscMalloc1(nf, &fidx);
250: for (f = 0, g = 0; f < onf && g < nf; ++f) {
251: if (ofld[f] == fld[g]) fidx[g++] = f;
252: }
253: ISRestoreIndices(dm->probs[0].fields, &ofld);
254: ISRestoreIndices(dsfields, &fld);
255: ISDestroy(&dsfields);
256: PetscDSSelectEquations(dm->probs[0].ds, nf, fidx, (*subdm)->probs[0].ds);
257: PetscFree(fidx);
258: } else {
259: PetscDSSelectEquations(dm->probs[0].ds, numFields, fields, (*subdm)->probs[0].ds);
260: }
261: }
262: if (dm->coarseMesh) {
263: DMCreateSubDM(dm->coarseMesh, numFields, fields, NULL, &(*subdm)->coarseMesh);
264: }
265: }
266: return(0);
267: }
269: /*@C
270: DMCreateSectionSuperDM - Returns an arrays of ISes and DM+Section encapsulating a superproblem defined by the DM+Sections passed in.
272: Not collective
274: Input Parameter:
275: + dms - The DM objects
276: - len - The number of DMs
278: Output Parameters:
279: + is - The global indices for the subproblem, or NULL
280: - superdm - The DM for the superproblem, which must already have be cloned
282: Note: This handles all information in the DM class and the PetscSection. This is used as the basis for creating subDMs in specialized classes,
283: such as Plex and Forest.
285: Level: intermediate
287: .seealso DMCreateSuperDM(), DMGetLocalSection(), DMPlexSetMigrationSF(), DMView()
288: @*/
289: PetscErrorCode DMCreateSectionSuperDM(DM dms[], PetscInt len, IS **is, DM *superdm)
290: {
291: MPI_Comm comm;
292: PetscSection supersection, *sections, *sectionGlobals;
293: PetscInt *Nfs, Nf = 0, f, supf, nullf = -1, i;
294: PetscBool haveNull = PETSC_FALSE;
298: PetscObjectGetComm((PetscObject)dms[0], &comm);
299: /* Pull out local and global sections */
300: PetscMalloc3(len, &Nfs, len, §ions, len, §ionGlobals);
301: for (i = 0 ; i < len; ++i) {
302: DMGetLocalSection(dms[i], §ions[i]);
303: DMGetGlobalSection(dms[i], §ionGlobals[i]);
304: if (!sections[i]) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields");
305: if (!sectionGlobals[i]) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields");
306: PetscSectionGetNumFields(sections[i], &Nfs[i]);
307: Nf += Nfs[i];
308: }
309: /* Create the supersection */
310: PetscSectionCreateSupersection(sections, len, &supersection);
311: DMSetLocalSection(*superdm, supersection);
312: /* Create ISes */
313: if (is) {
314: PetscSection supersectionGlobal;
315: PetscInt bs = -1, startf = 0;
317: PetscMalloc1(len, is);
318: DMGetGlobalSection(*superdm, &supersectionGlobal);
319: for (i = 0 ; i < len; startf += Nfs[i], ++i) {
320: PetscInt *subIndices;
321: PetscInt subSize, subOff, pStart, pEnd, p, start, end, dummy;
323: PetscSectionGetChart(sectionGlobals[i], &pStart, &pEnd);
324: PetscSectionGetConstrainedStorageSize(sectionGlobals[i], &subSize);
325: PetscMalloc1(subSize, &subIndices);
326: for (p = pStart, subOff = 0; p < pEnd; ++p) {
327: PetscInt gdof, gcdof, gtdof, d;
329: PetscSectionGetDof(sectionGlobals[i], p, &gdof);
330: PetscSectionGetConstraintDof(sections[i], p, &gcdof);
331: gtdof = gdof - gcdof;
332: if (gdof > 0 && gtdof) {
333: if (bs < 0) {bs = gtdof;}
334: else if (bs != gtdof) {bs = 1;}
335: DMGetGlobalFieldOffset_Private(*superdm, p, startf, &start, &dummy);
336: DMGetGlobalFieldOffset_Private(*superdm, p, startf+Nfs[i]-1, &dummy, &end);
337: if (end-start != gtdof) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid number of global dofs %D != %D for dm %D on point %D", end-start, gtdof, i, p);
338: for (d = start; d < end; ++d, ++subOff) subIndices[subOff] = d;
339: }
340: }
341: ISCreateGeneral(comm, subSize, subIndices, PETSC_OWN_POINTER, &(*is)[i]);
342: /* Must have same blocksize on all procs (some might have no points) */
343: {
344: PetscInt bs = -1, bsLocal[2], bsMinMax[2];
346: bsLocal[0] = bs < 0 ? PETSC_MAX_INT : bs; bsLocal[1] = bs;
347: PetscGlobalMinMaxInt(comm, bsLocal, bsMinMax);
348: if (bsMinMax[0] != bsMinMax[1]) {bs = 1;}
349: else {bs = bsMinMax[0];}
350: ISSetBlockSize((*is)[i], bs);
351: }
352: }
353: }
354: /* Preserve discretizations */
355: if (len && dms[0]->probs) {
356: DMSetNumFields(*superdm, Nf);
357: for (i = 0, supf = 0; i < len; ++i) {
358: for (f = 0; f < Nfs[i]; ++f, ++supf) {
359: PetscObject disc;
361: DMGetField(dms[i], f, NULL, &disc);
362: DMSetField(*superdm, supf, NULL, disc);
363: }
364: }
365: DMCreateDS(*superdm);
366: }
367: /* Preserve nullspaces */
368: for (i = 0, supf = 0; i < len; ++i) {
369: for (f = 0; f < Nfs[i]; ++f, ++supf) {
370: (*superdm)->nullspaceConstructors[supf] = dms[i]->nullspaceConstructors[f];
371: if ((*superdm)->nullspaceConstructors[supf]) {
372: haveNull = PETSC_TRUE;
373: nullf = supf;
374: }
375: }
376: }
377: /* Attach nullspace to IS */
378: if (haveNull && is) {
379: MatNullSpace nullSpace;
381: (*(*superdm)->nullspaceConstructors[nullf])(*superdm, nullf, &nullSpace);
382: PetscObjectCompose((PetscObject) (*is)[nullf], "nullspace", (PetscObject) nullSpace);
383: MatNullSpaceDestroy(&nullSpace);
384: }
385: PetscSectionDestroy(&supersection);
386: PetscFree3(Nfs, sections, sectionGlobals);
387: return(0);
388: }