xref: /petsc/src/dm/impls/plex/plexnatural.c (revision a2aa6d814414e5e930288a4c53e1da56c681e2e0)
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));
146462c564dSBarry 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;
237*a2aa6d81SJames Wright   PetscSection oldGlobalSection, newGlobalSection;
23816635c05SJames Wright   PetscInt    *remoteOffsets;
23916635c05SJames Wright   PetscBool    debug = PETSC_FALSE;
24016635c05SJames Wright 
24116635c05SJames Wright   PetscFunctionBegin;
24216635c05SJames Wright   PetscCall(PetscObjectGetComm((PetscObject)dmNew, &comm));
24316635c05SJames Wright   if (!sfMigration) {
24416635c05SJames Wright     /* If sfMigration is missing, sfNatural cannot be computed and is set to NULL */
24516635c05SJames Wright     *sfNaturalNew = NULL;
24616635c05SJames Wright     PetscFunctionReturn(PETSC_SUCCESS);
24716635c05SJames Wright   }
24816635c05SJames Wright   if (debug) PetscCall(PetscSFView(sfMigration, NULL));
24916635c05SJames Wright 
25016635c05SJames Wright   { // Create oldGlobalSection and newGlobalSection *with* localOffsets
25116635c05SJames Wright     PetscSection oldLocalSection, newLocalSection;
25216635c05SJames Wright     PetscSF      pointSF;
25316635c05SJames Wright 
25416635c05SJames Wright     PetscCall(DMGetLocalSection(dmOld, &oldLocalSection));
25516635c05SJames Wright     PetscCall(DMGetPointSF(dmOld, &pointSF));
25616635c05SJames Wright     PetscCall(PetscSectionCreateGlobalSection(oldLocalSection, pointSF, PETSC_TRUE, PETSC_TRUE, PETSC_TRUE, &oldGlobalSection));
25716635c05SJames Wright 
25816635c05SJames Wright     PetscCall(PetscSectionCreate(comm, &newLocalSection));
25916635c05SJames Wright     PetscCall(PetscSFDistributeSection(sfMigration, oldLocalSection, NULL, newLocalSection));
26016635c05SJames Wright     PetscCall(DMSetLocalSection(dmNew, newLocalSection));
26116635c05SJames Wright 
26216635c05SJames Wright     PetscCall(PetscSectionCreate(comm, &newGlobalSection));
26316635c05SJames Wright     PetscCall(DMGetPointSF(dmNew, &pointSF));
26416635c05SJames Wright     PetscCall(PetscSectionCreateGlobalSection(newLocalSection, pointSF, PETSC_TRUE, PETSC_TRUE, PETSC_TRUE, &newGlobalSection));
26516635c05SJames Wright 
26616635c05SJames Wright     PetscCall(PetscObjectSetName((PetscObject)oldLocalSection, "Old Local Section"));
26716635c05SJames Wright     if (debug) PetscCall(PetscSectionView(oldLocalSection, NULL));
26816635c05SJames Wright     PetscCall(PetscObjectSetName((PetscObject)oldGlobalSection, "Old Global Section"));
26916635c05SJames Wright     if (debug) PetscCall(PetscSectionView(oldGlobalSection, NULL));
27016635c05SJames Wright     PetscCall(PetscObjectSetName((PetscObject)newLocalSection, "New Local Section"));
27116635c05SJames Wright     if (debug) PetscCall(PetscSectionView(newLocalSection, NULL));
27216635c05SJames Wright     PetscCall(PetscObjectSetName((PetscObject)newGlobalSection, "New Global Section"));
27316635c05SJames Wright     if (debug) PetscCall(PetscSectionView(newGlobalSection, NULL));
27416635c05SJames Wright     PetscCall(PetscSectionDestroy(&newLocalSection));
27516635c05SJames Wright   }
27616635c05SJames Wright 
27716635c05SJames Wright   { // Create remoteOffsets array, mapping the oldGlobalSection offsets to the local points (according to sfMigration)
27816635c05SJames Wright     PetscInt lpStart, lpEnd, rpStart, rpEnd;
27916635c05SJames Wright 
28016635c05SJames Wright     PetscCall(PetscSectionGetChart(oldGlobalSection, &rpStart, &rpEnd));
28116635c05SJames Wright     PetscCall(PetscSectionGetChart(newGlobalSection, &lpStart, &lpEnd));
28216635c05SJames Wright 
28316635c05SJames Wright     // in `PetscSFDistributeSection` (where this is taken from), it possibly makes a new embedded SF. Should possibly do that here?
28416635c05SJames Wright     PetscCall(PetscMalloc1(lpEnd - lpStart, &remoteOffsets));
28516635c05SJames Wright     PetscCall(PetscSFBcastBegin(sfMigration, MPIU_INT, PetscSafePointerPlusOffset(oldGlobalSection->atlasOff, -rpStart), PetscSafePointerPlusOffset(remoteOffsets, -lpStart), MPI_REPLACE));
28616635c05SJames Wright     PetscCall(PetscSFBcastEnd(sfMigration, MPIU_INT, PetscSafePointerPlusOffset(oldGlobalSection->atlasOff, -rpStart), PetscSafePointerPlusOffset(remoteOffsets, -lpStart), MPI_REPLACE));
28716635c05SJames Wright     if (debug) {
28816635c05SJames Wright       PetscViewer viewer;
28916635c05SJames Wright 
29016635c05SJames Wright       PetscCall(PetscPrintf(comm, "Remote Offsets:\n"));
29116635c05SJames Wright       PetscCall(PetscViewerASCIIGetStdout(comm, &viewer));
29216635c05SJames Wright       PetscCall(PetscIntView(lpEnd - lpStart, remoteOffsets, viewer));
29316635c05SJames Wright     }
29416635c05SJames Wright   }
29516635c05SJames Wright 
29616635c05SJames Wright   { // Create SF from oldGlobalSection to newGlobalSection and compose with sfNaturalOld
29716635c05SJames Wright     PetscSF oldglob_to_newglob_sf, newglob_to_oldglob_sf;
29816635c05SJames Wright 
29916635c05SJames Wright     PetscCall(PetscSFCreateSectionSF(sfMigration, oldGlobalSection, remoteOffsets, newGlobalSection, &oldglob_to_newglob_sf));
30016635c05SJames Wright     PetscCall(PetscObjectSetName((PetscObject)oldglob_to_newglob_sf, "OldGlobal-to-NewGlobal SF"));
30116635c05SJames Wright     if (debug) PetscCall(PetscSFView(oldglob_to_newglob_sf, NULL));
30216635c05SJames Wright 
30316635c05SJames Wright     PetscCall(PetscSFCreateInverseSF(oldglob_to_newglob_sf, &newglob_to_oldglob_sf));
30416635c05SJames Wright     PetscCall(PetscObjectSetName((PetscObject)newglob_to_oldglob_sf, "NewGlobal-to-OldGlobal SF"));
30516635c05SJames Wright     PetscCall(PetscObjectViewFromOptions((PetscObject)newglob_to_oldglob_sf, (PetscObject)dmOld, "-natural_migrate_sf_view"));
30616635c05SJames Wright     PetscCall(PetscSFCompose(newglob_to_oldglob_sf, sfNaturalOld, sfNaturalNew));
30716635c05SJames Wright 
30816635c05SJames Wright     PetscCall(PetscSFDestroy(&oldglob_to_newglob_sf));
30916635c05SJames Wright     PetscCall(PetscSFDestroy(&newglob_to_oldglob_sf));
31016635c05SJames Wright   }
31116635c05SJames Wright 
31216635c05SJames Wright   PetscCall(PetscSectionDestroy(&oldGlobalSection));
31316635c05SJames Wright   PetscCall(PetscSectionDestroy(&newGlobalSection));
31416635c05SJames Wright   PetscCall(PetscFree(remoteOffsets));
31516635c05SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
31616635c05SJames Wright }
31716635c05SJames Wright 
31816635c05SJames Wright /*@
319a1cb98faSBarry Smith   DMPlexGlobalToNaturalBegin - Rearranges a global `Vec` in the natural order.
320fa534816SMatthew G. Knepley 
32120f4b53cSBarry Smith   Collective
322fa534816SMatthew G. Knepley 
323fa534816SMatthew G. Knepley   Input Parameters:
324a1cb98faSBarry Smith + dm - The distributed `DMPLEX`
325a1cb98faSBarry Smith - gv - The global `Vec`
326fa534816SMatthew G. Knepley 
3272fe279fdSBarry Smith   Output Parameter:
32820f4b53cSBarry Smith . nv - `Vec` in the canonical ordering distributed over all processors associated with `gv`
329fa534816SMatthew G. Knepley 
330fa534816SMatthew G. Knepley   Level: intermediate
331fa534816SMatthew G. Knepley 
332a1cb98faSBarry Smith   Note:
333a1cb98faSBarry Smith   The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`.
334a1cb98faSBarry Smith 
3351cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalEnd()`
336fa534816SMatthew G. Knepley @*/
337d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGlobalToNaturalBegin(DM dm, Vec gv, Vec nv)
338d71ae5a4SJacob Faibussowitsch {
339fa534816SMatthew G. Knepley   const PetscScalar *inarray;
340fa534816SMatthew G. Knepley   PetscScalar       *outarray;
341f6c7c579SMatthew G. Knepley   MPI_Comm           comm;
342a6a55facSMatthew G. Knepley   PetscMPIInt        size;
343fa534816SMatthew G. Knepley 
344fa534816SMatthew G. Knepley   PetscFunctionBegin;
3459566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_GlobalToNaturalBegin, dm, 0, 0, 0));
346f6c7c579SMatthew G. Knepley   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
347f6c7c579SMatthew G. Knepley   PetscCallMPI(MPI_Comm_size(comm, &size));
348fa534816SMatthew G. Knepley   if (dm->sfNatural) {
349f6c7c579SMatthew G. Knepley     if (PetscDefined(USE_DEBUG)) {
350f6c7c579SMatthew G. Knepley       PetscSection gs;
351f6c7c579SMatthew G. Knepley       PetscInt     Nl, n;
352f6c7c579SMatthew G. Knepley 
353f6c7c579SMatthew G. Knepley       PetscCall(PetscSFGetGraph(dm->sfNatural, NULL, &Nl, NULL, NULL));
354f6c7c579SMatthew G. Knepley       PetscCall(VecGetLocalSize(nv, &n));
355f6c7c579SMatthew 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);
356f6c7c579SMatthew G. Knepley 
357f6c7c579SMatthew G. Knepley       PetscCall(DMGetGlobalSection(dm, &gs));
358f6c7c579SMatthew G. Knepley       PetscCall(PetscSectionGetConstrainedStorageSize(gs, &Nl));
359f6c7c579SMatthew G. Knepley       PetscCall(VecGetLocalSize(gv, &n));
360f6c7c579SMatthew 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);
361f6c7c579SMatthew G. Knepley     }
3629566063dSJacob Faibussowitsch     PetscCall(VecGetArray(nv, &outarray));
3639566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(gv, &inarray));
3649566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_REPLACE));
3659566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(gv, &inarray));
3669566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(nv, &outarray));
367a6a55facSMatthew G. Knepley   } else if (size == 1) {
3689566063dSJacob Faibussowitsch     PetscCall(VecCopy(gv, nv));
369f7d195e4SLawrence Mitchell   } else {
37000045ab3SPierre 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.");
37100045ab3SPierre Jolivet     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created. You must call DMSetUseNatural() before DMPlexDistribute().");
372f7d195e4SLawrence Mitchell   }
3739566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_GlobalToNaturalBegin, dm, 0, 0, 0));
3743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
375fa534816SMatthew G. Knepley }
376fa534816SMatthew G. Knepley 
377fa534816SMatthew G. Knepley /*@
378a1cb98faSBarry Smith   DMPlexGlobalToNaturalEnd - Rearranges a global `Vec` in the natural order.
379fa534816SMatthew G. Knepley 
38020f4b53cSBarry Smith   Collective
381fa534816SMatthew G. Knepley 
382fa534816SMatthew G. Knepley   Input Parameters:
383a1cb98faSBarry Smith + dm - The distributed `DMPLEX`
384a1cb98faSBarry Smith - gv - The global `Vec`
385fa534816SMatthew G. Knepley 
38697bb3fdcSJose E. Roman   Output Parameter:
387a1cb98faSBarry Smith . nv - The natural `Vec`
388fa534816SMatthew G. Knepley 
389fa534816SMatthew G. Knepley   Level: intermediate
390fa534816SMatthew G. Knepley 
391a1cb98faSBarry Smith   Note:
392a1cb98faSBarry Smith   The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`.
393a1cb98faSBarry Smith 
3941cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()`
395fa534816SMatthew G. Knepley  @*/
396d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGlobalToNaturalEnd(DM dm, Vec gv, Vec nv)
397d71ae5a4SJacob Faibussowitsch {
398fa534816SMatthew G. Knepley   const PetscScalar *inarray;
399fa534816SMatthew G. Knepley   PetscScalar       *outarray;
400a6a55facSMatthew G. Knepley   PetscMPIInt        size;
401fa534816SMatthew G. Knepley 
402fa534816SMatthew G. Knepley   PetscFunctionBegin;
4039566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_GlobalToNaturalEnd, dm, 0, 0, 0));
4049566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
405fa534816SMatthew G. Knepley   if (dm->sfNatural) {
4069566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(gv, &inarray));
4079566063dSJacob Faibussowitsch     PetscCall(VecGetArray(nv, &outarray));
4089566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_REPLACE));
4099566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(gv, &inarray));
4109566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(nv, &outarray));
411a6a55facSMatthew G. Knepley   } else if (size == 1) {
412f7d195e4SLawrence Mitchell   } else {
41300045ab3SPierre 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.");
41400045ab3SPierre Jolivet     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created. You must call DMSetUseNatural() before DMPlexDistribute().");
415f7d195e4SLawrence Mitchell   }
4169566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_GlobalToNaturalEnd, dm, 0, 0, 0));
4173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
418fa534816SMatthew G. Knepley }
419fa534816SMatthew G. Knepley 
420fa534816SMatthew G. Knepley /*@
421a1cb98faSBarry Smith   DMPlexNaturalToGlobalBegin - Rearranges a `Vec` in the natural order to the Global order.
422fa534816SMatthew G. Knepley 
42320f4b53cSBarry Smith   Collective
424fa534816SMatthew G. Knepley 
425fa534816SMatthew G. Knepley   Input Parameters:
426a1cb98faSBarry Smith + dm - The distributed `DMPLEX`
427a1cb98faSBarry Smith - nv - The natural `Vec`
428fa534816SMatthew G. Knepley 
4292fe279fdSBarry Smith   Output Parameter:
430a1cb98faSBarry Smith . gv - The global `Vec`
431fa534816SMatthew G. Knepley 
432fa534816SMatthew G. Knepley   Level: intermediate
433fa534816SMatthew G. Knepley 
434a1cb98faSBarry Smith   Note:
435a1cb98faSBarry Smith   The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`.
436a1cb98faSBarry Smith 
43742747ad1SJacob Faibussowitsch .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexGlobalToNaturalEnd()`
438fa534816SMatthew G. Knepley @*/
439d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexNaturalToGlobalBegin(DM dm, Vec nv, Vec gv)
440d71ae5a4SJacob Faibussowitsch {
441fa534816SMatthew G. Knepley   const PetscScalar *inarray;
442fa534816SMatthew G. Knepley   PetscScalar       *outarray;
443a6a55facSMatthew G. Knepley   PetscMPIInt        size;
444fa534816SMatthew G. Knepley 
445fa534816SMatthew G. Knepley   PetscFunctionBegin;
4469566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_NaturalToGlobalBegin, dm, 0, 0, 0));
4479566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
448fa534816SMatthew G. Knepley   if (dm->sfNatural) {
449a5b23f4aSJose E. Roman     /* We only have access to the SF that goes from Global to Natural.
450b64f75a9SBlaise Bourdin        Instead of inverting dm->sfNatural, we can call PetscSFReduceBegin/End with MPI_Op MPI_SUM.
451b64f75a9SBlaise Bourdin        Here the SUM really does nothing since sfNatural is one to one, as long as gV is set to zero first. */
4529566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(gv));
4539566063dSJacob Faibussowitsch     PetscCall(VecGetArray(gv, &outarray));
4549566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(nv, &inarray));
4559566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_SUM));
4569566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(nv, &inarray));
4579566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(gv, &outarray));
458a6a55facSMatthew G. Knepley   } else if (size == 1) {
4599566063dSJacob Faibussowitsch     PetscCall(VecCopy(nv, gv));
460f7d195e4SLawrence Mitchell   } else {
46100045ab3SPierre 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.");
46200045ab3SPierre Jolivet     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created. You must call DMSetUseNatural() before DMPlexDistribute().");
463f7d195e4SLawrence Mitchell   }
4649566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_NaturalToGlobalBegin, dm, 0, 0, 0));
4653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
466fa534816SMatthew G. Knepley }
467fa534816SMatthew G. Knepley 
468fa534816SMatthew G. Knepley /*@
469a1cb98faSBarry Smith   DMPlexNaturalToGlobalEnd - Rearranges a `Vec` in the natural order to the Global order.
470fa534816SMatthew G. Knepley 
47120f4b53cSBarry Smith   Collective
472fa534816SMatthew G. Knepley 
473fa534816SMatthew G. Knepley   Input Parameters:
474a1cb98faSBarry Smith + dm - The distributed `DMPLEX`
475a1cb98faSBarry Smith - nv - The natural `Vec`
476fa534816SMatthew G. Knepley 
4772fe279fdSBarry Smith   Output Parameter:
478a1cb98faSBarry Smith . gv - The global `Vec`
479fa534816SMatthew G. Knepley 
480fa534816SMatthew G. Knepley   Level: intermediate
481fa534816SMatthew G. Knepley 
482a1cb98faSBarry Smith   Note:
483a1cb98faSBarry Smith   The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`.
484a1cb98faSBarry Smith 
4851cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()`
486fa534816SMatthew G. Knepley  @*/
487d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexNaturalToGlobalEnd(DM dm, Vec nv, Vec gv)
488d71ae5a4SJacob Faibussowitsch {
489fa534816SMatthew G. Knepley   const PetscScalar *inarray;
490fa534816SMatthew G. Knepley   PetscScalar       *outarray;
491a6a55facSMatthew G. Knepley   PetscMPIInt        size;
492fa534816SMatthew G. Knepley 
493fa534816SMatthew G. Knepley   PetscFunctionBegin;
4949566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_NaturalToGlobalEnd, dm, 0, 0, 0));
4959566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
496fa534816SMatthew G. Knepley   if (dm->sfNatural) {
4979566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(nv, &inarray));
4989566063dSJacob Faibussowitsch     PetscCall(VecGetArray(gv, &outarray));
4999566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_SUM));
5009566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(nv, &inarray));
5019566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(gv, &outarray));
502a6a55facSMatthew G. Knepley   } else if (size == 1) {
503f7d195e4SLawrence Mitchell   } else {
50400045ab3SPierre 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.");
50500045ab3SPierre Jolivet     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created. You must call DMSetUseNatural() before DMPlexDistribute().");
506f7d195e4SLawrence Mitchell   }
5079566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_NaturalToGlobalEnd, dm, 0, 0, 0));
5083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
509fa534816SMatthew G. Knepley }
51009cf685aSAlexis Marboeuf 
51109cf685aSAlexis Marboeuf /*@
5125cb80ecdSBarry Smith   DMPlexCreateNaturalVector - Provide a `Vec` capable of holding the natural ordering and distribution.
51309cf685aSAlexis Marboeuf 
51420f4b53cSBarry Smith   Collective
51509cf685aSAlexis Marboeuf 
5165cb80ecdSBarry Smith   Input Parameter:
517a1cb98faSBarry Smith . dm - The distributed `DMPLEX`
51809cf685aSAlexis Marboeuf 
5195cb80ecdSBarry Smith   Output Parameter:
5205cb80ecdSBarry Smith . nv - The natural `Vec`
52109cf685aSAlexis Marboeuf 
52209cf685aSAlexis Marboeuf   Level: intermediate
52309cf685aSAlexis Marboeuf 
524a1cb98faSBarry Smith   Note:
525a1cb98faSBarry Smith   The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`.
526a1cb98faSBarry Smith 
5271cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()`
52809cf685aSAlexis Marboeuf  @*/
529d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreateNaturalVector(DM dm, Vec *nv)
530d71ae5a4SJacob Faibussowitsch {
53109cf685aSAlexis Marboeuf   PetscMPIInt size;
53209cf685aSAlexis Marboeuf 
53309cf685aSAlexis Marboeuf   PetscFunctionBegin;
53409cf685aSAlexis Marboeuf   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
53509cf685aSAlexis Marboeuf   if (dm->sfNatural) {
5366fba1fe0SBlaise Bourdin     PetscInt nleaves, bs, maxbs;
5372e8d78feSBlaise Bourdin     Vec      v;
5386fba1fe0SBlaise Bourdin 
5396fba1fe0SBlaise Bourdin     /*
5406fba1fe0SBlaise Bourdin       Setting the natural vector block size.
5416fba1fe0SBlaise Bourdin       We can't get it from a global vector because of constraints, and the block size in the local vector
5426fba1fe0SBlaise Bourdin       may be inconsistent across processes, typically when some local vectors have size 0, their block size is set to 1
5436fba1fe0SBlaise Bourdin     */
5442e8d78feSBlaise Bourdin     PetscCall(DMGetLocalVector(dm, &v));
5452e8d78feSBlaise Bourdin     PetscCall(VecGetBlockSize(v, &bs));
546462c564dSBarry Smith     PetscCallMPI(MPIU_Allreduce(&bs, &maxbs, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm)));
5476fba1fe0SBlaise Bourdin     if (bs == 1 && maxbs > 1) bs = maxbs;
5482e8d78feSBlaise Bourdin     PetscCall(DMRestoreLocalVector(dm, &v));
5492e8d78feSBlaise Bourdin 
55009cf685aSAlexis Marboeuf     PetscCall(PetscSFGetGraph(dm->sfNatural, NULL, &nleaves, NULL, NULL));
55109cf685aSAlexis Marboeuf     PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), nv));
55209cf685aSAlexis Marboeuf     PetscCall(VecSetSizes(*nv, nleaves, PETSC_DETERMINE));
55309cf685aSAlexis Marboeuf     PetscCall(VecSetBlockSize(*nv, bs));
55409cf685aSAlexis Marboeuf     PetscCall(VecSetType(*nv, dm->vectype));
55509cf685aSAlexis Marboeuf     PetscCall(VecSetDM(*nv, dm));
55609cf685aSAlexis Marboeuf   } else if (size == 1) {
5572e8d78feSBlaise Bourdin     PetscCall(DMCreateLocalVector(dm, nv));
55809cf685aSAlexis Marboeuf   } else {
55900045ab3SPierre 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.");
56000045ab3SPierre Jolivet     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created. You must call DMSetUseNatural() before DMPlexDistribute().");
56109cf685aSAlexis Marboeuf   }
5623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
56309cf685aSAlexis Marboeuf }
564