xref: /petsc/src/dm/impls/plex/plexnatural.c (revision e535cce4c56ea8fc2df340b7d14fa8f4f30804fb)
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
7 
8   Input Parameters:
9 + dm          - The `DM`
10 - migrationSF - 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: [](ch_unstructured), `DM`, `DMPLEX`, `PetscSF`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexCreateMigrationSF()`, `DMPlexGetMigrationSF()`
18 @*/
19 PetscErrorCode DMPlexSetMigrationSF(DM dm, PetscSF migrationSF)
20 {
21   PetscFunctionBegin;
22   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
23   if (migrationSF) PetscValidHeaderSpecific(migrationSF, PETSCSF_CLASSID, 2);
24   PetscCall(PetscObjectReference((PetscObject)migrationSF));
25   PetscCall(PetscSFDestroy(&dm->sfMigration));
26   dm->sfMigration = migrationSF;
27   PetscFunctionReturn(PETSC_SUCCESS);
28 }
29 
30 /*@
31   DMPlexGetMigrationSF - Gets the `PetscSF` for migrating from a parent `DM` into this `DM`
32 
33   Note Collective
34 
35   Input Parameter:
36 . dm - The `DM`
37 
38   Output Parameter:
39 . migrationSF - The `PetscSF`
40 
41   Level: intermediate
42 
43 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `PetscSF`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexCreateMigrationSF()`, `DMPlexSetMigrationSF`
44 @*/
45 PetscErrorCode DMPlexGetMigrationSF(DM dm, PetscSF *migrationSF)
46 {
47   PetscFunctionBegin;
48   *migrationSF = dm->sfMigration;
49   PetscFunctionReturn(PETSC_SUCCESS);
50 }
51 
52 /*@
53   DMPlexCreateGlobalToNaturalSF - Creates the `PetscSF` for mapping Global `Vec` to the Natural `Vec`
54 
55   Input Parameters:
56 + dm          - The redistributed `DM`
57 . section     - The local `PetscSection` describing the `Vec` before the mesh was distributed, or `NULL` if not available
58 - sfMigration - The `PetscSF` used to distribute the mesh, or `NULL` if it cannot be computed
59 
60   Output Parameter:
61 . sfNatural - `PetscSF` for mapping the `Vec` in PETSc ordering to the canonical ordering
62 
63   Level: intermediate
64 
65   Note:
66   This is not typically called by the user.
67 
68 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `PetscSF`, `PetscSection`, `DMPlexDistribute()`, `DMPlexDistributeField()`
69  @*/
70 PetscErrorCode DMPlexCreateGlobalToNaturalSF(DM dm, PetscSection section, PetscSF sfMigration, PetscSF *sfNatural)
71 {
72   MPI_Comm     comm;
73   PetscSF      sf, sfEmbed, sfField;
74   PetscSection gSection, sectionDist, gLocSection;
75   PetscInt    *spoints, *remoteOffsets;
76   PetscInt     ssize, pStart, pEnd, p, localSize, maxStorageSize;
77   PetscBool    destroyFlag = PETSC_FALSE, debug = PETSC_FALSE;
78 
79   PetscFunctionBegin;
80   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
81   if (!sfMigration) {
82     /* If sfMigration is missing, sfNatural cannot be computed and is set to NULL */
83     *sfNatural = NULL;
84     PetscFunctionReturn(PETSC_SUCCESS);
85   } else if (!section) {
86     /* If the sequential section is not provided (NULL), it is reconstructed from the parallel section */
87     PetscSF      sfMigrationInv;
88     PetscSection localSection;
89 
90     PetscCall(DMGetLocalSection(dm, &localSection));
91     PetscCall(PetscSFCreateInverseSF(sfMigration, &sfMigrationInv));
92     PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), &section));
93     PetscCall(PetscSFDistributeSection(sfMigrationInv, localSection, NULL, section));
94     PetscCall(PetscSFDestroy(&sfMigrationInv));
95     destroyFlag = PETSC_TRUE;
96   }
97   if (debug) PetscCall(PetscSFView(sfMigration, NULL));
98   /* Create a new section from distributing the original section */
99   PetscCall(PetscSectionCreate(comm, &sectionDist));
100   PetscCall(PetscSFDistributeSection(sfMigration, section, &remoteOffsets, sectionDist));
101   PetscCall(PetscObjectSetName((PetscObject)sectionDist, "Migrated Section"));
102   if (debug) PetscCall(PetscSectionView(sectionDist, NULL));
103   PetscCall(DMSetLocalSection(dm, sectionDist));
104   /* If a sequential section is provided but no dof is affected, sfNatural cannot be computed and is set to NULL */
105   PetscCall(PetscSectionGetStorageSize(sectionDist, &localSize));
106   PetscCallMPI(MPIU_Allreduce(&localSize, &maxStorageSize, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm)));
107   if (maxStorageSize) {
108     const PetscInt *leaves;
109     PetscInt       *sortleaves, *indices;
110     PetscInt        Nl;
111 
112     /* Get a pruned version of migration SF */
113     PetscCall(DMGetGlobalSection(dm, &gSection));
114     if (debug) PetscCall(PetscSectionView(gSection, NULL));
115     PetscCall(PetscSFGetGraph(sfMigration, NULL, &Nl, &leaves, NULL));
116     PetscCall(PetscSectionGetChart(gSection, &pStart, &pEnd));
117     for (p = pStart, ssize = 0; p < pEnd; ++p) {
118       PetscInt dof, off;
119 
120       PetscCall(PetscSectionGetDof(gSection, p, &dof));
121       PetscCall(PetscSectionGetOffset(gSection, p, &off));
122       if ((dof > 0) && (off >= 0)) ++ssize;
123     }
124     PetscCall(PetscMalloc3(ssize, &spoints, Nl, &sortleaves, Nl, &indices));
125     for (p = 0; p < Nl; ++p) {
126       sortleaves[p] = leaves ? leaves[p] : p;
127       indices[p]    = p;
128     }
129     PetscCall(PetscSortIntWithArray(Nl, sortleaves, indices));
130     for (p = pStart, ssize = 0; p < pEnd; ++p) {
131       PetscInt dof, off, loc;
132 
133       PetscCall(PetscSectionGetDof(gSection, p, &dof));
134       PetscCall(PetscSectionGetOffset(gSection, p, &off));
135       if ((dof > 0) && (off >= 0)) {
136         PetscCall(PetscFindInt(p, Nl, sortleaves, &loc));
137         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);
138         spoints[ssize++] = indices[loc];
139       }
140     }
141     PetscCall(PetscSFCreateEmbeddedLeafSF(sfMigration, ssize, spoints, &sfEmbed));
142     PetscCall(PetscObjectSetName((PetscObject)sfEmbed, "Embedded SF"));
143     PetscCall(PetscFree3(spoints, sortleaves, indices));
144     if (debug) PetscCall(PetscSFView(sfEmbed, NULL));
145     /* Create the SF associated with this section
146          Roots are natural dofs, leaves are global dofs */
147     PetscCall(DMGetPointSF(dm, &sf));
148     PetscCall(PetscSectionCreateGlobalSection(sectionDist, sf, PETSC_TRUE, PETSC_TRUE, PETSC_TRUE, &gLocSection));
149     PetscCall(PetscSFCreateSectionSF(sfEmbed, section, remoteOffsets, gLocSection, &sfField));
150     PetscCall(PetscSFDestroy(&sfEmbed));
151     PetscCall(PetscSectionDestroy(&gLocSection));
152     PetscCall(PetscObjectSetName((PetscObject)sfField, "Natural-to-Global SF"));
153     if (debug) PetscCall(PetscSFView(sfField, NULL));
154     /* Invert the field SF
155          Roots are global dofs, leaves are natural dofs */
156     PetscCall(PetscSFCreateInverseSF(sfField, sfNatural));
157     PetscCall(PetscObjectSetName((PetscObject)*sfNatural, "Global-to-Natural SF"));
158     PetscCall(PetscObjectViewFromOptions((PetscObject)*sfNatural, NULL, "-globaltonatural_sf_view"));
159     /* Clean up */
160     PetscCall(PetscSFDestroy(&sfField));
161   } else {
162     *sfNatural = NULL;
163   }
164   PetscCall(PetscSectionDestroy(&sectionDist));
165   PetscCall(PetscFree(remoteOffsets));
166   if (destroyFlag) PetscCall(PetscSectionDestroy(&section));
167   PetscFunctionReturn(PETSC_SUCCESS);
168 }
169 
170 /*@
171   DMPlexMigrateGlobalToNaturalSF - Migrates the input `sfNatural` based on sfMigration
172 
173   Input Parameters:
174 + dmOld        - The original `DM`
175 . dmNew        - The `DM` to be migrated to
176 . sfNaturalOld - The sfNatural for the `dmOld`
177 - sfMigration  - The `PetscSF` used to distribute the mesh, or `NULL` if it cannot be computed
178 
179   Output Parameter:
180 . sfNaturalNew - `PetscSF` for mapping the `Vec` in PETSc ordering to the canonical ordering
181 
182   Level: intermediate
183 
184   Notes:
185   `sfNaturalOld` maps from the old Global section (roots) to the natural Vec layout (leaves, may or may not be described by a PetscSection).
186   `DMPlexMigrateGlobalToNaturalSF` creates an SF to map from the old global section to the new global section (generated from `sfMigration`).
187   That SF is then composed with the `sfNaturalOld` to generate `sfNaturalNew`.
188   This also distributes and sets the local section for `dmNew`.
189 
190   This is not typically called by the user.
191 
192 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `PetscSF`, `PetscSection`, `DMPlexDistribute()`, `DMPlexDistributeField()`
193  @*/
194 PetscErrorCode DMPlexMigrateGlobalToNaturalSF(DM dmOld, DM dmNew, PetscSF sfNaturalOld, PetscSF sfMigration, PetscSF *sfNaturalNew)
195 {
196   MPI_Comm     comm;
197   PetscSection oldGlobalSection, newGlobalSection;
198   PetscInt    *remoteOffsets;
199   PetscBool    debug = PETSC_FALSE;
200 
201   PetscFunctionBegin;
202   PetscCall(PetscObjectGetComm((PetscObject)dmNew, &comm));
203   if (!sfMigration) {
204     /* If sfMigration is missing, sfNatural cannot be computed and is set to NULL */
205     *sfNaturalNew = NULL;
206     PetscFunctionReturn(PETSC_SUCCESS);
207   }
208   if (debug) PetscCall(PetscSFView(sfMigration, NULL));
209 
210   { // Create oldGlobalSection and newGlobalSection *with* localOffsets
211     PetscSection oldLocalSection, newLocalSection;
212     PetscSF      pointSF;
213 
214     PetscCall(DMGetLocalSection(dmOld, &oldLocalSection));
215     PetscCall(DMGetPointSF(dmOld, &pointSF));
216     PetscCall(PetscSectionCreateGlobalSection(oldLocalSection, pointSF, PETSC_TRUE, PETSC_TRUE, PETSC_TRUE, &oldGlobalSection));
217 
218     PetscCall(PetscSectionCreate(comm, &newLocalSection));
219     PetscCall(PetscSFDistributeSection(sfMigration, oldLocalSection, NULL, newLocalSection));
220     PetscCall(DMSetLocalSection(dmNew, newLocalSection));
221 
222     PetscCall(PetscSectionCreate(comm, &newGlobalSection));
223     PetscCall(DMGetPointSF(dmNew, &pointSF));
224     PetscCall(PetscSectionCreateGlobalSection(newLocalSection, pointSF, PETSC_TRUE, PETSC_TRUE, PETSC_TRUE, &newGlobalSection));
225 
226     PetscCall(PetscObjectSetName((PetscObject)oldLocalSection, "Old Local Section"));
227     if (debug) PetscCall(PetscSectionView(oldLocalSection, NULL));
228     PetscCall(PetscObjectSetName((PetscObject)oldGlobalSection, "Old Global Section"));
229     if (debug) PetscCall(PetscSectionView(oldGlobalSection, NULL));
230     PetscCall(PetscObjectSetName((PetscObject)newLocalSection, "New Local Section"));
231     if (debug) PetscCall(PetscSectionView(newLocalSection, NULL));
232     PetscCall(PetscObjectSetName((PetscObject)newGlobalSection, "New Global Section"));
233     if (debug) PetscCall(PetscSectionView(newGlobalSection, NULL));
234     PetscCall(PetscSectionDestroy(&newLocalSection));
235   }
236 
237   { // Create remoteOffsets array, mapping the oldGlobalSection offsets to the local points (according to sfMigration)
238     PetscInt lpStart, lpEnd, rpStart, rpEnd;
239 
240     PetscCall(PetscSectionGetChart(oldGlobalSection, &rpStart, &rpEnd));
241     PetscCall(PetscSectionGetChart(newGlobalSection, &lpStart, &lpEnd));
242 
243     // in `PetscSFDistributeSection` (where this is taken from), it possibly makes a new embedded SF. Should possibly do that here?
244     PetscCall(PetscMalloc1(lpEnd - lpStart, &remoteOffsets));
245     PetscCall(PetscSFBcastBegin(sfMigration, MPIU_INT, PetscSafePointerPlusOffset(oldGlobalSection->atlasOff, -rpStart), PetscSafePointerPlusOffset(remoteOffsets, -lpStart), MPI_REPLACE));
246     PetscCall(PetscSFBcastEnd(sfMigration, MPIU_INT, PetscSafePointerPlusOffset(oldGlobalSection->atlasOff, -rpStart), PetscSafePointerPlusOffset(remoteOffsets, -lpStart), MPI_REPLACE));
247     if (debug) {
248       PetscViewer viewer;
249 
250       PetscCall(PetscPrintf(comm, "Remote Offsets:\n"));
251       PetscCall(PetscViewerASCIIGetStdout(comm, &viewer));
252       PetscCall(PetscIntView(lpEnd - lpStart, remoteOffsets, viewer));
253     }
254   }
255 
256   { // Create SF from oldGlobalSection to newGlobalSection and compose with sfNaturalOld
257     PetscSF oldglob_to_newglob_sf, newglob_to_oldglob_sf;
258 
259     PetscCall(PetscSFCreateSectionSF(sfMigration, oldGlobalSection, remoteOffsets, newGlobalSection, &oldglob_to_newglob_sf));
260     PetscCall(PetscObjectSetName((PetscObject)oldglob_to_newglob_sf, "OldGlobal-to-NewGlobal SF"));
261     if (debug) PetscCall(PetscSFView(oldglob_to_newglob_sf, NULL));
262 
263     PetscCall(PetscSFCreateInverseSF(oldglob_to_newglob_sf, &newglob_to_oldglob_sf));
264     PetscCall(PetscObjectSetName((PetscObject)newglob_to_oldglob_sf, "NewGlobal-to-OldGlobal SF"));
265     PetscCall(PetscObjectViewFromOptions((PetscObject)newglob_to_oldglob_sf, (PetscObject)dmOld, "-natural_migrate_sf_view"));
266     PetscCall(PetscSFCompose(newglob_to_oldglob_sf, sfNaturalOld, sfNaturalNew));
267 
268     PetscCall(PetscSFDestroy(&oldglob_to_newglob_sf));
269     PetscCall(PetscSFDestroy(&newglob_to_oldglob_sf));
270   }
271 
272   PetscCall(PetscSectionDestroy(&oldGlobalSection));
273   PetscCall(PetscSectionDestroy(&newGlobalSection));
274   PetscCall(PetscFree(remoteOffsets));
275   PetscFunctionReturn(PETSC_SUCCESS);
276 }
277 
278 /*@
279   DMPlexGlobalToNaturalBegin - Rearranges a global `Vec` in the natural order.
280 
281   Collective
282 
283   Input Parameters:
284 + dm - The distributed `DMPLEX`
285 - gv - The global `Vec`
286 
287   Output Parameter:
288 . nv - `Vec` in the canonical ordering distributed over all processors associated with `gv`
289 
290   Level: intermediate
291 
292   Note:
293   The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`.
294 
295 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalEnd()`
296 @*/
297 PetscErrorCode DMPlexGlobalToNaturalBegin(DM dm, Vec gv, Vec nv)
298 {
299   const PetscScalar *inarray;
300   PetscScalar       *outarray;
301   MPI_Comm           comm;
302   PetscMPIInt        size;
303 
304   PetscFunctionBegin;
305   PetscCall(PetscLogEventBegin(DMPLEX_GlobalToNaturalBegin, dm, 0, 0, 0));
306   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
307   PetscCallMPI(MPI_Comm_size(comm, &size));
308   if (dm->sfNatural) {
309     if (PetscDefined(USE_DEBUG)) {
310       PetscSection gs;
311       PetscInt     Nl, n;
312 
313       PetscCall(PetscSFGetGraph(dm->sfNatural, NULL, &Nl, NULL, NULL));
314       PetscCall(VecGetLocalSize(nv, &n));
315       PetscCheck(n == Nl, comm, PETSC_ERR_ARG_INCOMP, "Natural vector local size %" PetscInt_FMT " != %" PetscInt_FMT " local size of natural section", n, Nl);
316 
317       PetscCall(DMGetGlobalSection(dm, &gs));
318       PetscCall(PetscSectionGetConstrainedStorageSize(gs, &Nl));
319       PetscCall(VecGetLocalSize(gv, &n));
320       PetscCheck(n == Nl, comm, PETSC_ERR_ARG_INCOMP, "Global vector local size %" PetscInt_FMT " != %" PetscInt_FMT " local size of global section", n, Nl);
321     }
322     PetscCall(VecGetArray(nv, &outarray));
323     PetscCall(VecGetArrayRead(gv, &inarray));
324     PetscCall(PetscSFBcastBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_REPLACE));
325     PetscCall(VecRestoreArrayRead(gv, &inarray));
326     PetscCall(VecRestoreArray(nv, &outarray));
327   } else if (size == 1) {
328     PetscCall(VecCopy(gv, nv));
329   } else {
330     PetscCheck(!dm->useNatural, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "DM global to natural SF not present. If DMPlexDistribute() was called and a section was defined, report to petsc-maint@mcs.anl.gov.");
331     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created. You must call DMSetUseNatural() before DMPlexDistribute().");
332   }
333   PetscCall(PetscLogEventEnd(DMPLEX_GlobalToNaturalBegin, dm, 0, 0, 0));
334   PetscFunctionReturn(PETSC_SUCCESS);
335 }
336 
337 /*@
338   DMPlexGlobalToNaturalEnd - Rearranges a global `Vec` in the natural order.
339 
340   Collective
341 
342   Input Parameters:
343 + dm - The distributed `DMPLEX`
344 - gv - The global `Vec`
345 
346   Output Parameter:
347 . nv - The natural `Vec`
348 
349   Level: intermediate
350 
351   Note:
352   The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`.
353 
354 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()`
355  @*/
356 PetscErrorCode DMPlexGlobalToNaturalEnd(DM dm, Vec gv, Vec nv)
357 {
358   const PetscScalar *inarray;
359   PetscScalar       *outarray;
360   PetscMPIInt        size;
361 
362   PetscFunctionBegin;
363   PetscCall(PetscLogEventBegin(DMPLEX_GlobalToNaturalEnd, dm, 0, 0, 0));
364   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
365   if (dm->sfNatural) {
366     PetscCall(VecGetArrayRead(gv, &inarray));
367     PetscCall(VecGetArray(nv, &outarray));
368     PetscCall(PetscSFBcastEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_REPLACE));
369     PetscCall(VecRestoreArrayRead(gv, &inarray));
370     PetscCall(VecRestoreArray(nv, &outarray));
371   } else if (size == 1) {
372   } else {
373     PetscCheck(!dm->useNatural, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "DM global to natural SF not present. If DMPlexDistribute() was called and a section was defined, report to petsc-maint@mcs.anl.gov.");
374     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created. You must call DMSetUseNatural() before DMPlexDistribute().");
375   }
376   PetscCall(PetscLogEventEnd(DMPLEX_GlobalToNaturalEnd, dm, 0, 0, 0));
377   PetscFunctionReturn(PETSC_SUCCESS);
378 }
379 
380 /*@
381   DMPlexNaturalToGlobalBegin - Rearranges a `Vec` in the natural order to the Global order.
382 
383   Collective
384 
385   Input Parameters:
386 + dm - The distributed `DMPLEX`
387 - nv - The natural `Vec`
388 
389   Output Parameter:
390 . gv - The global `Vec`
391 
392   Level: intermediate
393 
394   Note:
395   The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`.
396 
397 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexGlobalToNaturalEnd()`
398 @*/
399 PetscErrorCode DMPlexNaturalToGlobalBegin(DM dm, Vec nv, Vec gv)
400 {
401   const PetscScalar *inarray;
402   PetscScalar       *outarray;
403   PetscMPIInt        size;
404 
405   PetscFunctionBegin;
406   PetscCall(PetscLogEventBegin(DMPLEX_NaturalToGlobalBegin, dm, 0, 0, 0));
407   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
408   if (dm->sfNatural) {
409     /* We only have access to the SF that goes from Global to Natural.
410        Instead of inverting dm->sfNatural, we can call PetscSFReduceBegin/End with MPI_Op MPI_SUM.
411        Here the SUM really does nothing since sfNatural is one to one, as long as gV is set to zero first. */
412     PetscCall(VecZeroEntries(gv));
413     PetscCall(VecGetArray(gv, &outarray));
414     PetscCall(VecGetArrayRead(nv, &inarray));
415     PetscCall(PetscSFReduceBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_SUM));
416     PetscCall(VecRestoreArrayRead(nv, &inarray));
417     PetscCall(VecRestoreArray(gv, &outarray));
418   } else if (size == 1) {
419     PetscCall(VecCopy(nv, gv));
420   } else {
421     PetscCheck(!dm->useNatural, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "DM global to natural SF not present. If DMPlexDistribute() was called and a section was defined, report to petsc-maint@mcs.anl.gov.");
422     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created. You must call DMSetUseNatural() before DMPlexDistribute().");
423   }
424   PetscCall(PetscLogEventEnd(DMPLEX_NaturalToGlobalBegin, dm, 0, 0, 0));
425   PetscFunctionReturn(PETSC_SUCCESS);
426 }
427 
428 /*@
429   DMPlexNaturalToGlobalEnd - Rearranges a `Vec` in the natural order to the Global order.
430 
431   Collective
432 
433   Input Parameters:
434 + dm - The distributed `DMPLEX`
435 - nv - The natural `Vec`
436 
437   Output Parameter:
438 . gv - The global `Vec`
439 
440   Level: intermediate
441 
442   Note:
443   The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`.
444 
445 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()`
446  @*/
447 PetscErrorCode DMPlexNaturalToGlobalEnd(DM dm, Vec nv, Vec gv)
448 {
449   const PetscScalar *inarray;
450   PetscScalar       *outarray;
451   PetscMPIInt        size;
452 
453   PetscFunctionBegin;
454   PetscCall(PetscLogEventBegin(DMPLEX_NaturalToGlobalEnd, dm, 0, 0, 0));
455   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
456   if (dm->sfNatural) {
457     PetscCall(VecGetArrayRead(nv, &inarray));
458     PetscCall(VecGetArray(gv, &outarray));
459     PetscCall(PetscSFReduceEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_SUM));
460     PetscCall(VecRestoreArrayRead(nv, &inarray));
461     PetscCall(VecRestoreArray(gv, &outarray));
462   } else if (size == 1) {
463   } else {
464     PetscCheck(!dm->useNatural, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "DM global to natural SF not present. If DMPlexDistribute() was called and a section was defined, report to petsc-maint@mcs.anl.gov.");
465     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created. You must call DMSetUseNatural() before DMPlexDistribute().");
466   }
467   PetscCall(PetscLogEventEnd(DMPLEX_NaturalToGlobalEnd, dm, 0, 0, 0));
468   PetscFunctionReturn(PETSC_SUCCESS);
469 }
470 
471 /*@
472   DMPlexCreateNaturalVector - Provide a `Vec` capable of holding the natural ordering and distribution.
473 
474   Collective
475 
476   Input Parameter:
477 . dm - The distributed `DMPLEX`
478 
479   Output Parameter:
480 . nv - The natural `Vec`
481 
482   Level: intermediate
483 
484   Note:
485   The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`.
486 
487 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()`
488  @*/
489 PetscErrorCode DMPlexCreateNaturalVector(DM dm, Vec *nv)
490 {
491   PetscMPIInt size;
492 
493   PetscFunctionBegin;
494   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
495   if (dm->sfNatural) {
496     PetscInt nleaves, bs, maxbs;
497     Vec      v;
498 
499     /*
500       Setting the natural vector block size.
501       We can't get it from a global vector because of constraints, and the block size in the local vector
502       may be inconsistent across processes, typically when some local vectors have size 0, their block size is set to 1
503     */
504     PetscCall(DMGetLocalVector(dm, &v));
505     PetscCall(VecGetBlockSize(v, &bs));
506     PetscCallMPI(MPIU_Allreduce(&bs, &maxbs, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm)));
507     if (bs == 1 && maxbs > 1) bs = maxbs;
508     PetscCall(DMRestoreLocalVector(dm, &v));
509 
510     PetscCall(PetscSFGetGraph(dm->sfNatural, NULL, &nleaves, NULL, NULL));
511     PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), nv));
512     PetscCall(VecSetSizes(*nv, nleaves, PETSC_DETERMINE));
513     PetscCall(VecSetBlockSize(*nv, bs));
514     PetscCall(VecSetType(*nv, dm->vectype));
515     PetscCall(VecSetDM(*nv, dm));
516   } else if (size == 1) {
517     PetscCall(DMCreateLocalVector(dm, nv));
518   } else {
519     PetscCheck(!dm->useNatural, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "DM global to natural SF not present. If DMPlexDistribute() was called and a section was defined, report to petsc-maint@mcs.anl.gov.");
520     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created. You must call DMSetUseNatural() before DMPlexDistribute().");
521   }
522   PetscFunctionReturn(PETSC_SUCCESS);
523 }
524