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