xref: /petsc/src/dm/impls/plex/plexnatural.c (revision b5f0bcd6e9e8ed97648738542f5163d94f7b1782)
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 @*/
DMPlexSetMigrationSF(DM dm,PetscSF migrationSF)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 
33*c039638dSMatthew G. Knepley   Not 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 @*/
DMPlexGetMigrationSF(DM dm,PetscSF * migrationSF)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   DMPlexCreateGlobalToNaturalSF - Creates the `PetscSF` for mapping Global `Vec` to the Natural `Vec`
54fa534816SMatthew G. Knepley 
55fa534816SMatthew G. Knepley   Input Parameters:
56a1cb98faSBarry Smith + dm          - The redistributed `DM`
5720f4b53cSBarry Smith . section     - The local `PetscSection` describing the `Vec` before the mesh was distributed, or `NULL` if not available
5820f4b53cSBarry Smith - sfMigration - The `PetscSF` used to distribute the mesh, or `NULL` if it cannot be computed
59fa534816SMatthew G. Knepley 
605d3b26e6SMatthew G. Knepley   Output Parameter:
61a1cb98faSBarry Smith . sfNatural - `PetscSF` for mapping the `Vec` in PETSc ordering to the canonical ordering
625d3b26e6SMatthew G. Knepley 
63fa534816SMatthew G. Knepley   Level: intermediate
64fa534816SMatthew G. Knepley 
65a1cb98faSBarry Smith   Note:
66a1cb98faSBarry Smith   This is not typically called by the user.
67a1cb98faSBarry Smith 
681cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `PetscSF`, `PetscSection`, `DMPlexDistribute()`, `DMPlexDistributeField()`
69fa534816SMatthew G. Knepley  @*/
DMPlexCreateGlobalToNaturalSF(DM dm,PetscSection section,PetscSF sfMigration,PetscSF * sfNatural)70d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreateGlobalToNaturalSF(DM dm, PetscSection section, PetscSF sfMigration, PetscSF *sfNatural)
71d71ae5a4SJacob Faibussowitsch {
72fa534816SMatthew G. Knepley   MPI_Comm     comm;
73e7858701SMatthew G. Knepley   PetscSF      sf, sfEmbed, sfField;
74e5b44f4fSMatthew G. Knepley   PetscSection gSection, sectionDist, gLocSection;
75fa534816SMatthew G. Knepley   PetscInt    *spoints, *remoteOffsets;
7657a2b1f8SBlaise Bourdin   PetscInt     ssize, pStart, pEnd, p, localSize, maxStorageSize;
77e7858701SMatthew G. Knepley   PetscBool    destroyFlag = PETSC_FALSE, debug = PETSC_FALSE;
78fa534816SMatthew G. Knepley 
79fa534816SMatthew G. Knepley   PetscFunctionBegin;
809566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
815c3f5608SAlexisMarb   if (!sfMigration) {
82e7858701SMatthew G. Knepley     /* If sfMigration is missing, sfNatural cannot be computed and is set to NULL */
835c3f5608SAlexisMarb     *sfNatural = NULL;
843ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
855c3f5608SAlexisMarb   } else if (!section) {
86e7858701SMatthew G. Knepley     /* If the sequential section is not provided (NULL), it is reconstructed from the parallel section */
875c3f5608SAlexisMarb     PetscSF      sfMigrationInv;
885c3f5608SAlexisMarb     PetscSection localSection;
895c3f5608SAlexisMarb 
909566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dm, &localSection));
919566063dSJacob Faibussowitsch     PetscCall(PetscSFCreateInverseSF(sfMigration, &sfMigrationInv));
929566063dSJacob Faibussowitsch     PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), &section));
939566063dSJacob Faibussowitsch     PetscCall(PetscSFDistributeSection(sfMigrationInv, localSection, NULL, section));
949566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&sfMigrationInv));
955c3f5608SAlexisMarb     destroyFlag = PETSC_TRUE;
965c3f5608SAlexisMarb   }
97e7858701SMatthew G. Knepley   if (debug) PetscCall(PetscSFView(sfMigration, NULL));
98e5b44f4fSMatthew G. Knepley   /* Create a new section from distributing the original section */
999566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(comm, &sectionDist));
1009566063dSJacob Faibussowitsch   PetscCall(PetscSFDistributeSection(sfMigration, section, &remoteOffsets, sectionDist));
101e7858701SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)sectionDist, "Migrated Section"));
102e7858701SMatthew G. Knepley   if (debug) PetscCall(PetscSectionView(sectionDist, NULL));
1039566063dSJacob Faibussowitsch   PetscCall(DMSetLocalSection(dm, sectionDist));
104e7858701SMatthew G. Knepley   /* If a sequential section is provided but no dof is affected, sfNatural cannot be computed and is set to NULL */
10557a2b1f8SBlaise Bourdin   PetscCall(PetscSectionGetStorageSize(sectionDist, &localSize));
106462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(&localSize, &maxStorageSize, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm)));
10757a2b1f8SBlaise Bourdin   if (maxStorageSize) {
108e7858701SMatthew G. Knepley     const PetscInt *leaves;
109e7858701SMatthew G. Knepley     PetscInt       *sortleaves, *indices;
110e7858701SMatthew G. Knepley     PetscInt        Nl;
111e7858701SMatthew G. Knepley 
112fa534816SMatthew G. Knepley     /* Get a pruned version of migration SF */
1139566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalSection(dm, &gSection));
114e7858701SMatthew G. Knepley     if (debug) PetscCall(PetscSectionView(gSection, NULL));
115e7858701SMatthew G. Knepley     PetscCall(PetscSFGetGraph(sfMigration, NULL, &Nl, &leaves, NULL));
1169566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetChart(gSection, &pStart, &pEnd));
117fa534816SMatthew G. Knepley     for (p = pStart, ssize = 0; p < pEnd; ++p) {
118fa534816SMatthew G. Knepley       PetscInt dof, off;
119fa534816SMatthew G. Knepley 
1209566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(gSection, p, &dof));
1219566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(gSection, p, &off));
122fa534816SMatthew G. Knepley       if ((dof > 0) && (off >= 0)) ++ssize;
123fa534816SMatthew G. Knepley     }
124e7858701SMatthew G. Knepley     PetscCall(PetscMalloc3(ssize, &spoints, Nl, &sortleaves, Nl, &indices));
1259371c9d4SSatish Balay     for (p = 0; p < Nl; ++p) {
1269371c9d4SSatish Balay       sortleaves[p] = leaves ? leaves[p] : p;
1279371c9d4SSatish Balay       indices[p]    = p;
1289371c9d4SSatish Balay     }
129e7858701SMatthew G. Knepley     PetscCall(PetscSortIntWithArray(Nl, sortleaves, indices));
130fa534816SMatthew G. Knepley     for (p = pStart, ssize = 0; p < pEnd; ++p) {
131e7858701SMatthew G. Knepley       PetscInt dof, off, loc;
132fa534816SMatthew G. Knepley 
1339566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(gSection, p, &dof));
1349566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(gSection, p, &off));
135e7858701SMatthew G. Knepley       if ((dof > 0) && (off >= 0)) {
136e7858701SMatthew G. Knepley         PetscCall(PetscFindInt(p, Nl, sortleaves, &loc));
137e7858701SMatthew 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);
138e7858701SMatthew G. Knepley         spoints[ssize++] = indices[loc];
139e7858701SMatthew G. Knepley       }
140fa534816SMatthew G. Knepley     }
1419566063dSJacob Faibussowitsch     PetscCall(PetscSFCreateEmbeddedLeafSF(sfMigration, ssize, spoints, &sfEmbed));
142e7858701SMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)sfEmbed, "Embedded SF"));
143e7858701SMatthew G. Knepley     PetscCall(PetscFree3(spoints, sortleaves, indices));
144e7858701SMatthew G. Knepley     if (debug) PetscCall(PetscSFView(sfEmbed, NULL));
145e7858701SMatthew G. Knepley     /* Create the SF associated with this section
146e7858701SMatthew G. Knepley          Roots are natural dofs, leaves are global dofs */
1479566063dSJacob Faibussowitsch     PetscCall(DMGetPointSF(dm, &sf));
148eb9d3e4dSMatthew G. Knepley     PetscCall(PetscSectionCreateGlobalSection(sectionDist, sf, PETSC_TRUE, PETSC_TRUE, PETSC_TRUE, &gLocSection));
1499566063dSJacob Faibussowitsch     PetscCall(PetscSFCreateSectionSF(sfEmbed, section, remoteOffsets, gLocSection, &sfField));
1509566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&sfEmbed));
1519566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&gLocSection));
152e7858701SMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)sfField, "Natural-to-Global SF"));
153e7858701SMatthew G. Knepley     if (debug) PetscCall(PetscSFView(sfField, NULL));
154e7858701SMatthew G. Knepley     /* Invert the field SF
155e7858701SMatthew G. Knepley          Roots are global dofs, leaves are natural dofs */
156e7858701SMatthew G. Knepley     PetscCall(PetscSFCreateInverseSF(sfField, sfNatural));
157e7858701SMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)*sfNatural, "Global-to-Natural SF"));
1589566063dSJacob Faibussowitsch     PetscCall(PetscObjectViewFromOptions((PetscObject)*sfNatural, NULL, "-globaltonatural_sf_view"));
159fa534816SMatthew G. Knepley     /* Clean up */
160e7858701SMatthew G. Knepley     PetscCall(PetscSFDestroy(&sfField));
16114744f87SBlaise Bourdin   } else {
16214744f87SBlaise Bourdin     *sfNatural = NULL;
16314744f87SBlaise Bourdin   }
1649566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&sectionDist));
1659566063dSJacob Faibussowitsch   PetscCall(PetscFree(remoteOffsets));
1669566063dSJacob Faibussowitsch   if (destroyFlag) PetscCall(PetscSectionDestroy(&section));
1673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
168fa534816SMatthew G. Knepley }
169fa534816SMatthew G. Knepley 
170fa534816SMatthew G. Knepley /*@
17116635c05SJames Wright   DMPlexMigrateGlobalToNaturalSF - Migrates the input `sfNatural` based on sfMigration
17216635c05SJames Wright 
17316635c05SJames Wright   Input Parameters:
17416635c05SJames Wright + dmOld        - The original `DM`
17516635c05SJames Wright . dmNew        - The `DM` to be migrated to
17616635c05SJames Wright . sfNaturalOld - The sfNatural for the `dmOld`
17716635c05SJames Wright - sfMigration  - The `PetscSF` used to distribute the mesh, or `NULL` if it cannot be computed
17816635c05SJames Wright 
17916635c05SJames Wright   Output Parameter:
18016635c05SJames Wright . sfNaturalNew - `PetscSF` for mapping the `Vec` in PETSc ordering to the canonical ordering
18116635c05SJames Wright 
18216635c05SJames Wright   Level: intermediate
18316635c05SJames Wright 
18416635c05SJames Wright   Notes:
18516635c05SJames Wright   `sfNaturalOld` maps from the old Global section (roots) to the natural Vec layout (leaves, may or may not be described by a PetscSection).
18616635c05SJames Wright   `DMPlexMigrateGlobalToNaturalSF` creates an SF to map from the old global section to the new global section (generated from `sfMigration`).
18716635c05SJames Wright   That SF is then composed with the `sfNaturalOld` to generate `sfNaturalNew`.
18816635c05SJames Wright   This also distributes and sets the local section for `dmNew`.
18916635c05SJames Wright 
19016635c05SJames Wright   This is not typically called by the user.
19116635c05SJames Wright 
19216635c05SJames Wright .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `PetscSF`, `PetscSection`, `DMPlexDistribute()`, `DMPlexDistributeField()`
19316635c05SJames Wright  @*/
DMPlexMigrateGlobalToNaturalSF(DM dmOld,DM dmNew,PetscSF sfNaturalOld,PetscSF sfMigration,PetscSF * sfNaturalNew)19416635c05SJames Wright PetscErrorCode DMPlexMigrateGlobalToNaturalSF(DM dmOld, DM dmNew, PetscSF sfNaturalOld, PetscSF sfMigration, PetscSF *sfNaturalNew)
19516635c05SJames Wright {
19616635c05SJames Wright   MPI_Comm     comm;
197a2aa6d81SJames Wright   PetscSection oldGlobalSection, newGlobalSection;
19816635c05SJames Wright   PetscInt    *remoteOffsets;
19916635c05SJames Wright   PetscBool    debug = PETSC_FALSE;
20016635c05SJames Wright 
20116635c05SJames Wright   PetscFunctionBegin;
20216635c05SJames Wright   PetscCall(PetscObjectGetComm((PetscObject)dmNew, &comm));
20316635c05SJames Wright   if (!sfMigration) {
20416635c05SJames Wright     /* If sfMigration is missing, sfNatural cannot be computed and is set to NULL */
20516635c05SJames Wright     *sfNaturalNew = NULL;
20616635c05SJames Wright     PetscFunctionReturn(PETSC_SUCCESS);
20716635c05SJames Wright   }
20816635c05SJames Wright   if (debug) PetscCall(PetscSFView(sfMigration, NULL));
20916635c05SJames Wright 
21016635c05SJames Wright   { // Create oldGlobalSection and newGlobalSection *with* localOffsets
21116635c05SJames Wright     PetscSection oldLocalSection, newLocalSection;
21216635c05SJames Wright     PetscSF      pointSF;
21316635c05SJames Wright 
21416635c05SJames Wright     PetscCall(DMGetLocalSection(dmOld, &oldLocalSection));
21516635c05SJames Wright     PetscCall(DMGetPointSF(dmOld, &pointSF));
21616635c05SJames Wright     PetscCall(PetscSectionCreateGlobalSection(oldLocalSection, pointSF, PETSC_TRUE, PETSC_TRUE, PETSC_TRUE, &oldGlobalSection));
21716635c05SJames Wright 
21816635c05SJames Wright     PetscCall(PetscSectionCreate(comm, &newLocalSection));
21916635c05SJames Wright     PetscCall(PetscSFDistributeSection(sfMigration, oldLocalSection, NULL, newLocalSection));
22016635c05SJames Wright     PetscCall(DMSetLocalSection(dmNew, newLocalSection));
22116635c05SJames Wright 
22216635c05SJames Wright     PetscCall(DMGetPointSF(dmNew, &pointSF));
22316635c05SJames Wright     PetscCall(PetscSectionCreateGlobalSection(newLocalSection, pointSF, PETSC_TRUE, PETSC_TRUE, PETSC_TRUE, &newGlobalSection));
22416635c05SJames Wright 
22516635c05SJames Wright     PetscCall(PetscObjectSetName((PetscObject)oldLocalSection, "Old Local Section"));
22616635c05SJames Wright     if (debug) PetscCall(PetscSectionView(oldLocalSection, NULL));
22716635c05SJames Wright     PetscCall(PetscObjectSetName((PetscObject)oldGlobalSection, "Old Global Section"));
22816635c05SJames Wright     if (debug) PetscCall(PetscSectionView(oldGlobalSection, NULL));
22916635c05SJames Wright     PetscCall(PetscObjectSetName((PetscObject)newLocalSection, "New Local Section"));
23016635c05SJames Wright     if (debug) PetscCall(PetscSectionView(newLocalSection, NULL));
23116635c05SJames Wright     PetscCall(PetscObjectSetName((PetscObject)newGlobalSection, "New Global Section"));
23216635c05SJames Wright     if (debug) PetscCall(PetscSectionView(newGlobalSection, NULL));
23316635c05SJames Wright     PetscCall(PetscSectionDestroy(&newLocalSection));
23416635c05SJames Wright   }
23516635c05SJames Wright 
23616635c05SJames Wright   { // Create remoteOffsets array, mapping the oldGlobalSection offsets to the local points (according to sfMigration)
23716635c05SJames Wright     PetscInt lpStart, lpEnd, rpStart, rpEnd;
23816635c05SJames Wright 
23916635c05SJames Wright     PetscCall(PetscSectionGetChart(oldGlobalSection, &rpStart, &rpEnd));
24016635c05SJames Wright     PetscCall(PetscSectionGetChart(newGlobalSection, &lpStart, &lpEnd));
24116635c05SJames Wright 
24216635c05SJames Wright     // in `PetscSFDistributeSection` (where this is taken from), it possibly makes a new embedded SF. Should possibly do that here?
24316635c05SJames Wright     PetscCall(PetscMalloc1(lpEnd - lpStart, &remoteOffsets));
24416635c05SJames Wright     PetscCall(PetscSFBcastBegin(sfMigration, MPIU_INT, PetscSafePointerPlusOffset(oldGlobalSection->atlasOff, -rpStart), PetscSafePointerPlusOffset(remoteOffsets, -lpStart), MPI_REPLACE));
24516635c05SJames Wright     PetscCall(PetscSFBcastEnd(sfMigration, MPIU_INT, PetscSafePointerPlusOffset(oldGlobalSection->atlasOff, -rpStart), PetscSafePointerPlusOffset(remoteOffsets, -lpStart), MPI_REPLACE));
24616635c05SJames Wright     if (debug) {
24716635c05SJames Wright       PetscViewer viewer;
24816635c05SJames Wright 
24916635c05SJames Wright       PetscCall(PetscPrintf(comm, "Remote Offsets:\n"));
25016635c05SJames Wright       PetscCall(PetscViewerASCIIGetStdout(comm, &viewer));
25116635c05SJames Wright       PetscCall(PetscIntView(lpEnd - lpStart, remoteOffsets, viewer));
25216635c05SJames Wright     }
25316635c05SJames Wright   }
25416635c05SJames Wright 
25516635c05SJames Wright   { // Create SF from oldGlobalSection to newGlobalSection and compose with sfNaturalOld
25616635c05SJames Wright     PetscSF oldglob_to_newglob_sf, newglob_to_oldglob_sf;
25716635c05SJames Wright 
25816635c05SJames Wright     PetscCall(PetscSFCreateSectionSF(sfMigration, oldGlobalSection, remoteOffsets, newGlobalSection, &oldglob_to_newglob_sf));
25916635c05SJames Wright     PetscCall(PetscObjectSetName((PetscObject)oldglob_to_newglob_sf, "OldGlobal-to-NewGlobal SF"));
26016635c05SJames Wright     if (debug) PetscCall(PetscSFView(oldglob_to_newglob_sf, NULL));
26116635c05SJames Wright 
26216635c05SJames Wright     PetscCall(PetscSFCreateInverseSF(oldglob_to_newglob_sf, &newglob_to_oldglob_sf));
26316635c05SJames Wright     PetscCall(PetscObjectSetName((PetscObject)newglob_to_oldglob_sf, "NewGlobal-to-OldGlobal SF"));
26416635c05SJames Wright     PetscCall(PetscObjectViewFromOptions((PetscObject)newglob_to_oldglob_sf, (PetscObject)dmOld, "-natural_migrate_sf_view"));
26516635c05SJames Wright     PetscCall(PetscSFCompose(newglob_to_oldglob_sf, sfNaturalOld, sfNaturalNew));
26616635c05SJames Wright 
26716635c05SJames Wright     PetscCall(PetscSFDestroy(&oldglob_to_newglob_sf));
26816635c05SJames Wright     PetscCall(PetscSFDestroy(&newglob_to_oldglob_sf));
26916635c05SJames Wright   }
27016635c05SJames Wright 
27116635c05SJames Wright   PetscCall(PetscSectionDestroy(&oldGlobalSection));
27216635c05SJames Wright   PetscCall(PetscSectionDestroy(&newGlobalSection));
27316635c05SJames Wright   PetscCall(PetscFree(remoteOffsets));
27416635c05SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
27516635c05SJames Wright }
27616635c05SJames Wright 
27716635c05SJames Wright /*@
278a1cb98faSBarry Smith   DMPlexGlobalToNaturalBegin - Rearranges a global `Vec` in the natural order.
279fa534816SMatthew G. Knepley 
28020f4b53cSBarry Smith   Collective
281fa534816SMatthew G. Knepley 
282fa534816SMatthew G. Knepley   Input Parameters:
283a1cb98faSBarry Smith + dm - The distributed `DMPLEX`
284a1cb98faSBarry Smith - gv - The global `Vec`
285fa534816SMatthew G. Knepley 
2862fe279fdSBarry Smith   Output Parameter:
28720f4b53cSBarry Smith . nv - `Vec` in the canonical ordering distributed over all processors associated with `gv`
288fa534816SMatthew G. Knepley 
289fa534816SMatthew G. Knepley   Level: intermediate
290fa534816SMatthew G. Knepley 
291a1cb98faSBarry Smith   Note:
292a1cb98faSBarry Smith   The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`.
293a1cb98faSBarry Smith 
2941cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalEnd()`
295fa534816SMatthew G. Knepley @*/
DMPlexGlobalToNaturalBegin(DM dm,Vec gv,Vec nv)296d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGlobalToNaturalBegin(DM dm, Vec gv, Vec nv)
297d71ae5a4SJacob Faibussowitsch {
298fa534816SMatthew G. Knepley   const PetscScalar *inarray;
299fa534816SMatthew G. Knepley   PetscScalar       *outarray;
300f6c7c579SMatthew G. Knepley   MPI_Comm           comm;
301a6a55facSMatthew G. Knepley   PetscMPIInt        size;
302fa534816SMatthew G. Knepley 
303fa534816SMatthew G. Knepley   PetscFunctionBegin;
3049566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_GlobalToNaturalBegin, dm, 0, 0, 0));
305f6c7c579SMatthew G. Knepley   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
306f6c7c579SMatthew G. Knepley   PetscCallMPI(MPI_Comm_size(comm, &size));
307fa534816SMatthew G. Knepley   if (dm->sfNatural) {
308f6c7c579SMatthew G. Knepley     if (PetscDefined(USE_DEBUG)) {
309f6c7c579SMatthew G. Knepley       PetscSection gs;
310f6c7c579SMatthew G. Knepley       PetscInt     Nl, n;
311f6c7c579SMatthew G. Knepley 
312f6c7c579SMatthew G. Knepley       PetscCall(PetscSFGetGraph(dm->sfNatural, NULL, &Nl, NULL, NULL));
313f6c7c579SMatthew G. Knepley       PetscCall(VecGetLocalSize(nv, &n));
314f6c7c579SMatthew 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);
315f6c7c579SMatthew G. Knepley 
316f6c7c579SMatthew G. Knepley       PetscCall(DMGetGlobalSection(dm, &gs));
317f6c7c579SMatthew G. Knepley       PetscCall(PetscSectionGetConstrainedStorageSize(gs, &Nl));
318f6c7c579SMatthew G. Knepley       PetscCall(VecGetLocalSize(gv, &n));
319f6c7c579SMatthew 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);
320f6c7c579SMatthew G. Knepley     }
3219566063dSJacob Faibussowitsch     PetscCall(VecGetArray(nv, &outarray));
3229566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(gv, &inarray));
3239566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_REPLACE));
3249566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(gv, &inarray));
3259566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(nv, &outarray));
326a6a55facSMatthew G. Knepley   } else if (size == 1) {
3279566063dSJacob Faibussowitsch     PetscCall(VecCopy(gv, nv));
328f7d195e4SLawrence Mitchell   } else {
32900045ab3SPierre 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.");
33000045ab3SPierre Jolivet     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created. You must call DMSetUseNatural() before DMPlexDistribute().");
331f7d195e4SLawrence Mitchell   }
3329566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_GlobalToNaturalBegin, dm, 0, 0, 0));
3333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
334fa534816SMatthew G. Knepley }
335fa534816SMatthew G. Knepley 
336fa534816SMatthew G. Knepley /*@
337a1cb98faSBarry Smith   DMPlexGlobalToNaturalEnd - Rearranges a global `Vec` in the natural order.
338fa534816SMatthew G. Knepley 
33920f4b53cSBarry Smith   Collective
340fa534816SMatthew G. Knepley 
341fa534816SMatthew G. Knepley   Input Parameters:
342a1cb98faSBarry Smith + dm - The distributed `DMPLEX`
343a1cb98faSBarry Smith - gv - The global `Vec`
344fa534816SMatthew G. Knepley 
34597bb3fdcSJose E. Roman   Output Parameter:
346a1cb98faSBarry Smith . nv - The natural `Vec`
347fa534816SMatthew G. Knepley 
348fa534816SMatthew G. Knepley   Level: intermediate
349fa534816SMatthew G. Knepley 
350a1cb98faSBarry Smith   Note:
351a1cb98faSBarry Smith   The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`.
352a1cb98faSBarry Smith 
3531cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()`
354fa534816SMatthew G. Knepley  @*/
DMPlexGlobalToNaturalEnd(DM dm,Vec gv,Vec nv)355d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGlobalToNaturalEnd(DM dm, Vec gv, Vec nv)
356d71ae5a4SJacob Faibussowitsch {
357fa534816SMatthew G. Knepley   const PetscScalar *inarray;
358fa534816SMatthew G. Knepley   PetscScalar       *outarray;
359a6a55facSMatthew G. Knepley   PetscMPIInt        size;
360fa534816SMatthew G. Knepley 
361fa534816SMatthew G. Knepley   PetscFunctionBegin;
3629566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_GlobalToNaturalEnd, dm, 0, 0, 0));
3639566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
364fa534816SMatthew G. Knepley   if (dm->sfNatural) {
3659566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(gv, &inarray));
3669566063dSJacob Faibussowitsch     PetscCall(VecGetArray(nv, &outarray));
3679566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_REPLACE));
3689566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(gv, &inarray));
3699566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(nv, &outarray));
370a6a55facSMatthew G. Knepley   } else if (size == 1) {
371f7d195e4SLawrence Mitchell   } else {
37200045ab3SPierre 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.");
37300045ab3SPierre Jolivet     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created. You must call DMSetUseNatural() before DMPlexDistribute().");
374f7d195e4SLawrence Mitchell   }
3759566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_GlobalToNaturalEnd, dm, 0, 0, 0));
3763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
377fa534816SMatthew G. Knepley }
378fa534816SMatthew G. Knepley 
379fa534816SMatthew G. Knepley /*@
380a1cb98faSBarry Smith   DMPlexNaturalToGlobalBegin - Rearranges a `Vec` in the natural order to the Global order.
381fa534816SMatthew G. Knepley 
38220f4b53cSBarry Smith   Collective
383fa534816SMatthew G. Knepley 
384fa534816SMatthew G. Knepley   Input Parameters:
385a1cb98faSBarry Smith + dm - The distributed `DMPLEX`
386a1cb98faSBarry Smith - nv - The natural `Vec`
387fa534816SMatthew G. Knepley 
3882fe279fdSBarry Smith   Output Parameter:
389a1cb98faSBarry Smith . gv - The global `Vec`
390fa534816SMatthew G. Knepley 
391fa534816SMatthew G. Knepley   Level: intermediate
392fa534816SMatthew G. Knepley 
393a1cb98faSBarry Smith   Note:
394a1cb98faSBarry Smith   The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`.
395a1cb98faSBarry Smith 
39642747ad1SJacob Faibussowitsch .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexGlobalToNaturalEnd()`
397fa534816SMatthew G. Knepley @*/
DMPlexNaturalToGlobalBegin(DM dm,Vec nv,Vec gv)398d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexNaturalToGlobalBegin(DM dm, Vec nv, Vec gv)
399d71ae5a4SJacob Faibussowitsch {
400fa534816SMatthew G. Knepley   const PetscScalar *inarray;
401fa534816SMatthew G. Knepley   PetscScalar       *outarray;
402a6a55facSMatthew G. Knepley   PetscMPIInt        size;
403fa534816SMatthew G. Knepley 
404fa534816SMatthew G. Knepley   PetscFunctionBegin;
4059566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_NaturalToGlobalBegin, dm, 0, 0, 0));
4069566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
407fa534816SMatthew G. Knepley   if (dm->sfNatural) {
408a5b23f4aSJose E. Roman     /* We only have access to the SF that goes from Global to Natural.
409b64f75a9SBlaise Bourdin        Instead of inverting dm->sfNatural, we can call PetscSFReduceBegin/End with MPI_Op MPI_SUM.
410b64f75a9SBlaise Bourdin        Here the SUM really does nothing since sfNatural is one to one, as long as gV is set to zero first. */
4119566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(gv));
4129566063dSJacob Faibussowitsch     PetscCall(VecGetArray(gv, &outarray));
4139566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(nv, &inarray));
4149566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_SUM));
4159566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(nv, &inarray));
4169566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(gv, &outarray));
417a6a55facSMatthew G. Knepley   } else if (size == 1) {
4189566063dSJacob Faibussowitsch     PetscCall(VecCopy(nv, gv));
419f7d195e4SLawrence Mitchell   } else {
42000045ab3SPierre 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.");
42100045ab3SPierre Jolivet     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created. You must call DMSetUseNatural() before DMPlexDistribute().");
422f7d195e4SLawrence Mitchell   }
4239566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_NaturalToGlobalBegin, dm, 0, 0, 0));
4243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
425fa534816SMatthew G. Knepley }
426fa534816SMatthew G. Knepley 
427fa534816SMatthew G. Knepley /*@
428a1cb98faSBarry Smith   DMPlexNaturalToGlobalEnd - Rearranges a `Vec` in the natural order to the Global order.
429fa534816SMatthew G. Knepley 
43020f4b53cSBarry Smith   Collective
431fa534816SMatthew G. Knepley 
432fa534816SMatthew G. Knepley   Input Parameters:
433a1cb98faSBarry Smith + dm - The distributed `DMPLEX`
434a1cb98faSBarry Smith - nv - The natural `Vec`
435fa534816SMatthew G. Knepley 
4362fe279fdSBarry Smith   Output Parameter:
437a1cb98faSBarry Smith . gv - The global `Vec`
438fa534816SMatthew G. Knepley 
439fa534816SMatthew G. Knepley   Level: intermediate
440fa534816SMatthew G. Knepley 
441a1cb98faSBarry Smith   Note:
442a1cb98faSBarry Smith   The user must call `DMSetUseNatural`(dm, `PETSC_TRUE`) before `DMPlexDistribute()`.
443a1cb98faSBarry Smith 
4441cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `Vec`, `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()`
445fa534816SMatthew G. Knepley  @*/
DMPlexNaturalToGlobalEnd(DM dm,Vec nv,Vec gv)446d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexNaturalToGlobalEnd(DM dm, Vec nv, Vec gv)
447d71ae5a4SJacob Faibussowitsch {
448fa534816SMatthew G. Knepley   const PetscScalar *inarray;
449fa534816SMatthew G. Knepley   PetscScalar       *outarray;
450a6a55facSMatthew G. Knepley   PetscMPIInt        size;
451fa534816SMatthew G. Knepley 
452fa534816SMatthew G. Knepley   PetscFunctionBegin;
4539566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_NaturalToGlobalEnd, dm, 0, 0, 0));
4549566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
455fa534816SMatthew G. Knepley   if (dm->sfNatural) {
4569566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(nv, &inarray));
4579566063dSJacob Faibussowitsch     PetscCall(VecGetArray(gv, &outarray));
4589566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_SUM));
4599566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(nv, &inarray));
4609566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(gv, &outarray));
461a6a55facSMatthew G. Knepley   } else if (size == 1) {
462f7d195e4SLawrence Mitchell   } else {
46300045ab3SPierre 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.");
46400045ab3SPierre Jolivet     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created. You must call DMSetUseNatural() before DMPlexDistribute().");
465f7d195e4SLawrence Mitchell   }
4669566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_NaturalToGlobalEnd, dm, 0, 0, 0));
4673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
468fa534816SMatthew G. Knepley }
46909cf685aSAlexis Marboeuf 
47009cf685aSAlexis Marboeuf /*@
4715cb80ecdSBarry Smith   DMPlexCreateNaturalVector - Provide a `Vec` capable of holding the natural ordering and distribution.
47209cf685aSAlexis Marboeuf 
47320f4b53cSBarry Smith   Collective
47409cf685aSAlexis Marboeuf 
4755cb80ecdSBarry Smith   Input Parameter:
476a1cb98faSBarry Smith . dm - The distributed `DMPLEX`
47709cf685aSAlexis Marboeuf 
4785cb80ecdSBarry Smith   Output Parameter:
4795cb80ecdSBarry Smith . nv - The natural `Vec`
48009cf685aSAlexis Marboeuf 
48109cf685aSAlexis Marboeuf   Level: intermediate
48209cf685aSAlexis Marboeuf 
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()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()`
48709cf685aSAlexis Marboeuf  @*/
DMPlexCreateNaturalVector(DM dm,Vec * nv)488d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreateNaturalVector(DM dm, Vec *nv)
489d71ae5a4SJacob Faibussowitsch {
49009cf685aSAlexis Marboeuf   PetscMPIInt size;
49109cf685aSAlexis Marboeuf 
49209cf685aSAlexis Marboeuf   PetscFunctionBegin;
49309cf685aSAlexis Marboeuf   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
49409cf685aSAlexis Marboeuf   if (dm->sfNatural) {
4956fba1fe0SBlaise Bourdin     PetscInt nleaves, bs, maxbs;
4962e8d78feSBlaise Bourdin     Vec      v;
4976fba1fe0SBlaise Bourdin 
4986fba1fe0SBlaise Bourdin     /*
4996fba1fe0SBlaise Bourdin       Setting the natural vector block size.
5006fba1fe0SBlaise Bourdin       We can't get it from a global vector because of constraints, and the block size in the local vector
5016fba1fe0SBlaise Bourdin       may be inconsistent across processes, typically when some local vectors have size 0, their block size is set to 1
5026fba1fe0SBlaise Bourdin     */
5032e8d78feSBlaise Bourdin     PetscCall(DMGetLocalVector(dm, &v));
5042e8d78feSBlaise Bourdin     PetscCall(VecGetBlockSize(v, &bs));
505462c564dSBarry Smith     PetscCallMPI(MPIU_Allreduce(&bs, &maxbs, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm)));
5066fba1fe0SBlaise Bourdin     if (bs == 1 && maxbs > 1) bs = maxbs;
5072e8d78feSBlaise Bourdin     PetscCall(DMRestoreLocalVector(dm, &v));
5082e8d78feSBlaise Bourdin 
50909cf685aSAlexis Marboeuf     PetscCall(PetscSFGetGraph(dm->sfNatural, NULL, &nleaves, NULL, NULL));
51009cf685aSAlexis Marboeuf     PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), nv));
51109cf685aSAlexis Marboeuf     PetscCall(VecSetSizes(*nv, nleaves, PETSC_DETERMINE));
51209cf685aSAlexis Marboeuf     PetscCall(VecSetBlockSize(*nv, bs));
51309cf685aSAlexis Marboeuf     PetscCall(VecSetType(*nv, dm->vectype));
51409cf685aSAlexis Marboeuf     PetscCall(VecSetDM(*nv, dm));
51509cf685aSAlexis Marboeuf   } else if (size == 1) {
5162e8d78feSBlaise Bourdin     PetscCall(DMCreateLocalVector(dm, nv));
51709cf685aSAlexis Marboeuf   } else {
51800045ab3SPierre 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.");
51900045ab3SPierre Jolivet     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created. You must call DMSetUseNatural() before DMPlexDistribute().");
52009cf685aSAlexis Marboeuf   }
5213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
52209cf685aSAlexis Marboeuf }
523