xref: /petsc/src/dm/impls/plex/plexnatural.c (revision 1cc06b555e92f8ec64db10330b8bbd830e5bc876)
1fa534816SMatthew G. Knepley #include <petsc/private/dmpleximpl.h> /*I      "petscdmplex.h"   I*/
2fa534816SMatthew G. Knepley 
3fa534816SMatthew G. Knepley /*@
4a1cb98faSBarry Smith   DMPlexSetMigrationSF - Sets the `PetscSF` for migrating from a parent `DM` into this `DM`
5a1cb98faSBarry Smith 
620f4b53cSBarry Smith   Logically Collective
7f94b4a02SBlaise Bourdin 
85d3b26e6SMatthew G. Knepley   Input Parameters:
9a1cb98faSBarry Smith + dm        - The `DM`
10a1cb98faSBarry Smith - naturalSF - The `PetscSF`
115d3b26e6SMatthew G. Knepley 
12f94b4a02SBlaise Bourdin   Level: intermediate
13f94b4a02SBlaise Bourdin 
14a1cb98faSBarry Smith   Note:
15a1cb98faSBarry Smith   It is necessary to call this in order to have `DMCreateSubDM()` or `DMCreateSuperDM()` build the Global-To-Natural map
16a1cb98faSBarry Smith 
17*1cc06b55SBarry Smith .seealso: [](ch_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));
243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
25f94b4a02SBlaise Bourdin }
26f94b4a02SBlaise Bourdin 
27f94b4a02SBlaise Bourdin /*@
28a1cb98faSBarry Smith   DMPlexGetMigrationSF - Gets the `PetscSF` for migrating from a parent `DM` into this `DM`
29a1cb98faSBarry Smith 
30a1cb98faSBarry Smith   Note Collective
31f94b4a02SBlaise Bourdin 
325d3b26e6SMatthew G. Knepley   Input Parameter:
33a1cb98faSBarry Smith . dm          - The `DM`
345d3b26e6SMatthew G. Knepley 
355d3b26e6SMatthew G. Knepley   Output Parameter:
36a1cb98faSBarry Smith . migrationSF - The `PetscSF`
375d3b26e6SMatthew G. Knepley 
38f94b4a02SBlaise Bourdin   Level: intermediate
39f94b4a02SBlaise Bourdin 
40*1cc06b55SBarry Smith .seealso: [](ch_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;
463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
47f94b4a02SBlaise Bourdin }
48f94b4a02SBlaise Bourdin 
49f94b4a02SBlaise Bourdin /*@
50a1cb98faSBarry Smith   DMPlexSetGlobalToNaturalSF - Sets the `PetscSF` for mapping Global `Vec` to the Natural `Vec`
51f94b4a02SBlaise Bourdin 
525d3b26e6SMatthew G. Knepley   Input Parameters:
53a1cb98faSBarry Smith + dm          - The `DM`
54a1cb98faSBarry Smith - naturalSF   - The `PetscSF`
555d3b26e6SMatthew G. Knepley 
56f94b4a02SBlaise Bourdin   Level: intermediate
57f94b4a02SBlaise Bourdin 
58*1cc06b55SBarry Smith .seealso: [](ch_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;
663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
67f94b4a02SBlaise Bourdin }
68f94b4a02SBlaise Bourdin 
69f94b4a02SBlaise Bourdin /*@
70a1cb98faSBarry Smith   DMPlexGetGlobalToNaturalSF - Gets the `PetscSF` for mapping Global `Vec` to the Natural `Vec`
71f94b4a02SBlaise Bourdin 
725d3b26e6SMatthew G. Knepley   Input Parameter:
73a1cb98faSBarry Smith . dm          - The `DM`
745d3b26e6SMatthew G. Knepley 
755d3b26e6SMatthew G. Knepley   Output Parameter:
76a1cb98faSBarry Smith . naturalSF   - The `PetscSF`
775d3b26e6SMatthew G. Knepley 
78f94b4a02SBlaise Bourdin   Level: intermediate
79f94b4a02SBlaise Bourdin 
80*1cc06b55SBarry Smith .seealso: [](ch_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;
863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
87f94b4a02SBlaise Bourdin }
88f94b4a02SBlaise Bourdin 
89f94b4a02SBlaise Bourdin /*@
90a1cb98faSBarry Smith   DMPlexCreateGlobalToNaturalSF - Creates the `PetscSF` for mapping Global `Vec` to the Natural `Vec`
91fa534816SMatthew G. Knepley 
92fa534816SMatthew G. Knepley   Input Parameters:
93a1cb98faSBarry Smith + dm          - The redistributed `DM`
9420f4b53cSBarry Smith . section     - The local `PetscSection` describing the `Vec` before the mesh was distributed, or `NULL` if not available
9520f4b53cSBarry Smith - sfMigration - The `PetscSF` used to distribute the mesh, or `NULL` if it cannot be computed
96fa534816SMatthew G. Knepley 
975d3b26e6SMatthew G. Knepley   Output Parameter:
98a1cb98faSBarry 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 
102a1cb98faSBarry Smith   Note:
103a1cb98faSBarry Smith   This is not typically called by the user.
104a1cb98faSBarry Smith 
105*1cc06b55SBarry Smith .seealso: [](ch_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;
1213ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
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));
143ab9fe21fSBlaise Bourdin   PetscCall(MPIU_Allreduce(&localSize, &maxStorageSize, 1, MPIU_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));
2043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
205fa534816SMatthew G. Knepley }
206fa534816SMatthew G. Knepley 
207fa534816SMatthew G. Knepley /*@
208a1cb98faSBarry Smith   DMPlexGlobalToNaturalBegin - Rearranges a global `Vec` in the natural order.
209fa534816SMatthew G. Knepley 
21020f4b53cSBarry Smith   Collective
211fa534816SMatthew G. Knepley 
212fa534816SMatthew G. Knepley   Input Parameters:
213a1cb98faSBarry Smith + dm - The distributed `DMPLEX`
214a1cb98faSBarry Smith - gv - The global `Vec`
215fa534816SMatthew G. Knepley 
2162fe279fdSBarry Smith   Output Parameter:
21720f4b53cSBarry 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 
221a1cb98faSBarry Smith   Note:
222a1cb98faSBarry Smith   The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`.
223a1cb98faSBarry Smith 
224*1cc06b55SBarry Smith .seealso: [](ch_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;
230f6c7c579SMatthew G. Knepley   MPI_Comm           comm;
231a6a55facSMatthew G. Knepley   PetscMPIInt        size;
232fa534816SMatthew G. Knepley 
233fa534816SMatthew G. Knepley   PetscFunctionBegin;
2349566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_GlobalToNaturalBegin, dm, 0, 0, 0));
235f6c7c579SMatthew G. Knepley   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
236f6c7c579SMatthew G. Knepley   PetscCallMPI(MPI_Comm_size(comm, &size));
237fa534816SMatthew G. Knepley   if (dm->sfNatural) {
238f6c7c579SMatthew G. Knepley     if (PetscDefined(USE_DEBUG)) {
239f6c7c579SMatthew G. Knepley       PetscSection gs;
240f6c7c579SMatthew G. Knepley       PetscInt     Nl, n;
241f6c7c579SMatthew G. Knepley 
242f6c7c579SMatthew G. Knepley       PetscCall(PetscSFGetGraph(dm->sfNatural, NULL, &Nl, NULL, NULL));
243f6c7c579SMatthew G. Knepley       PetscCall(VecGetLocalSize(nv, &n));
244f6c7c579SMatthew G. Knepley       PetscCheck(n == Nl, comm, PETSC_ERR_ARG_INCOMP, "Natural vector local size %" PetscInt_FMT " != %" PetscInt_FMT " local size of natural section", n, Nl);
245f6c7c579SMatthew G. Knepley 
246f6c7c579SMatthew G. Knepley       PetscCall(DMGetGlobalSection(dm, &gs));
247f6c7c579SMatthew G. Knepley       PetscCall(PetscSectionGetConstrainedStorageSize(gs, &Nl));
248f6c7c579SMatthew G. Knepley       PetscCall(VecGetLocalSize(gv, &n));
249f6c7c579SMatthew G. Knepley       PetscCheck(n == Nl, comm, PETSC_ERR_ARG_INCOMP, "Global vector local size %" PetscInt_FMT " != %" PetscInt_FMT " local size of global section", n, Nl);
250f6c7c579SMatthew G. Knepley     }
2519566063dSJacob Faibussowitsch     PetscCall(VecGetArray(nv, &outarray));
2529566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(gv, &inarray));
2539566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_REPLACE));
2549566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(gv, &inarray));
2559566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(nv, &outarray));
256a6a55facSMatthew G. Knepley   } else if (size == 1) {
2579566063dSJacob Faibussowitsch     PetscCall(VecCopy(gv, nv));
258f7d195e4SLawrence Mitchell   } else {
259f7d195e4SLawrence 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.");
260f7d195e4SLawrence Mitchell     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().");
261f7d195e4SLawrence Mitchell   }
2629566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_GlobalToNaturalBegin, dm, 0, 0, 0));
2633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
264fa534816SMatthew G. Knepley }
265fa534816SMatthew G. Knepley 
266fa534816SMatthew G. Knepley /*@
267a1cb98faSBarry Smith   DMPlexGlobalToNaturalEnd - Rearranges a global `Vec` in the natural order.
268fa534816SMatthew G. Knepley 
26920f4b53cSBarry Smith   Collective
270fa534816SMatthew G. Knepley 
271fa534816SMatthew G. Knepley   Input Parameters:
272a1cb98faSBarry Smith + dm - The distributed `DMPLEX`
273a1cb98faSBarry Smith - gv - The global `Vec`
274fa534816SMatthew G. Knepley 
27597bb3fdcSJose E. Roman   Output Parameter:
276a1cb98faSBarry Smith . nv - The natural `Vec`
277fa534816SMatthew G. Knepley 
278fa534816SMatthew G. Knepley   Level: intermediate
279fa534816SMatthew G. Knepley 
280a1cb98faSBarry Smith   Note:
281a1cb98faSBarry Smith   The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`.
282a1cb98faSBarry Smith 
283*1cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()`
284fa534816SMatthew G. Knepley  @*/
285d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGlobalToNaturalEnd(DM dm, Vec gv, Vec nv)
286d71ae5a4SJacob Faibussowitsch {
287fa534816SMatthew G. Knepley   const PetscScalar *inarray;
288fa534816SMatthew G. Knepley   PetscScalar       *outarray;
289a6a55facSMatthew G. Knepley   PetscMPIInt        size;
290fa534816SMatthew G. Knepley 
291fa534816SMatthew G. Knepley   PetscFunctionBegin;
2929566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_GlobalToNaturalEnd, dm, 0, 0, 0));
2939566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
294fa534816SMatthew G. Knepley   if (dm->sfNatural) {
2959566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(gv, &inarray));
2969566063dSJacob Faibussowitsch     PetscCall(VecGetArray(nv, &outarray));
2979566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_REPLACE));
2989566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(gv, &inarray));
2999566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(nv, &outarray));
300a6a55facSMatthew G. Knepley   } else if (size == 1) {
301f7d195e4SLawrence Mitchell   } else {
302f7d195e4SLawrence 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.");
303f7d195e4SLawrence Mitchell     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().");
304f7d195e4SLawrence Mitchell   }
3059566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_GlobalToNaturalEnd, dm, 0, 0, 0));
3063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
307fa534816SMatthew G. Knepley }
308fa534816SMatthew G. Knepley 
309fa534816SMatthew G. Knepley /*@
310a1cb98faSBarry Smith   DMPlexNaturalToGlobalBegin - Rearranges a `Vec` in the natural order to the Global order.
311fa534816SMatthew G. Knepley 
31220f4b53cSBarry Smith   Collective
313fa534816SMatthew G. Knepley 
314fa534816SMatthew G. Knepley   Input Parameters:
315a1cb98faSBarry Smith + dm - The distributed `DMPLEX`
316a1cb98faSBarry Smith - nv - The natural `Vec`
317fa534816SMatthew G. Knepley 
3182fe279fdSBarry Smith   Output Parameter:
319a1cb98faSBarry Smith . gv - The global `Vec`
320fa534816SMatthew G. Knepley 
321fa534816SMatthew G. Knepley   Level: intermediate
322fa534816SMatthew G. Knepley 
323a1cb98faSBarry Smith   Note:
324a1cb98faSBarry Smith   The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`.
325a1cb98faSBarry Smith 
326*1cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalEnd()`
327fa534816SMatthew G. Knepley @*/
328d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexNaturalToGlobalBegin(DM dm, Vec nv, Vec gv)
329d71ae5a4SJacob Faibussowitsch {
330fa534816SMatthew G. Knepley   const PetscScalar *inarray;
331fa534816SMatthew G. Knepley   PetscScalar       *outarray;
332a6a55facSMatthew G. Knepley   PetscMPIInt        size;
333fa534816SMatthew G. Knepley 
334fa534816SMatthew G. Knepley   PetscFunctionBegin;
3359566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_NaturalToGlobalBegin, dm, 0, 0, 0));
3369566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
337fa534816SMatthew G. Knepley   if (dm->sfNatural) {
338a5b23f4aSJose E. Roman     /* We only have access to the SF that goes from Global to Natural.
339b64f75a9SBlaise Bourdin        Instead of inverting dm->sfNatural, we can call PetscSFReduceBegin/End with MPI_Op MPI_SUM.
340b64f75a9SBlaise Bourdin        Here the SUM really does nothing since sfNatural is one to one, as long as gV is set to zero first. */
3419566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(gv));
3429566063dSJacob Faibussowitsch     PetscCall(VecGetArray(gv, &outarray));
3439566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(nv, &inarray));
3449566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_SUM));
3459566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(nv, &inarray));
3469566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(gv, &outarray));
347a6a55facSMatthew G. Knepley   } else if (size == 1) {
3489566063dSJacob Faibussowitsch     PetscCall(VecCopy(nv, gv));
349f7d195e4SLawrence Mitchell   } else {
350f7d195e4SLawrence 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.");
351f7d195e4SLawrence Mitchell     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().");
352f7d195e4SLawrence Mitchell   }
3539566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_NaturalToGlobalBegin, dm, 0, 0, 0));
3543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
355fa534816SMatthew G. Knepley }
356fa534816SMatthew G. Knepley 
357fa534816SMatthew G. Knepley /*@
358a1cb98faSBarry Smith   DMPlexNaturalToGlobalEnd - Rearranges a `Vec` in the natural order to the Global order.
359fa534816SMatthew G. Knepley 
36020f4b53cSBarry Smith   Collective
361fa534816SMatthew G. Knepley 
362fa534816SMatthew G. Knepley   Input Parameters:
363a1cb98faSBarry Smith + dm - The distributed `DMPLEX`
364a1cb98faSBarry Smith - nv - The natural `Vec`
365fa534816SMatthew G. Knepley 
3662fe279fdSBarry Smith   Output Parameter:
367a1cb98faSBarry Smith . gv - The global `Vec`
368fa534816SMatthew G. Knepley 
369fa534816SMatthew G. Knepley   Level: intermediate
370fa534816SMatthew G. Knepley 
371a1cb98faSBarry Smith   Note:
372a1cb98faSBarry Smith   The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`.
373a1cb98faSBarry Smith 
374*1cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()`
375fa534816SMatthew G. Knepley  @*/
376d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexNaturalToGlobalEnd(DM dm, Vec nv, Vec gv)
377d71ae5a4SJacob Faibussowitsch {
378fa534816SMatthew G. Knepley   const PetscScalar *inarray;
379fa534816SMatthew G. Knepley   PetscScalar       *outarray;
380a6a55facSMatthew G. Knepley   PetscMPIInt        size;
381fa534816SMatthew G. Knepley 
382fa534816SMatthew G. Knepley   PetscFunctionBegin;
3839566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_NaturalToGlobalEnd, dm, 0, 0, 0));
3849566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
385fa534816SMatthew G. Knepley   if (dm->sfNatural) {
3869566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(nv, &inarray));
3879566063dSJacob Faibussowitsch     PetscCall(VecGetArray(gv, &outarray));
3889566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_SUM));
3899566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(nv, &inarray));
3909566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(gv, &outarray));
391a6a55facSMatthew G. Knepley   } else if (size == 1) {
392f7d195e4SLawrence Mitchell   } else {
393f7d195e4SLawrence 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.");
394f7d195e4SLawrence Mitchell     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().");
395f7d195e4SLawrence Mitchell   }
3969566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_NaturalToGlobalEnd, dm, 0, 0, 0));
3973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
398fa534816SMatthew G. Knepley }
39909cf685aSAlexis Marboeuf 
40009cf685aSAlexis Marboeuf /*@
4015cb80ecdSBarry Smith   DMPlexCreateNaturalVector - Provide a `Vec` capable of holding the natural ordering and distribution.
40209cf685aSAlexis Marboeuf 
40320f4b53cSBarry Smith   Collective
40409cf685aSAlexis Marboeuf 
4055cb80ecdSBarry Smith   Input Parameter:
406a1cb98faSBarry Smith . dm - The distributed `DMPLEX`
40709cf685aSAlexis Marboeuf 
4085cb80ecdSBarry Smith   Output Parameter:
4095cb80ecdSBarry Smith . nv - The natural `Vec`
41009cf685aSAlexis Marboeuf 
41109cf685aSAlexis Marboeuf   Level: intermediate
41209cf685aSAlexis Marboeuf 
413a1cb98faSBarry Smith   Note:
414a1cb98faSBarry Smith   The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`.
415a1cb98faSBarry Smith 
416*1cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()`
41709cf685aSAlexis Marboeuf  @*/
418d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreateNaturalVector(DM dm, Vec *nv)
419d71ae5a4SJacob Faibussowitsch {
42009cf685aSAlexis Marboeuf   PetscMPIInt size;
42109cf685aSAlexis Marboeuf 
42209cf685aSAlexis Marboeuf   PetscFunctionBegin;
42309cf685aSAlexis Marboeuf   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
42409cf685aSAlexis Marboeuf   if (dm->sfNatural) {
4256fba1fe0SBlaise Bourdin     PetscInt nleaves, bs, maxbs;
4262e8d78feSBlaise Bourdin     Vec      v;
4276fba1fe0SBlaise Bourdin 
4286fba1fe0SBlaise Bourdin     /*
4296fba1fe0SBlaise Bourdin       Setting the natural vector block size.
4306fba1fe0SBlaise Bourdin       We can't get it from a global vector because of constraints, and the block size in the local vector
4316fba1fe0SBlaise Bourdin       may be inconsistent across processes, typically when some local vectors have size 0, their block size is set to 1
4326fba1fe0SBlaise Bourdin     */
4332e8d78feSBlaise Bourdin     PetscCall(DMGetLocalVector(dm, &v));
4342e8d78feSBlaise Bourdin     PetscCall(VecGetBlockSize(v, &bs));
435712fec58SPierre Jolivet     PetscCall(MPIU_Allreduce(&bs, &maxbs, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm)));
4366fba1fe0SBlaise Bourdin     if (bs == 1 && maxbs > 1) bs = maxbs;
4372e8d78feSBlaise Bourdin     PetscCall(DMRestoreLocalVector(dm, &v));
4382e8d78feSBlaise Bourdin 
43909cf685aSAlexis Marboeuf     PetscCall(PetscSFGetGraph(dm->sfNatural, NULL, &nleaves, NULL, NULL));
44009cf685aSAlexis Marboeuf     PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), nv));
44109cf685aSAlexis Marboeuf     PetscCall(VecSetSizes(*nv, nleaves, PETSC_DETERMINE));
44209cf685aSAlexis Marboeuf     PetscCall(VecSetBlockSize(*nv, bs));
44309cf685aSAlexis Marboeuf     PetscCall(VecSetType(*nv, dm->vectype));
44409cf685aSAlexis Marboeuf     PetscCall(VecSetDM(*nv, dm));
44509cf685aSAlexis Marboeuf   } else if (size == 1) {
4462e8d78feSBlaise Bourdin     PetscCall(DMCreateLocalVector(dm, nv));
44709cf685aSAlexis Marboeuf   } else {
44809cf685aSAlexis 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.");
44909cf685aSAlexis Marboeuf     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().");
45009cf685aSAlexis Marboeuf   }
4513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
45209cf685aSAlexis Marboeuf }
453