xref: /petsc/src/dm/impls/plex/plexnatural.c (revision eb9d3e4d1b0d07ce7eca2be9429d4698ffa7ae7f)
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`
1060225df5SJacob Faibussowitsch - migrationSF - 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 
171cc06b55SBarry 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;
22a987ce93SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
23a987ce93SMatthew G. Knepley   if (migrationSF) PetscValidHeaderSpecific(migrationSF, PETSCSF_CLASSID, 2);
249566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)migrationSF));
25a987ce93SMatthew G. Knepley   PetscCall(PetscSFDestroy(&dm->sfMigration));
26a987ce93SMatthew G. Knepley   dm->sfMigration = migrationSF;
273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
28f94b4a02SBlaise Bourdin }
29f94b4a02SBlaise Bourdin 
30f94b4a02SBlaise Bourdin /*@
31a1cb98faSBarry Smith   DMPlexGetMigrationSF - Gets the `PetscSF` for migrating from a parent `DM` into this `DM`
32a1cb98faSBarry Smith 
33a1cb98faSBarry Smith   Note Collective
34f94b4a02SBlaise Bourdin 
355d3b26e6SMatthew G. Knepley   Input Parameter:
36a1cb98faSBarry Smith . dm - The `DM`
375d3b26e6SMatthew G. Knepley 
385d3b26e6SMatthew G. Knepley   Output Parameter:
39a1cb98faSBarry Smith . migrationSF - The `PetscSF`
405d3b26e6SMatthew G. Knepley 
41f94b4a02SBlaise Bourdin   Level: intermediate
42f94b4a02SBlaise Bourdin 
431cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `PetscSF`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexCreateMigrationSF()`, `DMPlexSetMigrationSF`
44f94b4a02SBlaise Bourdin @*/
45d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGetMigrationSF(DM dm, PetscSF *migrationSF)
46d71ae5a4SJacob Faibussowitsch {
47f94b4a02SBlaise Bourdin   PetscFunctionBegin;
48f94b4a02SBlaise Bourdin   *migrationSF = dm->sfMigration;
493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
50f94b4a02SBlaise Bourdin }
51f94b4a02SBlaise Bourdin 
52f94b4a02SBlaise Bourdin /*@
53a1cb98faSBarry Smith   DMPlexSetGlobalToNaturalSF - Sets the `PetscSF` for mapping Global `Vec` to the Natural `Vec`
54f94b4a02SBlaise Bourdin 
555d3b26e6SMatthew G. Knepley   Input Parameters:
56a1cb98faSBarry Smith + dm        - The `DM`
57a1cb98faSBarry Smith - naturalSF - The `PetscSF`
585d3b26e6SMatthew G. Knepley 
59f94b4a02SBlaise Bourdin   Level: intermediate
60f94b4a02SBlaise Bourdin 
611cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `PetscSF`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexCreateGlobalToNaturalSF()`, `DMPlexGetGlobaltoNaturalSF()`
62f94b4a02SBlaise Bourdin @*/
63d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexSetGlobalToNaturalSF(DM dm, PetscSF naturalSF)
64d71ae5a4SJacob Faibussowitsch {
65f94b4a02SBlaise Bourdin   PetscFunctionBegin;
66f94b4a02SBlaise Bourdin   dm->sfNatural = naturalSF;
679566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)naturalSF));
686925d1c4SMatthew G. Knepley   dm->useNatural = naturalSF ? PETSC_TRUE : PETSC_FALSE;
693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
70f94b4a02SBlaise Bourdin }
71f94b4a02SBlaise Bourdin 
72f94b4a02SBlaise Bourdin /*@
73a1cb98faSBarry Smith   DMPlexGetGlobalToNaturalSF - Gets the `PetscSF` for mapping Global `Vec` to the Natural `Vec`
74f94b4a02SBlaise Bourdin 
755d3b26e6SMatthew G. Knepley   Input Parameter:
76a1cb98faSBarry Smith . dm - The `DM`
775d3b26e6SMatthew G. Knepley 
785d3b26e6SMatthew G. Knepley   Output Parameter:
79a1cb98faSBarry Smith . naturalSF - The `PetscSF`
805d3b26e6SMatthew G. Knepley 
81f94b4a02SBlaise Bourdin   Level: intermediate
82f94b4a02SBlaise Bourdin 
831cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `PetscSF`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexCreateGlobalToNaturalSF()`, `DMPlexSetGlobaltoNaturalSF`
84f94b4a02SBlaise Bourdin @*/
85d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGetGlobalToNaturalSF(DM dm, PetscSF *naturalSF)
86d71ae5a4SJacob Faibussowitsch {
87f94b4a02SBlaise Bourdin   PetscFunctionBegin;
88f94b4a02SBlaise Bourdin   *naturalSF = dm->sfNatural;
893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
90f94b4a02SBlaise Bourdin }
91f94b4a02SBlaise Bourdin 
92f94b4a02SBlaise Bourdin /*@
93a1cb98faSBarry Smith   DMPlexCreateGlobalToNaturalSF - Creates the `PetscSF` for mapping Global `Vec` to the Natural `Vec`
94fa534816SMatthew G. Knepley 
95fa534816SMatthew G. Knepley   Input Parameters:
96a1cb98faSBarry Smith + dm          - The redistributed `DM`
9720f4b53cSBarry Smith . section     - The local `PetscSection` describing the `Vec` before the mesh was distributed, or `NULL` if not available
9820f4b53cSBarry Smith - sfMigration - The `PetscSF` used to distribute the mesh, or `NULL` if it cannot be computed
99fa534816SMatthew G. Knepley 
1005d3b26e6SMatthew G. Knepley   Output Parameter:
101a1cb98faSBarry Smith . sfNatural - `PetscSF` for mapping the `Vec` in PETSc ordering to the canonical ordering
1025d3b26e6SMatthew G. Knepley 
103fa534816SMatthew G. Knepley   Level: intermediate
104fa534816SMatthew G. Knepley 
105a1cb98faSBarry Smith   Note:
106a1cb98faSBarry Smith   This is not typically called by the user.
107a1cb98faSBarry Smith 
1081cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `PetscSF`, `PetscSection`, `DMPlexDistribute()`, `DMPlexDistributeField()`
109fa534816SMatthew G. Knepley  @*/
110d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreateGlobalToNaturalSF(DM dm, PetscSection section, PetscSF sfMigration, PetscSF *sfNatural)
111d71ae5a4SJacob Faibussowitsch {
112fa534816SMatthew G. Knepley   MPI_Comm     comm;
113e7858701SMatthew G. Knepley   PetscSF      sf, sfEmbed, sfField;
114e5b44f4fSMatthew G. Knepley   PetscSection gSection, sectionDist, gLocSection;
115fa534816SMatthew G. Knepley   PetscInt    *spoints, *remoteOffsets;
11657a2b1f8SBlaise Bourdin   PetscInt     ssize, pStart, pEnd, p, localSize, maxStorageSize;
117e7858701SMatthew G. Knepley   PetscBool    destroyFlag = PETSC_FALSE, debug = PETSC_FALSE;
118fa534816SMatthew G. Knepley 
119fa534816SMatthew G. Knepley   PetscFunctionBegin;
1209566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
1215c3f5608SAlexisMarb   if (!sfMigration) {
122e7858701SMatthew G. Knepley     /* If sfMigration is missing, sfNatural cannot be computed and is set to NULL */
1235c3f5608SAlexisMarb     *sfNatural = NULL;
1243ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
1255c3f5608SAlexisMarb   } else if (!section) {
126e7858701SMatthew G. Knepley     /* If the sequential section is not provided (NULL), it is reconstructed from the parallel section */
1275c3f5608SAlexisMarb     PetscSF      sfMigrationInv;
1285c3f5608SAlexisMarb     PetscSection localSection;
1295c3f5608SAlexisMarb 
1309566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dm, &localSection));
1319566063dSJacob Faibussowitsch     PetscCall(PetscSFCreateInverseSF(sfMigration, &sfMigrationInv));
1329566063dSJacob Faibussowitsch     PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), &section));
1339566063dSJacob Faibussowitsch     PetscCall(PetscSFDistributeSection(sfMigrationInv, localSection, NULL, section));
1349566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&sfMigrationInv));
1355c3f5608SAlexisMarb     destroyFlag = PETSC_TRUE;
1365c3f5608SAlexisMarb   }
137e7858701SMatthew G. Knepley   if (debug) PetscCall(PetscSFView(sfMigration, NULL));
138e5b44f4fSMatthew G. Knepley   /* Create a new section from distributing the original section */
1399566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(comm, &sectionDist));
1409566063dSJacob Faibussowitsch   PetscCall(PetscSFDistributeSection(sfMigration, section, &remoteOffsets, sectionDist));
141e7858701SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)sectionDist, "Migrated Section"));
142e7858701SMatthew G. Knepley   if (debug) PetscCall(PetscSectionView(sectionDist, NULL));
1439566063dSJacob Faibussowitsch   PetscCall(DMSetLocalSection(dm, sectionDist));
144e7858701SMatthew G. Knepley   /* If a sequential section is provided but no dof is affected, sfNatural cannot be computed and is set to NULL */
14557a2b1f8SBlaise Bourdin   PetscCall(PetscSectionGetStorageSize(sectionDist, &localSize));
146ab9fe21fSBlaise Bourdin   PetscCall(MPIU_Allreduce(&localSize, &maxStorageSize, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm)));
14757a2b1f8SBlaise Bourdin   if (maxStorageSize) {
148e7858701SMatthew G. Knepley     const PetscInt *leaves;
149e7858701SMatthew G. Knepley     PetscInt       *sortleaves, *indices;
150e7858701SMatthew G. Knepley     PetscInt        Nl;
151e7858701SMatthew G. Knepley 
152fa534816SMatthew G. Knepley     /* Get a pruned version of migration SF */
1539566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalSection(dm, &gSection));
154e7858701SMatthew G. Knepley     if (debug) PetscCall(PetscSectionView(gSection, NULL));
155e7858701SMatthew G. Knepley     PetscCall(PetscSFGetGraph(sfMigration, NULL, &Nl, &leaves, NULL));
1569566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetChart(gSection, &pStart, &pEnd));
157fa534816SMatthew G. Knepley     for (p = pStart, ssize = 0; p < pEnd; ++p) {
158fa534816SMatthew G. Knepley       PetscInt dof, off;
159fa534816SMatthew G. Knepley 
1609566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(gSection, p, &dof));
1619566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(gSection, p, &off));
162fa534816SMatthew G. Knepley       if ((dof > 0) && (off >= 0)) ++ssize;
163fa534816SMatthew G. Knepley     }
164e7858701SMatthew G. Knepley     PetscCall(PetscMalloc3(ssize, &spoints, Nl, &sortleaves, Nl, &indices));
1659371c9d4SSatish Balay     for (p = 0; p < Nl; ++p) {
1669371c9d4SSatish Balay       sortleaves[p] = leaves ? leaves[p] : p;
1679371c9d4SSatish Balay       indices[p]    = p;
1689371c9d4SSatish Balay     }
169e7858701SMatthew G. Knepley     PetscCall(PetscSortIntWithArray(Nl, sortleaves, indices));
170fa534816SMatthew G. Knepley     for (p = pStart, ssize = 0; p < pEnd; ++p) {
171e7858701SMatthew G. Knepley       PetscInt dof, off, loc;
172fa534816SMatthew G. Knepley 
1739566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(gSection, p, &dof));
1749566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(gSection, p, &off));
175e7858701SMatthew G. Knepley       if ((dof > 0) && (off >= 0)) {
176e7858701SMatthew G. Knepley         PetscCall(PetscFindInt(p, Nl, sortleaves, &loc));
177e7858701SMatthew 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);
178e7858701SMatthew G. Knepley         spoints[ssize++] = indices[loc];
179e7858701SMatthew G. Knepley       }
180fa534816SMatthew G. Knepley     }
1819566063dSJacob Faibussowitsch     PetscCall(PetscSFCreateEmbeddedLeafSF(sfMigration, ssize, spoints, &sfEmbed));
182e7858701SMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)sfEmbed, "Embedded SF"));
183e7858701SMatthew G. Knepley     PetscCall(PetscFree3(spoints, sortleaves, indices));
184e7858701SMatthew G. Knepley     if (debug) PetscCall(PetscSFView(sfEmbed, NULL));
185e7858701SMatthew G. Knepley     /* Create the SF associated with this section
186e7858701SMatthew G. Knepley          Roots are natural dofs, leaves are global dofs */
1879566063dSJacob Faibussowitsch     PetscCall(DMGetPointSF(dm, &sf));
188*eb9d3e4dSMatthew G. Knepley     PetscCall(PetscSectionCreateGlobalSection(sectionDist, sf, PETSC_TRUE, PETSC_TRUE, PETSC_TRUE, &gLocSection));
1899566063dSJacob Faibussowitsch     PetscCall(PetscSFCreateSectionSF(sfEmbed, section, remoteOffsets, gLocSection, &sfField));
1909566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&sfEmbed));
1919566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&gLocSection));
192e7858701SMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)sfField, "Natural-to-Global SF"));
193e7858701SMatthew G. Knepley     if (debug) PetscCall(PetscSFView(sfField, NULL));
194e7858701SMatthew G. Knepley     /* Invert the field SF
195e7858701SMatthew G. Knepley          Roots are global dofs, leaves are natural dofs */
196e7858701SMatthew G. Knepley     PetscCall(PetscSFCreateInverseSF(sfField, sfNatural));
197e7858701SMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)*sfNatural, "Global-to-Natural SF"));
1989566063dSJacob Faibussowitsch     PetscCall(PetscObjectViewFromOptions((PetscObject)*sfNatural, NULL, "-globaltonatural_sf_view"));
199fa534816SMatthew G. Knepley     /* Clean up */
200e7858701SMatthew G. Knepley     PetscCall(PetscSFDestroy(&sfField));
20114744f87SBlaise Bourdin   } else {
20214744f87SBlaise Bourdin     *sfNatural = NULL;
20314744f87SBlaise Bourdin   }
2049566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&sectionDist));
2059566063dSJacob Faibussowitsch   PetscCall(PetscFree(remoteOffsets));
2069566063dSJacob Faibussowitsch   if (destroyFlag) PetscCall(PetscSectionDestroy(&section));
2073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
208fa534816SMatthew G. Knepley }
209fa534816SMatthew G. Knepley 
210fa534816SMatthew G. Knepley /*@
211a1cb98faSBarry Smith   DMPlexGlobalToNaturalBegin - Rearranges a global `Vec` in the natural order.
212fa534816SMatthew G. Knepley 
21320f4b53cSBarry Smith   Collective
214fa534816SMatthew G. Knepley 
215fa534816SMatthew G. Knepley   Input Parameters:
216a1cb98faSBarry Smith + dm - The distributed `DMPLEX`
217a1cb98faSBarry Smith - gv - The global `Vec`
218fa534816SMatthew G. Knepley 
2192fe279fdSBarry Smith   Output Parameter:
22020f4b53cSBarry Smith . nv - `Vec` in the canonical ordering distributed over all processors associated with `gv`
221fa534816SMatthew G. Knepley 
222fa534816SMatthew G. Knepley   Level: intermediate
223fa534816SMatthew G. Knepley 
224a1cb98faSBarry Smith   Note:
225a1cb98faSBarry Smith   The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`.
226a1cb98faSBarry Smith 
2271cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalEnd()`
228fa534816SMatthew G. Knepley @*/
229d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGlobalToNaturalBegin(DM dm, Vec gv, Vec nv)
230d71ae5a4SJacob Faibussowitsch {
231fa534816SMatthew G. Knepley   const PetscScalar *inarray;
232fa534816SMatthew G. Knepley   PetscScalar       *outarray;
233f6c7c579SMatthew G. Knepley   MPI_Comm           comm;
234a6a55facSMatthew G. Knepley   PetscMPIInt        size;
235fa534816SMatthew G. Knepley 
236fa534816SMatthew G. Knepley   PetscFunctionBegin;
2379566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_GlobalToNaturalBegin, dm, 0, 0, 0));
238f6c7c579SMatthew G. Knepley   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
239f6c7c579SMatthew G. Knepley   PetscCallMPI(MPI_Comm_size(comm, &size));
240fa534816SMatthew G. Knepley   if (dm->sfNatural) {
241f6c7c579SMatthew G. Knepley     if (PetscDefined(USE_DEBUG)) {
242f6c7c579SMatthew G. Knepley       PetscSection gs;
243f6c7c579SMatthew G. Knepley       PetscInt     Nl, n;
244f6c7c579SMatthew G. Knepley 
245f6c7c579SMatthew G. Knepley       PetscCall(PetscSFGetGraph(dm->sfNatural, NULL, &Nl, NULL, NULL));
246f6c7c579SMatthew G. Knepley       PetscCall(VecGetLocalSize(nv, &n));
247f6c7c579SMatthew 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);
248f6c7c579SMatthew G. Knepley 
249f6c7c579SMatthew G. Knepley       PetscCall(DMGetGlobalSection(dm, &gs));
250f6c7c579SMatthew G. Knepley       PetscCall(PetscSectionGetConstrainedStorageSize(gs, &Nl));
251f6c7c579SMatthew G. Knepley       PetscCall(VecGetLocalSize(gv, &n));
252f6c7c579SMatthew 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);
253f6c7c579SMatthew G. Knepley     }
2549566063dSJacob Faibussowitsch     PetscCall(VecGetArray(nv, &outarray));
2559566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(gv, &inarray));
2569566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_REPLACE));
2579566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(gv, &inarray));
2589566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(nv, &outarray));
259a6a55facSMatthew G. Knepley   } else if (size == 1) {
2609566063dSJacob Faibussowitsch     PetscCall(VecCopy(gv, nv));
261f7d195e4SLawrence Mitchell   } else {
262f7d195e4SLawrence 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.");
263f7d195e4SLawrence Mitchell     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().");
264f7d195e4SLawrence Mitchell   }
2659566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_GlobalToNaturalBegin, dm, 0, 0, 0));
2663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
267fa534816SMatthew G. Knepley }
268fa534816SMatthew G. Knepley 
269fa534816SMatthew G. Knepley /*@
270a1cb98faSBarry Smith   DMPlexGlobalToNaturalEnd - Rearranges a global `Vec` in the natural order.
271fa534816SMatthew G. Knepley 
27220f4b53cSBarry Smith   Collective
273fa534816SMatthew G. Knepley 
274fa534816SMatthew G. Knepley   Input Parameters:
275a1cb98faSBarry Smith + dm - The distributed `DMPLEX`
276a1cb98faSBarry Smith - gv - The global `Vec`
277fa534816SMatthew G. Knepley 
27897bb3fdcSJose E. Roman   Output Parameter:
279a1cb98faSBarry Smith . nv - The natural `Vec`
280fa534816SMatthew G. Knepley 
281fa534816SMatthew G. Knepley   Level: intermediate
282fa534816SMatthew G. Knepley 
283a1cb98faSBarry Smith   Note:
284a1cb98faSBarry Smith   The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`.
285a1cb98faSBarry Smith 
2861cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()`
287fa534816SMatthew G. Knepley  @*/
288d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGlobalToNaturalEnd(DM dm, Vec gv, Vec nv)
289d71ae5a4SJacob Faibussowitsch {
290fa534816SMatthew G. Knepley   const PetscScalar *inarray;
291fa534816SMatthew G. Knepley   PetscScalar       *outarray;
292a6a55facSMatthew G. Knepley   PetscMPIInt        size;
293fa534816SMatthew G. Knepley 
294fa534816SMatthew G. Knepley   PetscFunctionBegin;
2959566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_GlobalToNaturalEnd, dm, 0, 0, 0));
2969566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
297fa534816SMatthew G. Knepley   if (dm->sfNatural) {
2989566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(gv, &inarray));
2999566063dSJacob Faibussowitsch     PetscCall(VecGetArray(nv, &outarray));
3009566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_REPLACE));
3019566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(gv, &inarray));
3029566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(nv, &outarray));
303a6a55facSMatthew G. Knepley   } else if (size == 1) {
304f7d195e4SLawrence Mitchell   } else {
305f7d195e4SLawrence 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.");
306f7d195e4SLawrence Mitchell     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().");
307f7d195e4SLawrence Mitchell   }
3089566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_GlobalToNaturalEnd, dm, 0, 0, 0));
3093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
310fa534816SMatthew G. Knepley }
311fa534816SMatthew G. Knepley 
312fa534816SMatthew G. Knepley /*@
313a1cb98faSBarry Smith   DMPlexNaturalToGlobalBegin - Rearranges a `Vec` in the natural order to the Global order.
314fa534816SMatthew G. Knepley 
31520f4b53cSBarry Smith   Collective
316fa534816SMatthew G. Knepley 
317fa534816SMatthew G. Knepley   Input Parameters:
318a1cb98faSBarry Smith + dm - The distributed `DMPLEX`
319a1cb98faSBarry Smith - nv - The natural `Vec`
320fa534816SMatthew G. Knepley 
3212fe279fdSBarry Smith   Output Parameter:
322a1cb98faSBarry Smith . gv - The global `Vec`
323fa534816SMatthew G. Knepley 
324fa534816SMatthew G. Knepley   Level: intermediate
325fa534816SMatthew G. Knepley 
326a1cb98faSBarry Smith   Note:
327a1cb98faSBarry Smith   The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`.
328a1cb98faSBarry Smith 
32942747ad1SJacob Faibussowitsch .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexGlobalToNaturalEnd()`
330fa534816SMatthew G. Knepley @*/
331d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexNaturalToGlobalBegin(DM dm, Vec nv, Vec gv)
332d71ae5a4SJacob Faibussowitsch {
333fa534816SMatthew G. Knepley   const PetscScalar *inarray;
334fa534816SMatthew G. Knepley   PetscScalar       *outarray;
335a6a55facSMatthew G. Knepley   PetscMPIInt        size;
336fa534816SMatthew G. Knepley 
337fa534816SMatthew G. Knepley   PetscFunctionBegin;
3389566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_NaturalToGlobalBegin, dm, 0, 0, 0));
3399566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
340fa534816SMatthew G. Knepley   if (dm->sfNatural) {
341a5b23f4aSJose E. Roman     /* We only have access to the SF that goes from Global to Natural.
342b64f75a9SBlaise Bourdin        Instead of inverting dm->sfNatural, we can call PetscSFReduceBegin/End with MPI_Op MPI_SUM.
343b64f75a9SBlaise Bourdin        Here the SUM really does nothing since sfNatural is one to one, as long as gV is set to zero first. */
3449566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(gv));
3459566063dSJacob Faibussowitsch     PetscCall(VecGetArray(gv, &outarray));
3469566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(nv, &inarray));
3479566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_SUM));
3489566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(nv, &inarray));
3499566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(gv, &outarray));
350a6a55facSMatthew G. Knepley   } else if (size == 1) {
3519566063dSJacob Faibussowitsch     PetscCall(VecCopy(nv, gv));
352f7d195e4SLawrence Mitchell   } else {
353f7d195e4SLawrence 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.");
354f7d195e4SLawrence Mitchell     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().");
355f7d195e4SLawrence Mitchell   }
3569566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_NaturalToGlobalBegin, dm, 0, 0, 0));
3573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
358fa534816SMatthew G. Knepley }
359fa534816SMatthew G. Knepley 
360fa534816SMatthew G. Knepley /*@
361a1cb98faSBarry Smith   DMPlexNaturalToGlobalEnd - Rearranges a `Vec` in the natural order to the Global order.
362fa534816SMatthew G. Knepley 
36320f4b53cSBarry Smith   Collective
364fa534816SMatthew G. Knepley 
365fa534816SMatthew G. Knepley   Input Parameters:
366a1cb98faSBarry Smith + dm - The distributed `DMPLEX`
367a1cb98faSBarry Smith - nv - The natural `Vec`
368fa534816SMatthew G. Knepley 
3692fe279fdSBarry Smith   Output Parameter:
370a1cb98faSBarry Smith . gv - The global `Vec`
371fa534816SMatthew G. Knepley 
372fa534816SMatthew G. Knepley   Level: intermediate
373fa534816SMatthew G. Knepley 
374a1cb98faSBarry Smith   Note:
375a1cb98faSBarry Smith   The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`.
376a1cb98faSBarry Smith 
3771cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()`
378fa534816SMatthew G. Knepley  @*/
379d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexNaturalToGlobalEnd(DM dm, Vec nv, Vec gv)
380d71ae5a4SJacob Faibussowitsch {
381fa534816SMatthew G. Knepley   const PetscScalar *inarray;
382fa534816SMatthew G. Knepley   PetscScalar       *outarray;
383a6a55facSMatthew G. Knepley   PetscMPIInt        size;
384fa534816SMatthew G. Knepley 
385fa534816SMatthew G. Knepley   PetscFunctionBegin;
3869566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_NaturalToGlobalEnd, dm, 0, 0, 0));
3879566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
388fa534816SMatthew G. Knepley   if (dm->sfNatural) {
3899566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(nv, &inarray));
3909566063dSJacob Faibussowitsch     PetscCall(VecGetArray(gv, &outarray));
3919566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_SUM));
3929566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(nv, &inarray));
3939566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(gv, &outarray));
394a6a55facSMatthew G. Knepley   } else if (size == 1) {
395f7d195e4SLawrence Mitchell   } else {
396f7d195e4SLawrence 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.");
397f7d195e4SLawrence Mitchell     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().");
398f7d195e4SLawrence Mitchell   }
3999566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_NaturalToGlobalEnd, dm, 0, 0, 0));
4003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
401fa534816SMatthew G. Knepley }
40209cf685aSAlexis Marboeuf 
40309cf685aSAlexis Marboeuf /*@
4045cb80ecdSBarry Smith   DMPlexCreateNaturalVector - Provide a `Vec` capable of holding the natural ordering and distribution.
40509cf685aSAlexis Marboeuf 
40620f4b53cSBarry Smith   Collective
40709cf685aSAlexis Marboeuf 
4085cb80ecdSBarry Smith   Input Parameter:
409a1cb98faSBarry Smith . dm - The distributed `DMPLEX`
41009cf685aSAlexis Marboeuf 
4115cb80ecdSBarry Smith   Output Parameter:
4125cb80ecdSBarry Smith . nv - The natural `Vec`
41309cf685aSAlexis Marboeuf 
41409cf685aSAlexis Marboeuf   Level: intermediate
41509cf685aSAlexis Marboeuf 
416a1cb98faSBarry Smith   Note:
417a1cb98faSBarry Smith   The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`.
418a1cb98faSBarry Smith 
4191cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()`
42009cf685aSAlexis Marboeuf  @*/
421d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreateNaturalVector(DM dm, Vec *nv)
422d71ae5a4SJacob Faibussowitsch {
42309cf685aSAlexis Marboeuf   PetscMPIInt size;
42409cf685aSAlexis Marboeuf 
42509cf685aSAlexis Marboeuf   PetscFunctionBegin;
42609cf685aSAlexis Marboeuf   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
42709cf685aSAlexis Marboeuf   if (dm->sfNatural) {
4286fba1fe0SBlaise Bourdin     PetscInt nleaves, bs, maxbs;
4292e8d78feSBlaise Bourdin     Vec      v;
4306fba1fe0SBlaise Bourdin 
4316fba1fe0SBlaise Bourdin     /*
4326fba1fe0SBlaise Bourdin       Setting the natural vector block size.
4336fba1fe0SBlaise Bourdin       We can't get it from a global vector because of constraints, and the block size in the local vector
4346fba1fe0SBlaise Bourdin       may be inconsistent across processes, typically when some local vectors have size 0, their block size is set to 1
4356fba1fe0SBlaise Bourdin     */
4362e8d78feSBlaise Bourdin     PetscCall(DMGetLocalVector(dm, &v));
4372e8d78feSBlaise Bourdin     PetscCall(VecGetBlockSize(v, &bs));
438712fec58SPierre Jolivet     PetscCall(MPIU_Allreduce(&bs, &maxbs, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm)));
4396fba1fe0SBlaise Bourdin     if (bs == 1 && maxbs > 1) bs = maxbs;
4402e8d78feSBlaise Bourdin     PetscCall(DMRestoreLocalVector(dm, &v));
4412e8d78feSBlaise Bourdin 
44209cf685aSAlexis Marboeuf     PetscCall(PetscSFGetGraph(dm->sfNatural, NULL, &nleaves, NULL, NULL));
44309cf685aSAlexis Marboeuf     PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), nv));
44409cf685aSAlexis Marboeuf     PetscCall(VecSetSizes(*nv, nleaves, PETSC_DETERMINE));
44509cf685aSAlexis Marboeuf     PetscCall(VecSetBlockSize(*nv, bs));
44609cf685aSAlexis Marboeuf     PetscCall(VecSetType(*nv, dm->vectype));
44709cf685aSAlexis Marboeuf     PetscCall(VecSetDM(*nv, dm));
44809cf685aSAlexis Marboeuf   } else if (size == 1) {
4492e8d78feSBlaise Bourdin     PetscCall(DMCreateLocalVector(dm, nv));
45009cf685aSAlexis Marboeuf   } else {
45109cf685aSAlexis 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.");
45209cf685aSAlexis Marboeuf     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().");
45309cf685aSAlexis Marboeuf   }
4543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
45509cf685aSAlexis Marboeuf }
456