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