Actual source code: plexnatural.c
1: #include <petsc/private/dmpleximpl.h>
3: /*@
4: DMPlexSetMigrationSF - Sets the SF for migrating from a parent DM into this DM
6: Input Parameters:
7: + dm - The DM
8: - naturalSF - The PetscSF
10: Note: It is necessary to call this in order to have DMCreateSubDM() or DMCreateSuperDM() build the Global-To-Natural map
12: Level: intermediate
14: .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexCreateMigrationSF()`, `DMPlexGetMigrationSF()`
15: @*/
16: PetscErrorCode DMPlexSetMigrationSF(DM dm, PetscSF migrationSF)
17: {
18: dm->sfMigration = migrationSF;
19: PetscObjectReference((PetscObject)migrationSF);
20: return 0;
21: }
23: /*@
24: DMPlexGetMigrationSF - Gets the SF for migrating from a parent DM into this DM
26: Input Parameter:
27: . dm - The DM
29: Output Parameter:
30: . migrationSF - The PetscSF
32: Level: intermediate
34: .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexCreateMigrationSF()`, `DMPlexSetMigrationSF`
35: @*/
36: PetscErrorCode DMPlexGetMigrationSF(DM dm, PetscSF *migrationSF)
37: {
38: *migrationSF = dm->sfMigration;
39: return 0;
40: }
42: /*@
43: DMPlexSetGlobalToNaturalSF - Sets the SF for mapping Global Vec to the Natural Vec
45: Input Parameters:
46: + dm - The DM
47: - naturalSF - The PetscSF
49: Level: intermediate
51: .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexCreateGlobalToNaturalSF()`, `DMPlexGetGlobaltoNaturalSF()`
52: @*/
53: PetscErrorCode DMPlexSetGlobalToNaturalSF(DM dm, PetscSF naturalSF)
54: {
55: dm->sfNatural = naturalSF;
56: PetscObjectReference((PetscObject)naturalSF);
57: dm->useNatural = PETSC_TRUE;
58: return 0;
59: }
61: /*@
62: DMPlexGetGlobalToNaturalSF - Gets the SF for mapping Global Vec to the Natural Vec
64: Input Parameter:
65: . dm - The DM
67: Output Parameter:
68: . naturalSF - The PetscSF
70: Level: intermediate
72: .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexCreateGlobalToNaturalSF()`, `DMPlexSetGlobaltoNaturalSF`
73: @*/
74: PetscErrorCode DMPlexGetGlobalToNaturalSF(DM dm, PetscSF *naturalSF)
75: {
76: *naturalSF = dm->sfNatural;
77: return 0;
78: }
80: /*@
81: DMPlexCreateGlobalToNaturalSF - Creates the SF for mapping Global Vec to the Natural Vec
83: Input Parameters:
84: + dm - The redistributed DM
85: . section - The local PetscSection describing the Vec before the mesh was distributed, or NULL if not available
86: - sfMigration - The PetscSF used to distribute the mesh, or NULL if it cannot be computed
88: Output Parameter:
89: . sfNatural - PetscSF for mapping the Vec in PETSc ordering to the canonical ordering
91: Note: This is not typically called by the user.
93: Level: intermediate
95: .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`
96: @*/
97: PetscErrorCode DMPlexCreateGlobalToNaturalSF(DM dm, PetscSection section, PetscSF sfMigration, PetscSF *sfNatural)
98: {
99: MPI_Comm comm;
100: PetscSF sf, sfEmbed, sfField;
101: PetscSection gSection, sectionDist, gLocSection;
102: PetscInt *spoints, *remoteOffsets;
103: PetscInt ssize, pStart, pEnd, p, globalSize;
104: PetscBool destroyFlag = PETSC_FALSE, debug = PETSC_FALSE;
106: PetscObjectGetComm((PetscObject)dm, &comm);
107: if (!sfMigration) {
108: /* If sfMigration is missing, sfNatural cannot be computed and is set to NULL */
109: *sfNatural = NULL;
110: return 0;
111: } else if (!section) {
112: /* If the sequential section is not provided (NULL), it is reconstructed from the parallel section */
113: PetscSF sfMigrationInv;
114: PetscSection localSection;
116: DMGetLocalSection(dm, &localSection);
117: PetscSFCreateInverseSF(sfMigration, &sfMigrationInv);
118: PetscSectionCreate(PetscObjectComm((PetscObject)dm), §ion);
119: PetscSFDistributeSection(sfMigrationInv, localSection, NULL, section);
120: PetscSFDestroy(&sfMigrationInv);
121: destroyFlag = PETSC_TRUE;
122: }
123: if (debug) PetscSFView(sfMigration, NULL);
124: /* Create a new section from distributing the original section */
125: PetscSectionCreate(comm, §ionDist);
126: PetscSFDistributeSection(sfMigration, section, &remoteOffsets, sectionDist);
127: PetscObjectSetName((PetscObject)sectionDist, "Migrated Section");
128: if (debug) PetscSectionView(sectionDist, NULL);
129: DMSetLocalSection(dm, sectionDist);
130: /* If a sequential section is provided but no dof is affected, sfNatural cannot be computed and is set to NULL */
131: PetscSectionGetStorageSize(sectionDist, &globalSize);
132: if (globalSize) {
133: const PetscInt *leaves;
134: PetscInt *sortleaves, *indices;
135: PetscInt Nl;
137: /* Get a pruned version of migration SF */
138: DMGetGlobalSection(dm, &gSection);
139: if (debug) PetscSectionView(gSection, NULL);
140: PetscSFGetGraph(sfMigration, NULL, &Nl, &leaves, NULL);
141: PetscSectionGetChart(gSection, &pStart, &pEnd);
142: for (p = pStart, ssize = 0; p < pEnd; ++p) {
143: PetscInt dof, off;
145: PetscSectionGetDof(gSection, p, &dof);
146: PetscSectionGetOffset(gSection, p, &off);
147: if ((dof > 0) && (off >= 0)) ++ssize;
148: }
149: PetscMalloc3(ssize, &spoints, Nl, &sortleaves, Nl, &indices);
150: for (p = 0; p < Nl; ++p) {
151: sortleaves[p] = leaves ? leaves[p] : p;
152: indices[p] = p;
153: }
154: PetscSortIntWithArray(Nl, sortleaves, indices);
155: for (p = pStart, ssize = 0; p < pEnd; ++p) {
156: PetscInt dof, off, loc;
158: PetscSectionGetDof(gSection, p, &dof);
159: PetscSectionGetOffset(gSection, p, &off);
160: if ((dof > 0) && (off >= 0)) {
161: PetscFindInt(p, Nl, sortleaves, &loc);
163: spoints[ssize++] = indices[loc];
164: }
165: }
166: PetscSFCreateEmbeddedLeafSF(sfMigration, ssize, spoints, &sfEmbed);
167: PetscObjectSetName((PetscObject)sfEmbed, "Embedded SF");
168: PetscFree3(spoints, sortleaves, indices);
169: if (debug) PetscSFView(sfEmbed, NULL);
170: /* Create the SF associated with this section
171: Roots are natural dofs, leaves are global dofs */
172: DMGetPointSF(dm, &sf);
173: PetscSectionCreateGlobalSection(sectionDist, sf, PETSC_TRUE, PETSC_TRUE, &gLocSection);
174: PetscSFCreateSectionSF(sfEmbed, section, remoteOffsets, gLocSection, &sfField);
175: PetscSFDestroy(&sfEmbed);
176: PetscSectionDestroy(&gLocSection);
177: PetscObjectSetName((PetscObject)sfField, "Natural-to-Global SF");
178: if (debug) PetscSFView(sfField, NULL);
179: /* Invert the field SF
180: Roots are global dofs, leaves are natural dofs */
181: PetscSFCreateInverseSF(sfField, sfNatural);
182: PetscObjectSetName((PetscObject)*sfNatural, "Global-to-Natural SF");
183: PetscObjectViewFromOptions((PetscObject)*sfNatural, NULL, "-globaltonatural_sf_view");
184: /* Clean up */
185: PetscSFDestroy(&sfField);
186: } else {
187: *sfNatural = NULL;
188: }
189: PetscSectionDestroy(§ionDist);
190: PetscFree(remoteOffsets);
191: if (destroyFlag) PetscSectionDestroy(§ion);
192: return 0;
193: }
195: /*@
196: DMPlexGlobalToNaturalBegin - Rearranges a global Vector in the natural order.
198: Collective on dm
200: Input Parameters:
201: + dm - The distributed DMPlex
202: - gv - The global Vec
204: Output Parameters:
205: . nv - Vec in the canonical ordering distributed over all processors associated with gv
207: Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute().
209: Level: intermediate
211: .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalEnd()`
212: @*/
213: PetscErrorCode DMPlexGlobalToNaturalBegin(DM dm, Vec gv, Vec nv)
214: {
215: const PetscScalar *inarray;
216: PetscScalar *outarray;
217: PetscMPIInt size;
219: PetscLogEventBegin(DMPLEX_GlobalToNaturalBegin, dm, 0, 0, 0);
220: MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size);
221: if (dm->sfNatural) {
222: VecGetArray(nv, &outarray);
223: VecGetArrayRead(gv, &inarray);
224: PetscSFBcastBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_REPLACE);
225: VecRestoreArrayRead(gv, &inarray);
226: VecRestoreArray(nv, &outarray);
227: } else if (size == 1) {
228: VecCopy(gv, nv);
229: } else {
231: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().");
232: }
233: PetscLogEventEnd(DMPLEX_GlobalToNaturalBegin, dm, 0, 0, 0);
234: return 0;
235: }
237: /*@
238: DMPlexGlobalToNaturalEnd - Rearranges a global Vector in the natural order.
240: Collective on dm
242: Input Parameters:
243: + dm - The distributed DMPlex
244: - gv - The global Vec
246: Output Parameter:
247: . nv - The natural Vec
249: Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute().
251: Level: intermediate
253: .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()`
254: @*/
255: PetscErrorCode DMPlexGlobalToNaturalEnd(DM dm, Vec gv, Vec nv)
256: {
257: const PetscScalar *inarray;
258: PetscScalar *outarray;
259: PetscMPIInt size;
261: PetscLogEventBegin(DMPLEX_GlobalToNaturalEnd, dm, 0, 0, 0);
262: MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size);
263: if (dm->sfNatural) {
264: VecGetArrayRead(gv, &inarray);
265: VecGetArray(nv, &outarray);
266: PetscSFBcastEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_REPLACE);
267: VecRestoreArrayRead(gv, &inarray);
268: VecRestoreArray(nv, &outarray);
269: } else if (size == 1) {
270: } else {
272: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().");
273: }
274: PetscLogEventEnd(DMPLEX_GlobalToNaturalEnd, dm, 0, 0, 0);
275: return 0;
276: }
278: /*@
279: DMPlexNaturalToGlobalBegin - Rearranges a Vector in the natural order to the Global order.
281: Collective on dm
283: Input Parameters:
284: + dm - The distributed DMPlex
285: - nv - The natural Vec
287: Output Parameters:
288: . gv - The global Vec
290: Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute().
292: Level: intermediate
294: .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalEnd()`
295: @*/
296: PetscErrorCode DMPlexNaturalToGlobalBegin(DM dm, Vec nv, Vec gv)
297: {
298: const PetscScalar *inarray;
299: PetscScalar *outarray;
300: PetscMPIInt size;
302: PetscLogEventBegin(DMPLEX_NaturalToGlobalBegin, dm, 0, 0, 0);
303: MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size);
304: if (dm->sfNatural) {
305: /* We only have access to the SF that goes from Global to Natural.
306: Instead of inverting dm->sfNatural, we can call PetscSFReduceBegin/End with MPI_Op MPI_SUM.
307: Here the SUM really does nothing since sfNatural is one to one, as long as gV is set to zero first. */
308: VecZeroEntries(gv);
309: VecGetArray(gv, &outarray);
310: VecGetArrayRead(nv, &inarray);
311: PetscSFReduceBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_SUM);
312: VecRestoreArrayRead(nv, &inarray);
313: VecRestoreArray(gv, &outarray);
314: } else if (size == 1) {
315: VecCopy(nv, gv);
316: } else {
318: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().");
319: }
320: PetscLogEventEnd(DMPLEX_NaturalToGlobalBegin, dm, 0, 0, 0);
321: return 0;
322: }
324: /*@
325: DMPlexNaturalToGlobalEnd - Rearranges a Vector in the natural order to the Global order.
327: Collective on dm
329: Input Parameters:
330: + dm - The distributed DMPlex
331: - nv - The natural Vec
333: Output Parameters:
334: . gv - The global Vec
336: Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute().
338: Level: intermediate
340: .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()`
341: @*/
342: PetscErrorCode DMPlexNaturalToGlobalEnd(DM dm, Vec nv, Vec gv)
343: {
344: const PetscScalar *inarray;
345: PetscScalar *outarray;
346: PetscMPIInt size;
348: PetscLogEventBegin(DMPLEX_NaturalToGlobalEnd, dm, 0, 0, 0);
349: MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size);
350: if (dm->sfNatural) {
351: VecGetArrayRead(nv, &inarray);
352: VecGetArray(gv, &outarray);
353: PetscSFReduceEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_SUM);
354: VecRestoreArrayRead(nv, &inarray);
355: VecRestoreArray(gv, &outarray);
356: } else if (size == 1) {
357: } else {
359: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().");
360: }
361: PetscLogEventEnd(DMPLEX_NaturalToGlobalEnd, dm, 0, 0, 0);
362: return 0;
363: }
365: /*@
366: DMPlexCreateNaturalVector - Provide a `Vec` capable of holding the natural ordering and distribution.
368: Collective on dm
370: Input Parameter:
371: . dm - The distributed `DMPlex`
373: Output Parameter:
374: . nv - The natural `Vec`
376: Note: The user must call `DMSetUseNatural`(dm, PETSC_TRUE) before `DMPlexDistribute()`.
378: Level: intermediate
380: .seealso: `DMPlexDistribute()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()`
381: @*/
382: PetscErrorCode DMPlexCreateNaturalVector(DM dm, Vec *nv)
383: {
384: PetscMPIInt size;
386: MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size);
387: if (dm->sfNatural) {
388: PetscInt nleaves, bs;
389: Vec v;
390: DMGetLocalVector(dm, &v);
391: VecGetBlockSize(v, &bs);
392: DMRestoreLocalVector(dm, &v);
394: PetscSFGetGraph(dm->sfNatural, NULL, &nleaves, NULL, NULL);
395: VecCreate(PetscObjectComm((PetscObject)dm), nv);
396: VecSetSizes(*nv, nleaves, PETSC_DETERMINE);
397: VecSetBlockSize(*nv, bs);
398: VecSetType(*nv, dm->vectype);
399: VecSetDM(*nv, dm);
400: } else if (size == 1) {
401: DMCreateLocalVector(dm, nv);
402: } else {
404: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().");
405: }
406: return 0;
407: }