xref: /petsc/src/dm/impls/plex/plexnatural.c (revision a1cb98fac0cdf0eb4d3e8a0c8b58f3fe8f800bc6)
1fa534816SMatthew G. Knepley #include <petsc/private/dmpleximpl.h> /*I      "petscdmplex.h"   I*/
2fa534816SMatthew G. Knepley 
3fa534816SMatthew G. Knepley /*@
4*a1cb98faSBarry Smith   DMPlexSetMigrationSF - Sets the `PetscSF` for migrating from a parent `DM` into this `DM`
5*a1cb98faSBarry Smith 
6*a1cb98faSBarry Smith   Logically Collective on dm
7f94b4a02SBlaise Bourdin 
85d3b26e6SMatthew G. Knepley   Input Parameters:
9*a1cb98faSBarry Smith + dm        - The `DM`
10*a1cb98faSBarry Smith - naturalSF - The `PetscSF`
115d3b26e6SMatthew G. Knepley 
12f94b4a02SBlaise Bourdin   Level: intermediate
13f94b4a02SBlaise Bourdin 
14*a1cb98faSBarry Smith   Note:
15*a1cb98faSBarry Smith   It is necessary to call this in order to have `DMCreateSubDM()` or `DMCreateSuperDM()` build the Global-To-Natural map
16*a1cb98faSBarry Smith 
17*a1cb98faSBarry Smith .seealso: [](chapter_unstructured), `DM`, `DMPLEX`, `PetscSF`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexCreateMigrationSF()`, `DMPlexGetMigrationSF()`
18f94b4a02SBlaise Bourdin @*/
19d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexSetMigrationSF(DM dm, PetscSF migrationSF)
20d71ae5a4SJacob Faibussowitsch {
21f94b4a02SBlaise Bourdin   PetscFunctionBegin;
22f94b4a02SBlaise Bourdin   dm->sfMigration = migrationSF;
239566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)migrationSF));
24f94b4a02SBlaise Bourdin   PetscFunctionReturn(0);
25f94b4a02SBlaise Bourdin }
26f94b4a02SBlaise Bourdin 
27f94b4a02SBlaise Bourdin /*@
28*a1cb98faSBarry Smith   DMPlexGetMigrationSF - Gets the `PetscSF` for migrating from a parent `DM` into this `DM`
29*a1cb98faSBarry Smith 
30*a1cb98faSBarry Smith   Note Collective
31f94b4a02SBlaise Bourdin 
325d3b26e6SMatthew G. Knepley   Input Parameter:
33*a1cb98faSBarry Smith . dm          - The `DM`
345d3b26e6SMatthew G. Knepley 
355d3b26e6SMatthew G. Knepley   Output Parameter:
36*a1cb98faSBarry Smith . migrationSF - The `PetscSF`
375d3b26e6SMatthew G. Knepley 
38f94b4a02SBlaise Bourdin   Level: intermediate
39f94b4a02SBlaise Bourdin 
40*a1cb98faSBarry Smith .seealso: [](chapter_unstructured), `DM`, `DMPLEX`, `PetscSF`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexCreateMigrationSF()`, `DMPlexSetMigrationSF`
41f94b4a02SBlaise Bourdin @*/
42d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGetMigrationSF(DM dm, PetscSF *migrationSF)
43d71ae5a4SJacob Faibussowitsch {
44f94b4a02SBlaise Bourdin   PetscFunctionBegin;
45f94b4a02SBlaise Bourdin   *migrationSF = dm->sfMigration;
46f94b4a02SBlaise Bourdin   PetscFunctionReturn(0);
47f94b4a02SBlaise Bourdin }
48f94b4a02SBlaise Bourdin 
49f94b4a02SBlaise Bourdin /*@
50*a1cb98faSBarry Smith   DMPlexSetGlobalToNaturalSF - Sets the `PetscSF` for mapping Global `Vec` to the Natural `Vec`
51f94b4a02SBlaise Bourdin 
525d3b26e6SMatthew G. Knepley   Input Parameters:
53*a1cb98faSBarry Smith + dm          - The `DM`
54*a1cb98faSBarry Smith - naturalSF   - The `PetscSF`
555d3b26e6SMatthew G. Knepley 
56f94b4a02SBlaise Bourdin   Level: intermediate
57f94b4a02SBlaise Bourdin 
58*a1cb98faSBarry Smith .seealso: [](chapter_unstructured), `DM`, `DMPLEX`, `PetscSF`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexCreateGlobalToNaturalSF()`, `DMPlexGetGlobaltoNaturalSF()`
59f94b4a02SBlaise Bourdin @*/
60d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexSetGlobalToNaturalSF(DM dm, PetscSF naturalSF)
61d71ae5a4SJacob Faibussowitsch {
62f94b4a02SBlaise Bourdin   PetscFunctionBegin;
63f94b4a02SBlaise Bourdin   dm->sfNatural = naturalSF;
649566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)naturalSF));
65f94b4a02SBlaise Bourdin   dm->useNatural = PETSC_TRUE;
66f94b4a02SBlaise Bourdin   PetscFunctionReturn(0);
67f94b4a02SBlaise Bourdin }
68f94b4a02SBlaise Bourdin 
69f94b4a02SBlaise Bourdin /*@
70*a1cb98faSBarry Smith   DMPlexGetGlobalToNaturalSF - Gets the `PetscSF` for mapping Global `Vec` to the Natural `Vec`
71f94b4a02SBlaise Bourdin 
725d3b26e6SMatthew G. Knepley   Input Parameter:
73*a1cb98faSBarry Smith . dm          - The `DM`
745d3b26e6SMatthew G. Knepley 
755d3b26e6SMatthew G. Knepley   Output Parameter:
76*a1cb98faSBarry Smith . naturalSF   - The `PetscSF`
775d3b26e6SMatthew G. Knepley 
78f94b4a02SBlaise Bourdin   Level: intermediate
79f94b4a02SBlaise Bourdin 
80*a1cb98faSBarry Smith .seealso: [](chapter_unstructured), `DM`, `DMPLEX`, `PetscSF`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexCreateGlobalToNaturalSF()`, `DMPlexSetGlobaltoNaturalSF`
81f94b4a02SBlaise Bourdin @*/
82d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGetGlobalToNaturalSF(DM dm, PetscSF *naturalSF)
83d71ae5a4SJacob Faibussowitsch {
84f94b4a02SBlaise Bourdin   PetscFunctionBegin;
85f94b4a02SBlaise Bourdin   *naturalSF = dm->sfNatural;
86f94b4a02SBlaise Bourdin   PetscFunctionReturn(0);
87f94b4a02SBlaise Bourdin }
88f94b4a02SBlaise Bourdin 
89f94b4a02SBlaise Bourdin /*@
90*a1cb98faSBarry Smith   DMPlexCreateGlobalToNaturalSF - Creates the `PetscSF` for mapping Global `Vec` to the Natural `Vec`
91fa534816SMatthew G. Knepley 
92fa534816SMatthew G. Knepley   Input Parameters:
93*a1cb98faSBarry Smith + dm          - The redistributed `DM`
94*a1cb98faSBarry Smith . section     - The local `PetscSection` describing the `Vec` before the mesh was distributed, or NULL if not available
95*a1cb98faSBarry Smith - sfMigration - The `PetscSF` used to distribute the mesh, or NULL if it cannot be computed
96fa534816SMatthew G. Knepley 
975d3b26e6SMatthew G. Knepley   Output Parameter:
98*a1cb98faSBarry Smith . sfNatural   - `PetscSF` for mapping the `Vec` in PETSc ordering to the canonical ordering
995d3b26e6SMatthew G. Knepley 
100fa534816SMatthew G. Knepley   Level: intermediate
101fa534816SMatthew G. Knepley 
102*a1cb98faSBarry Smith   Note:
103*a1cb98faSBarry Smith   This is not typically called by the user.
104*a1cb98faSBarry Smith 
105*a1cb98faSBarry Smith .seealso: [](chapter_unstructured), `DM`, `DMPLEX`, `PetscSF`, `PetscSection`, `DMPlexDistribute()`, `DMPlexDistributeField()`
106fa534816SMatthew G. Knepley  @*/
107d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreateGlobalToNaturalSF(DM dm, PetscSection section, PetscSF sfMigration, PetscSF *sfNatural)
108d71ae5a4SJacob Faibussowitsch {
109fa534816SMatthew G. Knepley   MPI_Comm     comm;
110e7858701SMatthew G. Knepley   PetscSF      sf, sfEmbed, sfField;
111e5b44f4fSMatthew G. Knepley   PetscSection gSection, sectionDist, gLocSection;
112fa534816SMatthew G. Knepley   PetscInt    *spoints, *remoteOffsets;
11357a2b1f8SBlaise Bourdin   PetscInt     ssize, pStart, pEnd, p, localSize, maxStorageSize;
114e7858701SMatthew G. Knepley   PetscBool    destroyFlag = PETSC_FALSE, debug = PETSC_FALSE;
115fa534816SMatthew G. Knepley 
116fa534816SMatthew G. Knepley   PetscFunctionBegin;
1179566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
1185c3f5608SAlexisMarb   if (!sfMigration) {
119e7858701SMatthew G. Knepley     /* If sfMigration is missing, sfNatural cannot be computed and is set to NULL */
1205c3f5608SAlexisMarb     *sfNatural = NULL;
1215c3f5608SAlexisMarb     PetscFunctionReturn(0);
1225c3f5608SAlexisMarb   } else if (!section) {
123e7858701SMatthew G. Knepley     /* If the sequential section is not provided (NULL), it is reconstructed from the parallel section */
1245c3f5608SAlexisMarb     PetscSF      sfMigrationInv;
1255c3f5608SAlexisMarb     PetscSection localSection;
1265c3f5608SAlexisMarb 
1279566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dm, &localSection));
1289566063dSJacob Faibussowitsch     PetscCall(PetscSFCreateInverseSF(sfMigration, &sfMigrationInv));
1299566063dSJacob Faibussowitsch     PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), &section));
1309566063dSJacob Faibussowitsch     PetscCall(PetscSFDistributeSection(sfMigrationInv, localSection, NULL, section));
1319566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&sfMigrationInv));
1325c3f5608SAlexisMarb     destroyFlag = PETSC_TRUE;
1335c3f5608SAlexisMarb   }
134e7858701SMatthew G. Knepley   if (debug) PetscCall(PetscSFView(sfMigration, NULL));
135e5b44f4fSMatthew G. Knepley   /* Create a new section from distributing the original section */
1369566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(comm, &sectionDist));
1379566063dSJacob Faibussowitsch   PetscCall(PetscSFDistributeSection(sfMigration, section, &remoteOffsets, sectionDist));
138e7858701SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)sectionDist, "Migrated Section"));
139e7858701SMatthew G. Knepley   if (debug) PetscCall(PetscSectionView(sectionDist, NULL));
1409566063dSJacob Faibussowitsch   PetscCall(DMSetLocalSection(dm, sectionDist));
141e7858701SMatthew G. Knepley   /* If a sequential section is provided but no dof is affected, sfNatural cannot be computed and is set to NULL */
14257a2b1f8SBlaise Bourdin   PetscCall(PetscSectionGetStorageSize(sectionDist, &localSize));
14357a2b1f8SBlaise Bourdin   PetscCallMPI(MPI_Allreduce(&localSize, &maxStorageSize, 1, MPI_INT, MPI_MAX, PetscObjectComm((PetscObject)dm)));
14457a2b1f8SBlaise Bourdin   if (maxStorageSize) {
145e7858701SMatthew G. Knepley     const PetscInt *leaves;
146e7858701SMatthew G. Knepley     PetscInt       *sortleaves, *indices;
147e7858701SMatthew G. Knepley     PetscInt        Nl;
148e7858701SMatthew G. Knepley 
149fa534816SMatthew G. Knepley     /* Get a pruned version of migration SF */
1509566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalSection(dm, &gSection));
151e7858701SMatthew G. Knepley     if (debug) PetscCall(PetscSectionView(gSection, NULL));
152e7858701SMatthew G. Knepley     PetscCall(PetscSFGetGraph(sfMigration, NULL, &Nl, &leaves, NULL));
1539566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetChart(gSection, &pStart, &pEnd));
154fa534816SMatthew G. Knepley     for (p = pStart, ssize = 0; p < pEnd; ++p) {
155fa534816SMatthew G. Knepley       PetscInt dof, off;
156fa534816SMatthew G. Knepley 
1579566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(gSection, p, &dof));
1589566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(gSection, p, &off));
159fa534816SMatthew G. Knepley       if ((dof > 0) && (off >= 0)) ++ssize;
160fa534816SMatthew G. Knepley     }
161e7858701SMatthew G. Knepley     PetscCall(PetscMalloc3(ssize, &spoints, Nl, &sortleaves, Nl, &indices));
1629371c9d4SSatish Balay     for (p = 0; p < Nl; ++p) {
1639371c9d4SSatish Balay       sortleaves[p] = leaves ? leaves[p] : p;
1649371c9d4SSatish Balay       indices[p]    = p;
1659371c9d4SSatish Balay     }
166e7858701SMatthew G. Knepley     PetscCall(PetscSortIntWithArray(Nl, sortleaves, indices));
167fa534816SMatthew G. Knepley     for (p = pStart, ssize = 0; p < pEnd; ++p) {
168e7858701SMatthew G. Knepley       PetscInt dof, off, loc;
169fa534816SMatthew G. Knepley 
1709566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(gSection, p, &dof));
1719566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(gSection, p, &off));
172e7858701SMatthew G. Knepley       if ((dof > 0) && (off >= 0)) {
173e7858701SMatthew G. Knepley         PetscCall(PetscFindInt(p, Nl, sortleaves, &loc));
174e7858701SMatthew G. Knepley         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);
175e7858701SMatthew G. Knepley         spoints[ssize++] = indices[loc];
176e7858701SMatthew G. Knepley       }
177fa534816SMatthew G. Knepley     }
1789566063dSJacob Faibussowitsch     PetscCall(PetscSFCreateEmbeddedLeafSF(sfMigration, ssize, spoints, &sfEmbed));
179e7858701SMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)sfEmbed, "Embedded SF"));
180e7858701SMatthew G. Knepley     PetscCall(PetscFree3(spoints, sortleaves, indices));
181e7858701SMatthew G. Knepley     if (debug) PetscCall(PetscSFView(sfEmbed, NULL));
182e7858701SMatthew G. Knepley     /* Create the SF associated with this section
183e7858701SMatthew G. Knepley          Roots are natural dofs, leaves are global dofs */
1849566063dSJacob Faibussowitsch     PetscCall(DMGetPointSF(dm, &sf));
18599a026b9SAlexis Marboeuf     PetscCall(PetscSectionCreateGlobalSection(sectionDist, sf, PETSC_TRUE, PETSC_TRUE, &gLocSection));
1869566063dSJacob Faibussowitsch     PetscCall(PetscSFCreateSectionSF(sfEmbed, section, remoteOffsets, gLocSection, &sfField));
1879566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&sfEmbed));
1889566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&gLocSection));
189e7858701SMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)sfField, "Natural-to-Global SF"));
190e7858701SMatthew G. Knepley     if (debug) PetscCall(PetscSFView(sfField, NULL));
191e7858701SMatthew G. Knepley     /* Invert the field SF
192e7858701SMatthew G. Knepley          Roots are global dofs, leaves are natural dofs */
193e7858701SMatthew G. Knepley     PetscCall(PetscSFCreateInverseSF(sfField, sfNatural));
194e7858701SMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)*sfNatural, "Global-to-Natural SF"));
1959566063dSJacob Faibussowitsch     PetscCall(PetscObjectViewFromOptions((PetscObject)*sfNatural, NULL, "-globaltonatural_sf_view"));
196fa534816SMatthew G. Knepley     /* Clean up */
197e7858701SMatthew G. Knepley     PetscCall(PetscSFDestroy(&sfField));
19814744f87SBlaise Bourdin   } else {
19914744f87SBlaise Bourdin     *sfNatural = NULL;
20014744f87SBlaise Bourdin   }
2019566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&sectionDist));
2029566063dSJacob Faibussowitsch   PetscCall(PetscFree(remoteOffsets));
2039566063dSJacob Faibussowitsch   if (destroyFlag) PetscCall(PetscSectionDestroy(&section));
204fa534816SMatthew G. Knepley   PetscFunctionReturn(0);
205fa534816SMatthew G. Knepley }
206fa534816SMatthew G. Knepley 
207fa534816SMatthew G. Knepley /*@
208*a1cb98faSBarry Smith   DMPlexGlobalToNaturalBegin - Rearranges a global `Vec` in the natural order.
209fa534816SMatthew G. Knepley 
210fa534816SMatthew G. Knepley   Collective on dm
211fa534816SMatthew G. Knepley 
212fa534816SMatthew G. Knepley   Input Parameters:
213*a1cb98faSBarry Smith + dm - The distributed `DMPLEX`
214*a1cb98faSBarry Smith - gv - The global `Vec`
215fa534816SMatthew G. Knepley 
216fa534816SMatthew G. Knepley   Output Parameters:
217*a1cb98faSBarry Smith . nv - `Vec` in the canonical ordering distributed over all processors associated with gv
218fa534816SMatthew G. Knepley 
219fa534816SMatthew G. Knepley   Level: intermediate
220fa534816SMatthew G. Knepley 
221*a1cb98faSBarry Smith   Note:
222*a1cb98faSBarry Smith   The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`.
223*a1cb98faSBarry Smith 
224*a1cb98faSBarry Smith .seealso: [](chapter_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalEnd()`
225fa534816SMatthew G. Knepley @*/
226d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGlobalToNaturalBegin(DM dm, Vec gv, Vec nv)
227d71ae5a4SJacob Faibussowitsch {
228fa534816SMatthew G. Knepley   const PetscScalar *inarray;
229fa534816SMatthew G. Knepley   PetscScalar       *outarray;
230a6a55facSMatthew G. Knepley   PetscMPIInt        size;
231fa534816SMatthew G. Knepley 
232fa534816SMatthew G. Knepley   PetscFunctionBegin;
2339566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_GlobalToNaturalBegin, dm, 0, 0, 0));
2349566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
235fa534816SMatthew G. Knepley   if (dm->sfNatural) {
2369566063dSJacob Faibussowitsch     PetscCall(VecGetArray(nv, &outarray));
2379566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(gv, &inarray));
2389566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_REPLACE));
2399566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(gv, &inarray));
2409566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(nv, &outarray));
241a6a55facSMatthew G. Knepley   } else if (size == 1) {
2429566063dSJacob Faibussowitsch     PetscCall(VecCopy(gv, nv));
243f7d195e4SLawrence Mitchell   } else {
244f7d195e4SLawrence Mitchell     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.");
245f7d195e4SLawrence Mitchell     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().");
246f7d195e4SLawrence Mitchell   }
2479566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_GlobalToNaturalBegin, dm, 0, 0, 0));
248fa534816SMatthew G. Knepley   PetscFunctionReturn(0);
249fa534816SMatthew G. Knepley }
250fa534816SMatthew G. Knepley 
251fa534816SMatthew G. Knepley /*@
252*a1cb98faSBarry Smith   DMPlexGlobalToNaturalEnd - Rearranges a global `Vec` in the natural order.
253fa534816SMatthew G. Knepley 
254fa534816SMatthew G. Knepley   Collective on dm
255fa534816SMatthew G. Knepley 
256fa534816SMatthew G. Knepley   Input Parameters:
257*a1cb98faSBarry Smith + dm - The distributed `DMPLEX`
258*a1cb98faSBarry Smith - gv - The global `Vec`
259fa534816SMatthew G. Knepley 
26097bb3fdcSJose E. Roman   Output Parameter:
261*a1cb98faSBarry Smith . nv - The natural `Vec`
262fa534816SMatthew G. Knepley 
263fa534816SMatthew G. Knepley   Level: intermediate
264fa534816SMatthew G. Knepley 
265*a1cb98faSBarry Smith   Note:
266*a1cb98faSBarry Smith   The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`.
267*a1cb98faSBarry Smith 
268*a1cb98faSBarry Smith  .seealso: [](chapter_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()`
269fa534816SMatthew G. Knepley  @*/
270d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGlobalToNaturalEnd(DM dm, Vec gv, Vec nv)
271d71ae5a4SJacob Faibussowitsch {
272fa534816SMatthew G. Knepley   const PetscScalar *inarray;
273fa534816SMatthew G. Knepley   PetscScalar       *outarray;
274a6a55facSMatthew G. Knepley   PetscMPIInt        size;
275fa534816SMatthew G. Knepley 
276fa534816SMatthew G. Knepley   PetscFunctionBegin;
2779566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_GlobalToNaturalEnd, dm, 0, 0, 0));
2789566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
279fa534816SMatthew G. Knepley   if (dm->sfNatural) {
2809566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(gv, &inarray));
2819566063dSJacob Faibussowitsch     PetscCall(VecGetArray(nv, &outarray));
2829566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_REPLACE));
2839566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(gv, &inarray));
2849566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(nv, &outarray));
285a6a55facSMatthew G. Knepley   } else if (size == 1) {
286f7d195e4SLawrence Mitchell   } else {
287f7d195e4SLawrence Mitchell     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.");
288f7d195e4SLawrence Mitchell     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().");
289f7d195e4SLawrence Mitchell   }
2909566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_GlobalToNaturalEnd, dm, 0, 0, 0));
291fa534816SMatthew G. Knepley   PetscFunctionReturn(0);
292fa534816SMatthew G. Knepley }
293fa534816SMatthew G. Knepley 
294fa534816SMatthew G. Knepley /*@
295*a1cb98faSBarry Smith   DMPlexNaturalToGlobalBegin - Rearranges a `Vec` in the natural order to the Global order.
296fa534816SMatthew G. Knepley 
297fa534816SMatthew G. Knepley   Collective on dm
298fa534816SMatthew G. Knepley 
299fa534816SMatthew G. Knepley   Input Parameters:
300*a1cb98faSBarry Smith + dm - The distributed `DMPLEX`
301*a1cb98faSBarry Smith - nv - The natural `Vec`
302fa534816SMatthew G. Knepley 
303fa534816SMatthew G. Knepley   Output Parameters:
304*a1cb98faSBarry Smith . gv - The global `Vec`
305fa534816SMatthew G. Knepley 
306fa534816SMatthew G. Knepley   Level: intermediate
307fa534816SMatthew G. Knepley 
308*a1cb98faSBarry Smith   Note:
309*a1cb98faSBarry Smith   The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`.
310*a1cb98faSBarry Smith 
311*a1cb98faSBarry Smith .seealso: [](chapter_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalEnd()`
312fa534816SMatthew G. Knepley @*/
313d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexNaturalToGlobalBegin(DM dm, Vec nv, Vec gv)
314d71ae5a4SJacob Faibussowitsch {
315fa534816SMatthew G. Knepley   const PetscScalar *inarray;
316fa534816SMatthew G. Knepley   PetscScalar       *outarray;
317a6a55facSMatthew G. Knepley   PetscMPIInt        size;
318fa534816SMatthew G. Knepley 
319fa534816SMatthew G. Knepley   PetscFunctionBegin;
3209566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_NaturalToGlobalBegin, dm, 0, 0, 0));
3219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
322fa534816SMatthew G. Knepley   if (dm->sfNatural) {
323a5b23f4aSJose E. Roman     /* We only have access to the SF that goes from Global to Natural.
324b64f75a9SBlaise Bourdin        Instead of inverting dm->sfNatural, we can call PetscSFReduceBegin/End with MPI_Op MPI_SUM.
325b64f75a9SBlaise Bourdin        Here the SUM really does nothing since sfNatural is one to one, as long as gV is set to zero first. */
3269566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(gv));
3279566063dSJacob Faibussowitsch     PetscCall(VecGetArray(gv, &outarray));
3289566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(nv, &inarray));
3299566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_SUM));
3309566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(nv, &inarray));
3319566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(gv, &outarray));
332a6a55facSMatthew G. Knepley   } else if (size == 1) {
3339566063dSJacob Faibussowitsch     PetscCall(VecCopy(nv, gv));
334f7d195e4SLawrence Mitchell   } else {
335f7d195e4SLawrence Mitchell     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.");
336f7d195e4SLawrence Mitchell     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().");
337f7d195e4SLawrence Mitchell   }
3389566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_NaturalToGlobalBegin, dm, 0, 0, 0));
339fa534816SMatthew G. Knepley   PetscFunctionReturn(0);
340fa534816SMatthew G. Knepley }
341fa534816SMatthew G. Knepley 
342fa534816SMatthew G. Knepley /*@
343*a1cb98faSBarry Smith   DMPlexNaturalToGlobalEnd - Rearranges a `Vec` in the natural order to the Global order.
344fa534816SMatthew G. Knepley 
345fa534816SMatthew G. Knepley   Collective on dm
346fa534816SMatthew G. Knepley 
347fa534816SMatthew G. Knepley   Input Parameters:
348*a1cb98faSBarry Smith + dm - The distributed `DMPLEX`
349*a1cb98faSBarry Smith - nv - The natural `Vec`
350fa534816SMatthew G. Knepley 
351fa534816SMatthew G. Knepley   Output Parameters:
352*a1cb98faSBarry Smith . gv - The global `Vec`
353fa534816SMatthew G. Knepley 
354fa534816SMatthew G. Knepley   Level: intermediate
355fa534816SMatthew G. Knepley 
356*a1cb98faSBarry Smith   Note:
357*a1cb98faSBarry Smith   The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`.
358*a1cb98faSBarry Smith 
359*a1cb98faSBarry Smith .seealso: [](chapter_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()`
360fa534816SMatthew G. Knepley  @*/
361d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexNaturalToGlobalEnd(DM dm, Vec nv, Vec gv)
362d71ae5a4SJacob Faibussowitsch {
363fa534816SMatthew G. Knepley   const PetscScalar *inarray;
364fa534816SMatthew G. Knepley   PetscScalar       *outarray;
365a6a55facSMatthew G. Knepley   PetscMPIInt        size;
366fa534816SMatthew G. Knepley 
367fa534816SMatthew G. Knepley   PetscFunctionBegin;
3689566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_NaturalToGlobalEnd, dm, 0, 0, 0));
3699566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
370fa534816SMatthew G. Knepley   if (dm->sfNatural) {
3719566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(nv, &inarray));
3729566063dSJacob Faibussowitsch     PetscCall(VecGetArray(gv, &outarray));
3739566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_SUM));
3749566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(nv, &inarray));
3759566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(gv, &outarray));
376a6a55facSMatthew G. Knepley   } else if (size == 1) {
377f7d195e4SLawrence Mitchell   } else {
378f7d195e4SLawrence Mitchell     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.");
379f7d195e4SLawrence Mitchell     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().");
380f7d195e4SLawrence Mitchell   }
3819566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_NaturalToGlobalEnd, dm, 0, 0, 0));
382fa534816SMatthew G. Knepley   PetscFunctionReturn(0);
383fa534816SMatthew G. Knepley }
38409cf685aSAlexis Marboeuf 
38509cf685aSAlexis Marboeuf /*@
3865cb80ecdSBarry Smith   DMPlexCreateNaturalVector - Provide a `Vec` capable of holding the natural ordering and distribution.
38709cf685aSAlexis Marboeuf 
38809cf685aSAlexis Marboeuf   Collective on dm
38909cf685aSAlexis Marboeuf 
3905cb80ecdSBarry Smith   Input Parameter:
391*a1cb98faSBarry Smith . dm - The distributed `DMPLEX`
39209cf685aSAlexis Marboeuf 
3935cb80ecdSBarry Smith   Output Parameter:
3945cb80ecdSBarry Smith . nv - The natural `Vec`
39509cf685aSAlexis Marboeuf 
39609cf685aSAlexis Marboeuf   Level: intermediate
39709cf685aSAlexis Marboeuf 
398*a1cb98faSBarry Smith   Note:
399*a1cb98faSBarry Smith   The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`.
400*a1cb98faSBarry Smith 
401*a1cb98faSBarry Smith .seealso: [](chapter_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()`
40209cf685aSAlexis Marboeuf  @*/
403d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreateNaturalVector(DM dm, Vec *nv)
404d71ae5a4SJacob Faibussowitsch {
40509cf685aSAlexis Marboeuf   PetscMPIInt size;
40609cf685aSAlexis Marboeuf 
40709cf685aSAlexis Marboeuf   PetscFunctionBegin;
40809cf685aSAlexis Marboeuf   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
40909cf685aSAlexis Marboeuf   if (dm->sfNatural) {
4102e8d78feSBlaise Bourdin     PetscInt nleaves, bs;
4112e8d78feSBlaise Bourdin     Vec      v;
4122e8d78feSBlaise Bourdin     PetscCall(DMGetLocalVector(dm, &v));
4132e8d78feSBlaise Bourdin     PetscCall(VecGetBlockSize(v, &bs));
4142e8d78feSBlaise Bourdin     PetscCall(DMRestoreLocalVector(dm, &v));
4152e8d78feSBlaise Bourdin 
41609cf685aSAlexis Marboeuf     PetscCall(PetscSFGetGraph(dm->sfNatural, NULL, &nleaves, NULL, NULL));
41709cf685aSAlexis Marboeuf     PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), nv));
41809cf685aSAlexis Marboeuf     PetscCall(VecSetSizes(*nv, nleaves, PETSC_DETERMINE));
41909cf685aSAlexis Marboeuf     PetscCall(VecSetBlockSize(*nv, bs));
42009cf685aSAlexis Marboeuf     PetscCall(VecSetType(*nv, dm->vectype));
42109cf685aSAlexis Marboeuf     PetscCall(VecSetDM(*nv, dm));
42209cf685aSAlexis Marboeuf   } else if (size == 1) {
4232e8d78feSBlaise Bourdin     PetscCall(DMCreateLocalVector(dm, nv));
42409cf685aSAlexis Marboeuf   } else {
42509cf685aSAlexis Marboeuf     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.");
42609cf685aSAlexis Marboeuf     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().");
42709cf685aSAlexis Marboeuf   }
42809cf685aSAlexis Marboeuf   PetscFunctionReturn(0);
42909cf685aSAlexis Marboeuf }
430