xref: /petsc/src/dm/impls/plex/plexnatural.c (revision 462c564db5c15796f1f07c2b71ff29812cdd1aad)
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));
146*462c564dSBarry Smith   PetscCallMPI(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));
188eb9d3e4dSMatthew 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 /*@
21116635c05SJames Wright   DMPlexMigrateGlobalToNaturalSF - Migrates the input `sfNatural` based on sfMigration
21216635c05SJames Wright 
21316635c05SJames Wright   Input Parameters:
21416635c05SJames Wright + dmOld        - The original `DM`
21516635c05SJames Wright . dmNew        - The `DM` to be migrated to
21616635c05SJames Wright . sfNaturalOld - The sfNatural for the `dmOld`
21716635c05SJames Wright - sfMigration  - The `PetscSF` used to distribute the mesh, or `NULL` if it cannot be computed
21816635c05SJames Wright 
21916635c05SJames Wright   Output Parameter:
22016635c05SJames Wright . sfNaturalNew - `PetscSF` for mapping the `Vec` in PETSc ordering to the canonical ordering
22116635c05SJames Wright 
22216635c05SJames Wright   Level: intermediate
22316635c05SJames Wright 
22416635c05SJames Wright   Notes:
22516635c05SJames Wright   `sfNaturalOld` maps from the old Global section (roots) to the natural Vec layout (leaves, may or may not be described by a PetscSection).
22616635c05SJames Wright   `DMPlexMigrateGlobalToNaturalSF` creates an SF to map from the old global section to the new global section (generated from `sfMigration`).
22716635c05SJames Wright   That SF is then composed with the `sfNaturalOld` to generate `sfNaturalNew`.
22816635c05SJames Wright   This also distributes and sets the local section for `dmNew`.
22916635c05SJames Wright 
23016635c05SJames Wright   This is not typically called by the user.
23116635c05SJames Wright 
23216635c05SJames Wright .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `PetscSF`, `PetscSection`, `DMPlexDistribute()`, `DMPlexDistributeField()`
23316635c05SJames Wright  @*/
23416635c05SJames Wright PetscErrorCode DMPlexMigrateGlobalToNaturalSF(DM dmOld, DM dmNew, PetscSF sfNaturalOld, PetscSF sfMigration, PetscSF *sfNaturalNew)
23516635c05SJames Wright {
23616635c05SJames Wright   MPI_Comm     comm;
23716635c05SJames Wright   PetscSection oldGlobalSection;
23816635c05SJames Wright   PetscSection newGlobalSection;
23916635c05SJames Wright   PetscInt    *remoteOffsets;
24016635c05SJames Wright   PetscBool    debug = PETSC_FALSE;
24116635c05SJames Wright 
24216635c05SJames Wright   PetscFunctionBegin;
24316635c05SJames Wright   PetscCall(PetscObjectGetComm((PetscObject)dmNew, &comm));
24416635c05SJames Wright   if (!sfMigration) {
24516635c05SJames Wright     /* If sfMigration is missing, sfNatural cannot be computed and is set to NULL */
24616635c05SJames Wright     *sfNaturalNew = NULL;
24716635c05SJames Wright     PetscFunctionReturn(PETSC_SUCCESS);
24816635c05SJames Wright   }
24916635c05SJames Wright   if (debug) PetscCall(PetscSFView(sfMigration, NULL));
25016635c05SJames Wright 
25116635c05SJames Wright   { // Create oldGlobalSection and newGlobalSection *with* localOffsets
25216635c05SJames Wright     PetscSection oldLocalSection, newLocalSection;
25316635c05SJames Wright     PetscSF      pointSF;
25416635c05SJames Wright 
25516635c05SJames Wright     PetscCall(DMGetLocalSection(dmOld, &oldLocalSection));
25616635c05SJames Wright     PetscCall(DMGetPointSF(dmOld, &pointSF));
25716635c05SJames Wright     PetscCall(PetscSectionCreateGlobalSection(oldLocalSection, pointSF, PETSC_TRUE, PETSC_TRUE, PETSC_TRUE, &oldGlobalSection));
25816635c05SJames Wright 
25916635c05SJames Wright     PetscCall(PetscSectionCreate(comm, &newLocalSection));
26016635c05SJames Wright     PetscCall(PetscSFDistributeSection(sfMigration, oldLocalSection, NULL, newLocalSection));
26116635c05SJames Wright     PetscCall(DMSetLocalSection(dmNew, newLocalSection));
26216635c05SJames Wright 
26316635c05SJames Wright     PetscCall(PetscSectionCreate(comm, &newGlobalSection));
26416635c05SJames Wright     PetscCall(DMGetPointSF(dmNew, &pointSF));
26516635c05SJames Wright     PetscCall(PetscSectionCreateGlobalSection(newLocalSection, pointSF, PETSC_TRUE, PETSC_TRUE, PETSC_TRUE, &newGlobalSection));
26616635c05SJames Wright 
26716635c05SJames Wright     PetscCall(PetscObjectSetName((PetscObject)oldLocalSection, "Old Local Section"));
26816635c05SJames Wright     if (debug) PetscCall(PetscSectionView(oldLocalSection, NULL));
26916635c05SJames Wright     PetscCall(PetscObjectSetName((PetscObject)oldGlobalSection, "Old Global Section"));
27016635c05SJames Wright     if (debug) PetscCall(PetscSectionView(oldGlobalSection, NULL));
27116635c05SJames Wright     PetscCall(PetscObjectSetName((PetscObject)newLocalSection, "New Local Section"));
27216635c05SJames Wright     if (debug) PetscCall(PetscSectionView(newLocalSection, NULL));
27316635c05SJames Wright     PetscCall(PetscObjectSetName((PetscObject)newGlobalSection, "New Global Section"));
27416635c05SJames Wright     if (debug) PetscCall(PetscSectionView(newGlobalSection, NULL));
27516635c05SJames Wright     PetscCall(PetscSectionDestroy(&newLocalSection));
27616635c05SJames Wright   }
27716635c05SJames Wright 
27816635c05SJames Wright   { // Create remoteOffsets array, mapping the oldGlobalSection offsets to the local points (according to sfMigration)
27916635c05SJames Wright     PetscInt lpStart, lpEnd, rpStart, rpEnd;
28016635c05SJames Wright 
28116635c05SJames Wright     PetscCall(PetscSectionGetChart(oldGlobalSection, &rpStart, &rpEnd));
28216635c05SJames Wright     PetscCall(PetscSectionGetChart(newGlobalSection, &lpStart, &lpEnd));
28316635c05SJames Wright 
28416635c05SJames Wright     // in `PetscSFDistributeSection` (where this is taken from), it possibly makes a new embedded SF. Should possibly do that here?
28516635c05SJames Wright     PetscCall(PetscMalloc1(lpEnd - lpStart, &remoteOffsets));
28616635c05SJames Wright     PetscCall(PetscSFBcastBegin(sfMigration, MPIU_INT, PetscSafePointerPlusOffset(oldGlobalSection->atlasOff, -rpStart), PetscSafePointerPlusOffset(remoteOffsets, -lpStart), MPI_REPLACE));
28716635c05SJames Wright     PetscCall(PetscSFBcastEnd(sfMigration, MPIU_INT, PetscSafePointerPlusOffset(oldGlobalSection->atlasOff, -rpStart), PetscSafePointerPlusOffset(remoteOffsets, -lpStart), MPI_REPLACE));
28816635c05SJames Wright     if (debug) {
28916635c05SJames Wright       PetscViewer viewer;
29016635c05SJames Wright 
29116635c05SJames Wright       PetscCall(PetscPrintf(comm, "Remote Offsets:\n"));
29216635c05SJames Wright       PetscCall(PetscViewerASCIIGetStdout(comm, &viewer));
29316635c05SJames Wright       PetscCall(PetscIntView(lpEnd - lpStart, remoteOffsets, viewer));
29416635c05SJames Wright     }
29516635c05SJames Wright   }
29616635c05SJames Wright 
29716635c05SJames Wright   { // Create SF from oldGlobalSection to newGlobalSection and compose with sfNaturalOld
29816635c05SJames Wright     PetscSF oldglob_to_newglob_sf, newglob_to_oldglob_sf;
29916635c05SJames Wright 
30016635c05SJames Wright     PetscCall(PetscSFCreateSectionSF(sfMigration, oldGlobalSection, remoteOffsets, newGlobalSection, &oldglob_to_newglob_sf));
30116635c05SJames Wright     PetscCall(PetscObjectSetName((PetscObject)oldglob_to_newglob_sf, "OldGlobal-to-NewGlobal SF"));
30216635c05SJames Wright     if (debug) PetscCall(PetscSFView(oldglob_to_newglob_sf, NULL));
30316635c05SJames Wright 
30416635c05SJames Wright     PetscCall(PetscSFCreateInverseSF(oldglob_to_newglob_sf, &newglob_to_oldglob_sf));
30516635c05SJames Wright     PetscCall(PetscObjectSetName((PetscObject)newglob_to_oldglob_sf, "NewGlobal-to-OldGlobal SF"));
30616635c05SJames Wright     PetscCall(PetscObjectViewFromOptions((PetscObject)newglob_to_oldglob_sf, (PetscObject)dmOld, "-natural_migrate_sf_view"));
30716635c05SJames Wright     PetscCall(PetscSFCompose(newglob_to_oldglob_sf, sfNaturalOld, sfNaturalNew));
30816635c05SJames Wright 
30916635c05SJames Wright     PetscCall(PetscSFDestroy(&oldglob_to_newglob_sf));
31016635c05SJames Wright     PetscCall(PetscSFDestroy(&newglob_to_oldglob_sf));
31116635c05SJames Wright   }
31216635c05SJames Wright 
31316635c05SJames Wright   PetscCall(PetscSectionDestroy(&oldGlobalSection));
31416635c05SJames Wright   PetscCall(PetscSectionDestroy(&newGlobalSection));
31516635c05SJames Wright   PetscCall(PetscFree(remoteOffsets));
31616635c05SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
31716635c05SJames Wright }
31816635c05SJames Wright 
31916635c05SJames Wright /*@
320a1cb98faSBarry Smith   DMPlexGlobalToNaturalBegin - Rearranges a global `Vec` in the natural order.
321fa534816SMatthew G. Knepley 
32220f4b53cSBarry Smith   Collective
323fa534816SMatthew G. Knepley 
324fa534816SMatthew G. Knepley   Input Parameters:
325a1cb98faSBarry Smith + dm - The distributed `DMPLEX`
326a1cb98faSBarry Smith - gv - The global `Vec`
327fa534816SMatthew G. Knepley 
3282fe279fdSBarry Smith   Output Parameter:
32920f4b53cSBarry Smith . nv - `Vec` in the canonical ordering distributed over all processors associated with `gv`
330fa534816SMatthew G. Knepley 
331fa534816SMatthew G. Knepley   Level: intermediate
332fa534816SMatthew G. Knepley 
333a1cb98faSBarry Smith   Note:
334a1cb98faSBarry Smith   The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`.
335a1cb98faSBarry Smith 
3361cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalEnd()`
337fa534816SMatthew G. Knepley @*/
338d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGlobalToNaturalBegin(DM dm, Vec gv, Vec nv)
339d71ae5a4SJacob Faibussowitsch {
340fa534816SMatthew G. Knepley   const PetscScalar *inarray;
341fa534816SMatthew G. Knepley   PetscScalar       *outarray;
342f6c7c579SMatthew G. Knepley   MPI_Comm           comm;
343a6a55facSMatthew G. Knepley   PetscMPIInt        size;
344fa534816SMatthew G. Knepley 
345fa534816SMatthew G. Knepley   PetscFunctionBegin;
3469566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_GlobalToNaturalBegin, dm, 0, 0, 0));
347f6c7c579SMatthew G. Knepley   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
348f6c7c579SMatthew G. Knepley   PetscCallMPI(MPI_Comm_size(comm, &size));
349fa534816SMatthew G. Knepley   if (dm->sfNatural) {
350f6c7c579SMatthew G. Knepley     if (PetscDefined(USE_DEBUG)) {
351f6c7c579SMatthew G. Knepley       PetscSection gs;
352f6c7c579SMatthew G. Knepley       PetscInt     Nl, n;
353f6c7c579SMatthew G. Knepley 
354f6c7c579SMatthew G. Knepley       PetscCall(PetscSFGetGraph(dm->sfNatural, NULL, &Nl, NULL, NULL));
355f6c7c579SMatthew G. Knepley       PetscCall(VecGetLocalSize(nv, &n));
356f6c7c579SMatthew 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);
357f6c7c579SMatthew G. Knepley 
358f6c7c579SMatthew G. Knepley       PetscCall(DMGetGlobalSection(dm, &gs));
359f6c7c579SMatthew G. Knepley       PetscCall(PetscSectionGetConstrainedStorageSize(gs, &Nl));
360f6c7c579SMatthew G. Knepley       PetscCall(VecGetLocalSize(gv, &n));
361f6c7c579SMatthew 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);
362f6c7c579SMatthew G. Knepley     }
3639566063dSJacob Faibussowitsch     PetscCall(VecGetArray(nv, &outarray));
3649566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(gv, &inarray));
3659566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_REPLACE));
3669566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(gv, &inarray));
3679566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(nv, &outarray));
368a6a55facSMatthew G. Knepley   } else if (size == 1) {
3699566063dSJacob Faibussowitsch     PetscCall(VecCopy(gv, nv));
370f7d195e4SLawrence Mitchell   } else {
37100045ab3SPierre Jolivet     PetscCheck(!dm->useNatural, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "DM global to natural SF not present. If DMPlexDistribute() was called and a section was defined, report to petsc-maint@mcs.anl.gov.");
37200045ab3SPierre Jolivet     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created. You must call DMSetUseNatural() before DMPlexDistribute().");
373f7d195e4SLawrence Mitchell   }
3749566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_GlobalToNaturalBegin, dm, 0, 0, 0));
3753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
376fa534816SMatthew G. Knepley }
377fa534816SMatthew G. Knepley 
378fa534816SMatthew G. Knepley /*@
379a1cb98faSBarry Smith   DMPlexGlobalToNaturalEnd - Rearranges a global `Vec` in the natural order.
380fa534816SMatthew G. Knepley 
38120f4b53cSBarry Smith   Collective
382fa534816SMatthew G. Knepley 
383fa534816SMatthew G. Knepley   Input Parameters:
384a1cb98faSBarry Smith + dm - The distributed `DMPLEX`
385a1cb98faSBarry Smith - gv - The global `Vec`
386fa534816SMatthew G. Knepley 
38797bb3fdcSJose E. Roman   Output Parameter:
388a1cb98faSBarry Smith . nv - The natural `Vec`
389fa534816SMatthew G. Knepley 
390fa534816SMatthew G. Knepley   Level: intermediate
391fa534816SMatthew G. Knepley 
392a1cb98faSBarry Smith   Note:
393a1cb98faSBarry Smith   The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`.
394a1cb98faSBarry Smith 
3951cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()`
396fa534816SMatthew G. Knepley  @*/
397d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGlobalToNaturalEnd(DM dm, Vec gv, Vec nv)
398d71ae5a4SJacob Faibussowitsch {
399fa534816SMatthew G. Knepley   const PetscScalar *inarray;
400fa534816SMatthew G. Knepley   PetscScalar       *outarray;
401a6a55facSMatthew G. Knepley   PetscMPIInt        size;
402fa534816SMatthew G. Knepley 
403fa534816SMatthew G. Knepley   PetscFunctionBegin;
4049566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_GlobalToNaturalEnd, dm, 0, 0, 0));
4059566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
406fa534816SMatthew G. Knepley   if (dm->sfNatural) {
4079566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(gv, &inarray));
4089566063dSJacob Faibussowitsch     PetscCall(VecGetArray(nv, &outarray));
4099566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_REPLACE));
4109566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(gv, &inarray));
4119566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(nv, &outarray));
412a6a55facSMatthew G. Knepley   } else if (size == 1) {
413f7d195e4SLawrence Mitchell   } else {
41400045ab3SPierre Jolivet     PetscCheck(!dm->useNatural, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "DM global to natural SF not present. If DMPlexDistribute() was called and a section was defined, report to petsc-maint@mcs.anl.gov.");
41500045ab3SPierre Jolivet     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created. You must call DMSetUseNatural() before DMPlexDistribute().");
416f7d195e4SLawrence Mitchell   }
4179566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_GlobalToNaturalEnd, dm, 0, 0, 0));
4183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
419fa534816SMatthew G. Knepley }
420fa534816SMatthew G. Knepley 
421fa534816SMatthew G. Knepley /*@
422a1cb98faSBarry Smith   DMPlexNaturalToGlobalBegin - Rearranges a `Vec` in the natural order to the Global order.
423fa534816SMatthew G. Knepley 
42420f4b53cSBarry Smith   Collective
425fa534816SMatthew G. Knepley 
426fa534816SMatthew G. Knepley   Input Parameters:
427a1cb98faSBarry Smith + dm - The distributed `DMPLEX`
428a1cb98faSBarry Smith - nv - The natural `Vec`
429fa534816SMatthew G. Knepley 
4302fe279fdSBarry Smith   Output Parameter:
431a1cb98faSBarry Smith . gv - The global `Vec`
432fa534816SMatthew G. Knepley 
433fa534816SMatthew G. Knepley   Level: intermediate
434fa534816SMatthew G. Knepley 
435a1cb98faSBarry Smith   Note:
436a1cb98faSBarry Smith   The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`.
437a1cb98faSBarry Smith 
43842747ad1SJacob Faibussowitsch .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexGlobalToNaturalEnd()`
439fa534816SMatthew G. Knepley @*/
440d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexNaturalToGlobalBegin(DM dm, Vec nv, Vec gv)
441d71ae5a4SJacob Faibussowitsch {
442fa534816SMatthew G. Knepley   const PetscScalar *inarray;
443fa534816SMatthew G. Knepley   PetscScalar       *outarray;
444a6a55facSMatthew G. Knepley   PetscMPIInt        size;
445fa534816SMatthew G. Knepley 
446fa534816SMatthew G. Knepley   PetscFunctionBegin;
4479566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_NaturalToGlobalBegin, dm, 0, 0, 0));
4489566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
449fa534816SMatthew G. Knepley   if (dm->sfNatural) {
450a5b23f4aSJose E. Roman     /* We only have access to the SF that goes from Global to Natural.
451b64f75a9SBlaise Bourdin        Instead of inverting dm->sfNatural, we can call PetscSFReduceBegin/End with MPI_Op MPI_SUM.
452b64f75a9SBlaise Bourdin        Here the SUM really does nothing since sfNatural is one to one, as long as gV is set to zero first. */
4539566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(gv));
4549566063dSJacob Faibussowitsch     PetscCall(VecGetArray(gv, &outarray));
4559566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(nv, &inarray));
4569566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_SUM));
4579566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(nv, &inarray));
4589566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(gv, &outarray));
459a6a55facSMatthew G. Knepley   } else if (size == 1) {
4609566063dSJacob Faibussowitsch     PetscCall(VecCopy(nv, gv));
461f7d195e4SLawrence Mitchell   } else {
46200045ab3SPierre Jolivet     PetscCheck(!dm->useNatural, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "DM global to natural SF not present. If DMPlexDistribute() was called and a section was defined, report to petsc-maint@mcs.anl.gov.");
46300045ab3SPierre Jolivet     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created. You must call DMSetUseNatural() before DMPlexDistribute().");
464f7d195e4SLawrence Mitchell   }
4659566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_NaturalToGlobalBegin, dm, 0, 0, 0));
4663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
467fa534816SMatthew G. Knepley }
468fa534816SMatthew G. Knepley 
469fa534816SMatthew G. Knepley /*@
470a1cb98faSBarry Smith   DMPlexNaturalToGlobalEnd - Rearranges a `Vec` in the natural order to the Global order.
471fa534816SMatthew G. Knepley 
47220f4b53cSBarry Smith   Collective
473fa534816SMatthew G. Knepley 
474fa534816SMatthew G. Knepley   Input Parameters:
475a1cb98faSBarry Smith + dm - The distributed `DMPLEX`
476a1cb98faSBarry Smith - nv - The natural `Vec`
477fa534816SMatthew G. Knepley 
4782fe279fdSBarry Smith   Output Parameter:
479a1cb98faSBarry Smith . gv - The global `Vec`
480fa534816SMatthew G. Knepley 
481fa534816SMatthew G. Knepley   Level: intermediate
482fa534816SMatthew G. Knepley 
483a1cb98faSBarry Smith   Note:
484a1cb98faSBarry Smith   The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`.
485a1cb98faSBarry Smith 
4861cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()`
487fa534816SMatthew G. Knepley  @*/
488d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexNaturalToGlobalEnd(DM dm, Vec nv, Vec gv)
489d71ae5a4SJacob Faibussowitsch {
490fa534816SMatthew G. Knepley   const PetscScalar *inarray;
491fa534816SMatthew G. Knepley   PetscScalar       *outarray;
492a6a55facSMatthew G. Knepley   PetscMPIInt        size;
493fa534816SMatthew G. Knepley 
494fa534816SMatthew G. Knepley   PetscFunctionBegin;
4959566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_NaturalToGlobalEnd, dm, 0, 0, 0));
4969566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
497fa534816SMatthew G. Knepley   if (dm->sfNatural) {
4989566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(nv, &inarray));
4999566063dSJacob Faibussowitsch     PetscCall(VecGetArray(gv, &outarray));
5009566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_SUM));
5019566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(nv, &inarray));
5029566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(gv, &outarray));
503a6a55facSMatthew G. Knepley   } else if (size == 1) {
504f7d195e4SLawrence Mitchell   } else {
50500045ab3SPierre Jolivet     PetscCheck(!dm->useNatural, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "DM global to natural SF not present. If DMPlexDistribute() was called and a section was defined, report to petsc-maint@mcs.anl.gov.");
50600045ab3SPierre Jolivet     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created. You must call DMSetUseNatural() before DMPlexDistribute().");
507f7d195e4SLawrence Mitchell   }
5089566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_NaturalToGlobalEnd, dm, 0, 0, 0));
5093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
510fa534816SMatthew G. Knepley }
51109cf685aSAlexis Marboeuf 
51209cf685aSAlexis Marboeuf /*@
5135cb80ecdSBarry Smith   DMPlexCreateNaturalVector - Provide a `Vec` capable of holding the natural ordering and distribution.
51409cf685aSAlexis Marboeuf 
51520f4b53cSBarry Smith   Collective
51609cf685aSAlexis Marboeuf 
5175cb80ecdSBarry Smith   Input Parameter:
518a1cb98faSBarry Smith . dm - The distributed `DMPLEX`
51909cf685aSAlexis Marboeuf 
5205cb80ecdSBarry Smith   Output Parameter:
5215cb80ecdSBarry Smith . nv - The natural `Vec`
52209cf685aSAlexis Marboeuf 
52309cf685aSAlexis Marboeuf   Level: intermediate
52409cf685aSAlexis Marboeuf 
525a1cb98faSBarry Smith   Note:
526a1cb98faSBarry Smith   The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`.
527a1cb98faSBarry Smith 
5281cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()`
52909cf685aSAlexis Marboeuf  @*/
530d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreateNaturalVector(DM dm, Vec *nv)
531d71ae5a4SJacob Faibussowitsch {
53209cf685aSAlexis Marboeuf   PetscMPIInt size;
53309cf685aSAlexis Marboeuf 
53409cf685aSAlexis Marboeuf   PetscFunctionBegin;
53509cf685aSAlexis Marboeuf   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
53609cf685aSAlexis Marboeuf   if (dm->sfNatural) {
5376fba1fe0SBlaise Bourdin     PetscInt nleaves, bs, maxbs;
5382e8d78feSBlaise Bourdin     Vec      v;
5396fba1fe0SBlaise Bourdin 
5406fba1fe0SBlaise Bourdin     /*
5416fba1fe0SBlaise Bourdin       Setting the natural vector block size.
5426fba1fe0SBlaise Bourdin       We can't get it from a global vector because of constraints, and the block size in the local vector
5436fba1fe0SBlaise Bourdin       may be inconsistent across processes, typically when some local vectors have size 0, their block size is set to 1
5446fba1fe0SBlaise Bourdin     */
5452e8d78feSBlaise Bourdin     PetscCall(DMGetLocalVector(dm, &v));
5462e8d78feSBlaise Bourdin     PetscCall(VecGetBlockSize(v, &bs));
547*462c564dSBarry Smith     PetscCallMPI(MPIU_Allreduce(&bs, &maxbs, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm)));
5486fba1fe0SBlaise Bourdin     if (bs == 1 && maxbs > 1) bs = maxbs;
5492e8d78feSBlaise Bourdin     PetscCall(DMRestoreLocalVector(dm, &v));
5502e8d78feSBlaise Bourdin 
55109cf685aSAlexis Marboeuf     PetscCall(PetscSFGetGraph(dm->sfNatural, NULL, &nleaves, NULL, NULL));
55209cf685aSAlexis Marboeuf     PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), nv));
55309cf685aSAlexis Marboeuf     PetscCall(VecSetSizes(*nv, nleaves, PETSC_DETERMINE));
55409cf685aSAlexis Marboeuf     PetscCall(VecSetBlockSize(*nv, bs));
55509cf685aSAlexis Marboeuf     PetscCall(VecSetType(*nv, dm->vectype));
55609cf685aSAlexis Marboeuf     PetscCall(VecSetDM(*nv, dm));
55709cf685aSAlexis Marboeuf   } else if (size == 1) {
5582e8d78feSBlaise Bourdin     PetscCall(DMCreateLocalVector(dm, nv));
55909cf685aSAlexis Marboeuf   } else {
56000045ab3SPierre Jolivet     PetscCheck(!dm->useNatural, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "DM global to natural SF not present. If DMPlexDistribute() was called and a section was defined, report to petsc-maint@mcs.anl.gov.");
56100045ab3SPierre Jolivet     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created. You must call DMSetUseNatural() before DMPlexDistribute().");
56209cf685aSAlexis Marboeuf   }
5633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
56409cf685aSAlexis Marboeuf }
565