xref: /petsc/src/dm/impls/plex/plexnatural.c (revision 04346f8ccbfaa9fb73fb7dce6ecc10b18d1c2572)
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, globalSize;
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, &globalSize));
137   if (globalSize) {
138     const PetscInt *leaves;
139     PetscInt       *sortleaves, *indices;
140     PetscInt        Nl;
141 
142     /* Get a pruned version of migration SF */
143     PetscCall(DMGetGlobalSection(dm, &gSection));
144     if (debug) PetscCall(PetscSectionView(gSection, NULL));
145     PetscCall(PetscSFGetGraph(sfMigration, NULL, &Nl, &leaves, NULL));
146     PetscCall(PetscSectionGetChart(gSection, &pStart, &pEnd));
147     for (p = pStart, ssize = 0; p < pEnd; ++p) {
148       PetscInt dof, off;
149 
150       PetscCall(PetscSectionGetDof(gSection, p, &dof));
151       PetscCall(PetscSectionGetOffset(gSection, p, &off));
152       if ((dof > 0) && (off >= 0)) ++ssize;
153     }
154     PetscCall(PetscMalloc3(ssize, &spoints, Nl, &sortleaves, Nl, &indices));
155     for (p = 0; p < Nl; ++p) {sortleaves[p] = leaves ? leaves[p] : p; indices[p] = p;}
156     PetscCall(PetscSortIntWithArray(Nl, sortleaves, indices));
157     for (p = pStart, ssize = 0; p < pEnd; ++p) {
158       PetscInt dof, off, loc;
159 
160       PetscCall(PetscSectionGetDof(gSection, p, &dof));
161       PetscCall(PetscSectionGetOffset(gSection, p, &off));
162       if ((dof > 0) && (off >= 0)) {
163         PetscCall(PetscFindInt(p, Nl, sortleaves, &loc));
164         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);
165         spoints[ssize++] = indices[loc];
166       }
167     }
168     PetscCall(PetscSFCreateEmbeddedLeafSF(sfMigration, ssize, spoints, &sfEmbed));
169     PetscCall(PetscObjectSetName((PetscObject) sfEmbed, "Embedded SF"));
170     PetscCall(PetscFree3(spoints, sortleaves, indices));
171     if (debug) PetscCall(PetscSFView(sfEmbed, NULL));
172     /* Create the SF associated with this section
173          Roots are natural dofs, leaves are global dofs */
174     PetscCall(DMGetPointSF(dm, &sf));
175     PetscCall(PetscSectionCreateGlobalSection(sectionDist, sf, PETSC_TRUE, PETSC_TRUE, &gLocSection));
176     PetscCall(PetscSFCreateSectionSF(sfEmbed, section, remoteOffsets, gLocSection, &sfField));
177     PetscCall(PetscSFDestroy(&sfEmbed));
178     PetscCall(PetscSectionDestroy(&gLocSection));
179     PetscCall(PetscObjectSetName((PetscObject) sfField, "Natural-to-Global SF"));
180     if (debug) PetscCall(PetscSFView(sfField, NULL));
181     /* Invert the field SF
182          Roots are global dofs, leaves are natural dofs */
183     PetscCall(PetscSFCreateInverseSF(sfField, sfNatural));
184     PetscCall(PetscObjectSetName((PetscObject) *sfNatural, "Global-to-Natural SF"));
185     PetscCall(PetscObjectViewFromOptions((PetscObject) *sfNatural, NULL, "-globaltonatural_sf_view"));
186     /* Clean up */
187     PetscCall(PetscSFDestroy(&sfField));
188   } else {
189     *sfNatural = NULL;
190   }
191   PetscCall(PetscSectionDestroy(&sectionDist));
192   PetscCall(PetscFree(remoteOffsets));
193   if (destroyFlag) PetscCall(PetscSectionDestroy(&section));
194   PetscFunctionReturn(0);
195 }
196 
197 /*@
198   DMPlexGlobalToNaturalBegin - Rearranges a global Vector in the natural order.
199 
200   Collective on dm
201 
202   Input Parameters:
203 + dm - The distributed DMPlex
204 - gv - The global Vec
205 
206   Output Parameters:
207 . nv - Vec in the canonical ordering distributed over all processors associated with gv
208 
209   Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute().
210 
211   Level: intermediate
212 
213 .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalEnd()`
214 @*/
215 PetscErrorCode DMPlexGlobalToNaturalBegin(DM dm, Vec gv, Vec nv)
216 {
217   const PetscScalar *inarray;
218   PetscScalar       *outarray;
219   PetscMPIInt        size;
220 
221   PetscFunctionBegin;
222   PetscCall(PetscLogEventBegin(DMPLEX_GlobalToNaturalBegin,dm,0,0,0));
223   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size));
224   if (dm->sfNatural) {
225     PetscCall(VecGetArray(nv, &outarray));
226     PetscCall(VecGetArrayRead(gv, &inarray));
227     PetscCall(PetscSFBcastBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray,MPI_REPLACE));
228     PetscCall(VecRestoreArrayRead(gv, &inarray));
229     PetscCall(VecRestoreArray(nv, &outarray));
230   } else if (size == 1) {
231     PetscCall(VecCopy(gv, nv));
232   } else {
233     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.");
234     SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().");
235   }
236   PetscCall(PetscLogEventEnd(DMPLEX_GlobalToNaturalBegin,dm,0,0,0));
237   PetscFunctionReturn(0);
238 }
239 
240 /*@
241   DMPlexGlobalToNaturalEnd - Rearranges a global Vector in the natural order.
242 
243   Collective on dm
244 
245   Input Parameters:
246 + dm - The distributed DMPlex
247 - gv - The global Vec
248 
249   Output Parameter:
250 . nv - The natural Vec
251 
252   Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute().
253 
254   Level: intermediate
255 
256  .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()`
257  @*/
258 PetscErrorCode DMPlexGlobalToNaturalEnd(DM dm, Vec gv, Vec nv)
259 {
260   const PetscScalar *inarray;
261   PetscScalar       *outarray;
262   PetscMPIInt        size;
263 
264   PetscFunctionBegin;
265   PetscCall(PetscLogEventBegin(DMPLEX_GlobalToNaturalEnd,dm,0,0,0));
266   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size));
267   if (dm->sfNatural) {
268     PetscCall(VecGetArrayRead(gv, &inarray));
269     PetscCall(VecGetArray(nv, &outarray));
270     PetscCall(PetscSFBcastEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray,MPI_REPLACE));
271     PetscCall(VecRestoreArrayRead(gv, &inarray));
272     PetscCall(VecRestoreArray(nv, &outarray));
273   } else if (size == 1) {
274   } else {
275     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.");
276     SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().");
277   }
278   PetscCall(PetscLogEventEnd(DMPLEX_GlobalToNaturalEnd,dm,0,0,0));
279   PetscFunctionReturn(0);
280 }
281 
282 /*@
283   DMPlexNaturalToGlobalBegin - Rearranges a Vector in the natural order to the Global order.
284 
285   Collective on dm
286 
287   Input Parameters:
288 + dm - The distributed DMPlex
289 - nv - The natural Vec
290 
291   Output Parameters:
292 . gv - The global Vec
293 
294   Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute().
295 
296   Level: intermediate
297 
298 .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalEnd()`
299 @*/
300 PetscErrorCode DMPlexNaturalToGlobalBegin(DM dm, Vec nv, Vec gv)
301 {
302   const PetscScalar *inarray;
303   PetscScalar       *outarray;
304   PetscMPIInt        size;
305 
306   PetscFunctionBegin;
307   PetscCall(PetscLogEventBegin(DMPLEX_NaturalToGlobalBegin,dm,0,0,0));
308   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size));
309   if (dm->sfNatural) {
310     /* We only have access to the SF that goes from Global to Natural.
311        Instead of inverting dm->sfNatural, we can call PetscSFReduceBegin/End with MPI_Op MPI_SUM.
312        Here the SUM really does nothing since sfNatural is one to one, as long as gV is set to zero first. */
313     PetscCall(VecZeroEntries(gv));
314     PetscCall(VecGetArray(gv, &outarray));
315     PetscCall(VecGetArrayRead(nv, &inarray));
316     PetscCall(PetscSFReduceBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray, MPI_SUM));
317     PetscCall(VecRestoreArrayRead(nv, &inarray));
318     PetscCall(VecRestoreArray(gv, &outarray));
319   } else if (size == 1) {
320     PetscCall(VecCopy(nv, gv));
321   } else {
322     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.");
323     SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().");
324   }
325   PetscCall(PetscLogEventEnd(DMPLEX_NaturalToGlobalBegin,dm,0,0,0));
326   PetscFunctionReturn(0);
327 }
328 
329 /*@
330   DMPlexNaturalToGlobalEnd - Rearranges a Vector in the natural order to the Global order.
331 
332   Collective on dm
333 
334   Input Parameters:
335 + dm - The distributed DMPlex
336 - nv - The natural Vec
337 
338   Output Parameters:
339 . gv - The global Vec
340 
341   Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute().
342 
343   Level: intermediate
344 
345 .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()`
346  @*/
347 PetscErrorCode DMPlexNaturalToGlobalEnd(DM dm, Vec nv, Vec gv)
348 {
349   const PetscScalar *inarray;
350   PetscScalar       *outarray;
351   PetscMPIInt        size;
352 
353   PetscFunctionBegin;
354   PetscCall(PetscLogEventBegin(DMPLEX_NaturalToGlobalEnd,dm,0,0,0));
355   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size));
356   if (dm->sfNatural) {
357     PetscCall(VecGetArrayRead(nv, &inarray));
358     PetscCall(VecGetArray(gv, &outarray));
359     PetscCall(PetscSFReduceEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray, MPI_SUM));
360     PetscCall(VecRestoreArrayRead(nv, &inarray));
361     PetscCall(VecRestoreArray(gv, &outarray));
362   } else if (size == 1) {
363   } else {
364     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.");
365     SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().");
366   }
367   PetscCall(PetscLogEventEnd(DMPLEX_NaturalToGlobalEnd,dm,0,0,0));
368   PetscFunctionReturn(0);
369 }
370 
371 /*@
372   DMPlexCreateNaturalVector - Provide a Vec capable of holding the natural ordering and distribution.
373 
374   Collective on dm
375 
376   Input Parameters:
377 + dm - The distributed DMPlex
378 
379   Output Parameters:
380 . nv - The natural Vec
381 
382   Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute().
383 
384   Level: intermediate
385 
386 .seealso: `DMPlexDistribute()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()`
387  @*/
388 PetscErrorCode DMPlexCreateNaturalVector(DM dm, Vec *nv)
389 {
390   PetscMPIInt size;
391 
392   PetscFunctionBegin;
393   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size));
394   if (dm->sfNatural) {
395     PetscInt    nleaves, bs;
396     Vec         v;
397     PetscCall(DMGetLocalVector(dm,&v));
398     PetscCall(VecGetBlockSize(v,&bs));
399     PetscCall(DMRestoreLocalVector(dm,&v));
400 
401     PetscCall(PetscSFGetGraph(dm->sfNatural,NULL,&nleaves,NULL,NULL));
402     PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), nv));
403     PetscCall(VecSetSizes(*nv, nleaves, PETSC_DETERMINE));
404     PetscCall(VecSetBlockSize(*nv, bs));
405     PetscCall(VecSetType(*nv,dm->vectype));
406     PetscCall(VecSetDM(*nv, dm));
407   } else if (size == 1) {
408     PetscCall(DMCreateLocalVector(dm,nv));
409   } else {
410     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.");
411     SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().");
412   }
413   PetscFunctionReturn(0);
414 }
415