xref: /petsc/src/dm/impls/plex/plexnatural.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1fa534816SMatthew G. Knepley #include <petsc/private/dmpleximpl.h> /*I      "petscdmplex.h"   I*/
2fa534816SMatthew G. Knepley 
3fa534816SMatthew G. Knepley /*@
4f94b4a02SBlaise Bourdin   DMPlexSetMigrationSF - Sets the SF for migrating from a parent DM into this DM
5f94b4a02SBlaise Bourdin 
65d3b26e6SMatthew G. Knepley   Input Parameters:
7f94b4a02SBlaise Bourdin + dm        - The DM
85d3b26e6SMatthew G. Knepley - naturalSF - The PetscSF
95d3b26e6SMatthew G. Knepley 
105d3b26e6SMatthew G. Knepley   Note: It is necessary to call this in order to have DMCreateSubDM() or DMCreateSuperDM() build the Global-To-Natural map
115d3b26e6SMatthew G. Knepley 
12f94b4a02SBlaise Bourdin   Level: intermediate
13f94b4a02SBlaise Bourdin 
14db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexCreateMigrationSF()`, `DMPlexGetMigrationSF()`
15f94b4a02SBlaise Bourdin @*/
16*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexSetMigrationSF(DM dm, PetscSF migrationSF)
17*d71ae5a4SJacob Faibussowitsch {
18f94b4a02SBlaise Bourdin   PetscFunctionBegin;
19f94b4a02SBlaise Bourdin   dm->sfMigration = migrationSF;
209566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)migrationSF));
21f94b4a02SBlaise Bourdin   PetscFunctionReturn(0);
22f94b4a02SBlaise Bourdin }
23f94b4a02SBlaise Bourdin 
24f94b4a02SBlaise Bourdin /*@
25f94b4a02SBlaise Bourdin   DMPlexGetMigrationSF - Gets the SF for migrating from a parent DM into this DM
26f94b4a02SBlaise Bourdin 
275d3b26e6SMatthew G. Knepley   Input Parameter:
285d3b26e6SMatthew G. Knepley . dm          - The DM
295d3b26e6SMatthew G. Knepley 
305d3b26e6SMatthew G. Knepley   Output Parameter:
315d3b26e6SMatthew G. Knepley . migrationSF - The PetscSF
325d3b26e6SMatthew G. Knepley 
33f94b4a02SBlaise Bourdin   Level: intermediate
34f94b4a02SBlaise Bourdin 
35db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexCreateMigrationSF()`, `DMPlexSetMigrationSF`
36f94b4a02SBlaise Bourdin @*/
37*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGetMigrationSF(DM dm, PetscSF *migrationSF)
38*d71ae5a4SJacob Faibussowitsch {
39f94b4a02SBlaise Bourdin   PetscFunctionBegin;
40f94b4a02SBlaise Bourdin   *migrationSF = dm->sfMigration;
41f94b4a02SBlaise Bourdin   PetscFunctionReturn(0);
42f94b4a02SBlaise Bourdin }
43f94b4a02SBlaise Bourdin 
44f94b4a02SBlaise Bourdin /*@
45f94b4a02SBlaise Bourdin   DMPlexSetGlobalToNaturalSF - Sets the SF for mapping Global Vec to the Natural Vec
46f94b4a02SBlaise Bourdin 
475d3b26e6SMatthew G. Knepley   Input Parameters:
48f94b4a02SBlaise Bourdin + dm          - The DM
495d3b26e6SMatthew G. Knepley - naturalSF   - The PetscSF
505d3b26e6SMatthew G. Knepley 
51f94b4a02SBlaise Bourdin   Level: intermediate
52f94b4a02SBlaise Bourdin 
53db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexCreateGlobalToNaturalSF()`, `DMPlexGetGlobaltoNaturalSF()`
54f94b4a02SBlaise Bourdin @*/
55*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexSetGlobalToNaturalSF(DM dm, PetscSF naturalSF)
56*d71ae5a4SJacob Faibussowitsch {
57f94b4a02SBlaise Bourdin   PetscFunctionBegin;
58f94b4a02SBlaise Bourdin   dm->sfNatural = naturalSF;
599566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)naturalSF));
60f94b4a02SBlaise Bourdin   dm->useNatural = PETSC_TRUE;
61f94b4a02SBlaise Bourdin   PetscFunctionReturn(0);
62f94b4a02SBlaise Bourdin }
63f94b4a02SBlaise Bourdin 
64f94b4a02SBlaise Bourdin /*@
65f94b4a02SBlaise Bourdin   DMPlexGetGlobalToNaturalSF - Gets the SF for mapping Global Vec to the Natural Vec
66f94b4a02SBlaise Bourdin 
675d3b26e6SMatthew G. Knepley   Input Parameter:
685d3b26e6SMatthew G. Knepley . dm          - The DM
695d3b26e6SMatthew G. Knepley 
705d3b26e6SMatthew G. Knepley   Output Parameter:
715d3b26e6SMatthew G. Knepley . naturalSF   - The PetscSF
725d3b26e6SMatthew G. Knepley 
73f94b4a02SBlaise Bourdin   Level: intermediate
74f94b4a02SBlaise Bourdin 
75db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexCreateGlobalToNaturalSF()`, `DMPlexSetGlobaltoNaturalSF`
76f94b4a02SBlaise Bourdin @*/
77*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGetGlobalToNaturalSF(DM dm, PetscSF *naturalSF)
78*d71ae5a4SJacob Faibussowitsch {
79f94b4a02SBlaise Bourdin   PetscFunctionBegin;
80f94b4a02SBlaise Bourdin   *naturalSF = dm->sfNatural;
81f94b4a02SBlaise Bourdin   PetscFunctionReturn(0);
82f94b4a02SBlaise Bourdin }
83f94b4a02SBlaise Bourdin 
84f94b4a02SBlaise Bourdin /*@
85fa534816SMatthew G. Knepley   DMPlexCreateGlobalToNaturalSF - Creates the SF for mapping Global Vec to the Natural Vec
86fa534816SMatthew G. Knepley 
87fa534816SMatthew G. Knepley   Input Parameters:
88e7858701SMatthew G. Knepley + dm          - The redistributed DM
89e7858701SMatthew G. Knepley . section     - The local PetscSection describing the Vec before the mesh was distributed, or NULL if not available
905c3f5608SAlexisMarb - sfMigration - The PetscSF used to distribute the mesh, or NULL if it cannot be computed
91fa534816SMatthew G. Knepley 
925d3b26e6SMatthew G. Knepley   Output Parameter:
93fa534816SMatthew G. Knepley . sfNatural   - PetscSF for mapping the Vec in PETSc ordering to the canonical ordering
94fa534816SMatthew G. Knepley 
955d3b26e6SMatthew G. Knepley   Note: This is not typically called by the user.
965d3b26e6SMatthew G. Knepley 
97fa534816SMatthew G. Knepley   Level: intermediate
98fa534816SMatthew G. Knepley 
99db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`
100fa534816SMatthew G. Knepley  @*/
101*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreateGlobalToNaturalSF(DM dm, PetscSection section, PetscSF sfMigration, PetscSF *sfNatural)
102*d71ae5a4SJacob Faibussowitsch {
103fa534816SMatthew G. Knepley   MPI_Comm     comm;
104e7858701SMatthew G. Knepley   PetscSF      sf, sfEmbed, sfField;
105e5b44f4fSMatthew G. Knepley   PetscSection gSection, sectionDist, gLocSection;
106fa534816SMatthew G. Knepley   PetscInt    *spoints, *remoteOffsets;
1075c3f5608SAlexisMarb   PetscInt     ssize, pStart, pEnd, p, globalSize;
108e7858701SMatthew G. Knepley   PetscBool    destroyFlag = PETSC_FALSE, debug = PETSC_FALSE;
109fa534816SMatthew G. Knepley 
110fa534816SMatthew G. Knepley   PetscFunctionBegin;
1119566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
1125c3f5608SAlexisMarb   if (!sfMigration) {
113e7858701SMatthew G. Knepley     /* If sfMigration is missing, sfNatural cannot be computed and is set to NULL */
1145c3f5608SAlexisMarb     *sfNatural = NULL;
1155c3f5608SAlexisMarb     PetscFunctionReturn(0);
1165c3f5608SAlexisMarb   } else if (!section) {
117e7858701SMatthew G. Knepley     /* If the sequential section is not provided (NULL), it is reconstructed from the parallel section */
1185c3f5608SAlexisMarb     PetscSF      sfMigrationInv;
1195c3f5608SAlexisMarb     PetscSection localSection;
1205c3f5608SAlexisMarb 
1219566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dm, &localSection));
1229566063dSJacob Faibussowitsch     PetscCall(PetscSFCreateInverseSF(sfMigration, &sfMigrationInv));
1239566063dSJacob Faibussowitsch     PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), &section));
1249566063dSJacob Faibussowitsch     PetscCall(PetscSFDistributeSection(sfMigrationInv, localSection, NULL, section));
1259566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&sfMigrationInv));
1265c3f5608SAlexisMarb     destroyFlag = PETSC_TRUE;
1275c3f5608SAlexisMarb   }
128e7858701SMatthew G. Knepley   if (debug) PetscCall(PetscSFView(sfMigration, NULL));
129e5b44f4fSMatthew G. Knepley   /* Create a new section from distributing the original section */
1309566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(comm, &sectionDist));
1319566063dSJacob Faibussowitsch   PetscCall(PetscSFDistributeSection(sfMigration, section, &remoteOffsets, sectionDist));
132e7858701SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)sectionDist, "Migrated Section"));
133e7858701SMatthew G. Knepley   if (debug) PetscCall(PetscSectionView(sectionDist, NULL));
1349566063dSJacob Faibussowitsch   PetscCall(DMSetLocalSection(dm, sectionDist));
135e7858701SMatthew G. Knepley   /* If a sequential section is provided but no dof is affected, sfNatural cannot be computed and is set to NULL */
13699a026b9SAlexis Marboeuf   PetscCall(PetscSectionGetStorageSize(sectionDist, &globalSize));
13714744f87SBlaise Bourdin   if (globalSize) {
138e7858701SMatthew G. Knepley     const PetscInt *leaves;
139e7858701SMatthew G. Knepley     PetscInt       *sortleaves, *indices;
140e7858701SMatthew G. Knepley     PetscInt        Nl;
141e7858701SMatthew G. Knepley 
142fa534816SMatthew G. Knepley     /* Get a pruned version of migration SF */
1439566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalSection(dm, &gSection));
144e7858701SMatthew G. Knepley     if (debug) PetscCall(PetscSectionView(gSection, NULL));
145e7858701SMatthew G. Knepley     PetscCall(PetscSFGetGraph(sfMigration, NULL, &Nl, &leaves, NULL));
1469566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetChart(gSection, &pStart, &pEnd));
147fa534816SMatthew G. Knepley     for (p = pStart, ssize = 0; p < pEnd; ++p) {
148fa534816SMatthew G. Knepley       PetscInt dof, off;
149fa534816SMatthew G. Knepley 
1509566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(gSection, p, &dof));
1519566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(gSection, p, &off));
152fa534816SMatthew G. Knepley       if ((dof > 0) && (off >= 0)) ++ssize;
153fa534816SMatthew G. Knepley     }
154e7858701SMatthew G. Knepley     PetscCall(PetscMalloc3(ssize, &spoints, Nl, &sortleaves, Nl, &indices));
1559371c9d4SSatish Balay     for (p = 0; p < Nl; ++p) {
1569371c9d4SSatish Balay       sortleaves[p] = leaves ? leaves[p] : p;
1579371c9d4SSatish Balay       indices[p]    = p;
1589371c9d4SSatish Balay     }
159e7858701SMatthew G. Knepley     PetscCall(PetscSortIntWithArray(Nl, sortleaves, indices));
160fa534816SMatthew G. Knepley     for (p = pStart, ssize = 0; p < pEnd; ++p) {
161e7858701SMatthew G. Knepley       PetscInt dof, off, loc;
162fa534816SMatthew G. Knepley 
1639566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(gSection, p, &dof));
1649566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(gSection, p, &off));
165e7858701SMatthew G. Knepley       if ((dof > 0) && (off >= 0)) {
166e7858701SMatthew G. Knepley         PetscCall(PetscFindInt(p, Nl, sortleaves, &loc));
167e7858701SMatthew 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);
168e7858701SMatthew G. Knepley         spoints[ssize++] = indices[loc];
169e7858701SMatthew G. Knepley       }
170fa534816SMatthew G. Knepley     }
1719566063dSJacob Faibussowitsch     PetscCall(PetscSFCreateEmbeddedLeafSF(sfMigration, ssize, spoints, &sfEmbed));
172e7858701SMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)sfEmbed, "Embedded SF"));
173e7858701SMatthew G. Knepley     PetscCall(PetscFree3(spoints, sortleaves, indices));
174e7858701SMatthew G. Knepley     if (debug) PetscCall(PetscSFView(sfEmbed, NULL));
175e7858701SMatthew G. Knepley     /* Create the SF associated with this section
176e7858701SMatthew G. Knepley          Roots are natural dofs, leaves are global dofs */
1779566063dSJacob Faibussowitsch     PetscCall(DMGetPointSF(dm, &sf));
17899a026b9SAlexis Marboeuf     PetscCall(PetscSectionCreateGlobalSection(sectionDist, sf, PETSC_TRUE, PETSC_TRUE, &gLocSection));
1799566063dSJacob Faibussowitsch     PetscCall(PetscSFCreateSectionSF(sfEmbed, section, remoteOffsets, gLocSection, &sfField));
1809566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&sfEmbed));
1819566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&gLocSection));
182e7858701SMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)sfField, "Natural-to-Global SF"));
183e7858701SMatthew G. Knepley     if (debug) PetscCall(PetscSFView(sfField, NULL));
184e7858701SMatthew G. Knepley     /* Invert the field SF
185e7858701SMatthew G. Knepley          Roots are global dofs, leaves are natural dofs */
186e7858701SMatthew G. Knepley     PetscCall(PetscSFCreateInverseSF(sfField, sfNatural));
187e7858701SMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)*sfNatural, "Global-to-Natural SF"));
1889566063dSJacob Faibussowitsch     PetscCall(PetscObjectViewFromOptions((PetscObject)*sfNatural, NULL, "-globaltonatural_sf_view"));
189fa534816SMatthew G. Knepley     /* Clean up */
190e7858701SMatthew G. Knepley     PetscCall(PetscSFDestroy(&sfField));
19114744f87SBlaise Bourdin   } else {
19214744f87SBlaise Bourdin     *sfNatural = NULL;
19314744f87SBlaise Bourdin   }
1949566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&sectionDist));
1959566063dSJacob Faibussowitsch   PetscCall(PetscFree(remoteOffsets));
1969566063dSJacob Faibussowitsch   if (destroyFlag) PetscCall(PetscSectionDestroy(&section));
197fa534816SMatthew G. Knepley   PetscFunctionReturn(0);
198fa534816SMatthew G. Knepley }
199fa534816SMatthew G. Knepley 
200fa534816SMatthew G. Knepley /*@
201fa534816SMatthew G. Knepley   DMPlexGlobalToNaturalBegin - Rearranges a global Vector in the natural order.
202fa534816SMatthew G. Knepley 
203fa534816SMatthew G. Knepley   Collective on dm
204fa534816SMatthew G. Knepley 
205fa534816SMatthew G. Knepley   Input Parameters:
206fa534816SMatthew G. Knepley + dm - The distributed DMPlex
207fa534816SMatthew G. Knepley - gv - The global Vec
208fa534816SMatthew G. Knepley 
209fa534816SMatthew G. Knepley   Output Parameters:
210fa534816SMatthew G. Knepley . nv - Vec in the canonical ordering distributed over all processors associated with gv
211fa534816SMatthew G. Knepley 
21242ea106eSTristan Konolige   Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute().
213fa534816SMatthew G. Knepley 
214fa534816SMatthew G. Knepley   Level: intermediate
215fa534816SMatthew G. Knepley 
216db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalEnd()`
217fa534816SMatthew G. Knepley @*/
218*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGlobalToNaturalBegin(DM dm, Vec gv, Vec nv)
219*d71ae5a4SJacob Faibussowitsch {
220fa534816SMatthew G. Knepley   const PetscScalar *inarray;
221fa534816SMatthew G. Knepley   PetscScalar       *outarray;
222a6a55facSMatthew G. Knepley   PetscMPIInt        size;
223fa534816SMatthew G. Knepley 
224fa534816SMatthew G. Knepley   PetscFunctionBegin;
2259566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_GlobalToNaturalBegin, dm, 0, 0, 0));
2269566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
227fa534816SMatthew G. Knepley   if (dm->sfNatural) {
2289566063dSJacob Faibussowitsch     PetscCall(VecGetArray(nv, &outarray));
2299566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(gv, &inarray));
2309566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_REPLACE));
2319566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(gv, &inarray));
2329566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(nv, &outarray));
233a6a55facSMatthew G. Knepley   } else if (size == 1) {
2349566063dSJacob Faibussowitsch     PetscCall(VecCopy(gv, nv));
235f7d195e4SLawrence Mitchell   } else {
236f7d195e4SLawrence Mitchell     PetscCheck(!dm->useNatural, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "DM global to natural SF not present.\nIf DMPlexDistribute() was called and a section was defined, report to petsc-maint@mcs.anl.gov.");
237f7d195e4SLawrence Mitchell     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().");
238f7d195e4SLawrence Mitchell   }
2399566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_GlobalToNaturalBegin, dm, 0, 0, 0));
240fa534816SMatthew G. Knepley   PetscFunctionReturn(0);
241fa534816SMatthew G. Knepley }
242fa534816SMatthew G. Knepley 
243fa534816SMatthew G. Knepley /*@
244fa534816SMatthew G. Knepley   DMPlexGlobalToNaturalEnd - Rearranges a global Vector in the natural order.
245fa534816SMatthew G. Knepley 
246fa534816SMatthew G. Knepley   Collective on dm
247fa534816SMatthew G. Knepley 
248fa534816SMatthew G. Knepley   Input Parameters:
249fa534816SMatthew G. Knepley + dm - The distributed DMPlex
250fa534816SMatthew G. Knepley - gv - The global Vec
251fa534816SMatthew G. Knepley 
25297bb3fdcSJose E. Roman   Output Parameter:
253fa534816SMatthew G. Knepley . nv - The natural Vec
254fa534816SMatthew G. Knepley 
25542ea106eSTristan Konolige   Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute().
256fa534816SMatthew G. Knepley 
257fa534816SMatthew G. Knepley   Level: intermediate
258fa534816SMatthew G. Knepley 
259db781477SPatrick Sanan  .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()`
260fa534816SMatthew G. Knepley  @*/
261*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGlobalToNaturalEnd(DM dm, Vec gv, Vec nv)
262*d71ae5a4SJacob Faibussowitsch {
263fa534816SMatthew G. Knepley   const PetscScalar *inarray;
264fa534816SMatthew G. Knepley   PetscScalar       *outarray;
265a6a55facSMatthew G. Knepley   PetscMPIInt        size;
266fa534816SMatthew G. Knepley 
267fa534816SMatthew G. Knepley   PetscFunctionBegin;
2689566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_GlobalToNaturalEnd, dm, 0, 0, 0));
2699566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
270fa534816SMatthew G. Knepley   if (dm->sfNatural) {
2719566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(gv, &inarray));
2729566063dSJacob Faibussowitsch     PetscCall(VecGetArray(nv, &outarray));
2739566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_REPLACE));
2749566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(gv, &inarray));
2759566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(nv, &outarray));
276a6a55facSMatthew G. Knepley   } else if (size == 1) {
277f7d195e4SLawrence Mitchell   } else {
278f7d195e4SLawrence Mitchell     PetscCheck(!dm->useNatural, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "DM global to natural SF not present.\nIf DMPlexDistribute() was called and a section was defined, report to petsc-maint@mcs.anl.gov.");
279f7d195e4SLawrence Mitchell     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().");
280f7d195e4SLawrence Mitchell   }
2819566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_GlobalToNaturalEnd, dm, 0, 0, 0));
282fa534816SMatthew G. Knepley   PetscFunctionReturn(0);
283fa534816SMatthew G. Knepley }
284fa534816SMatthew G. Knepley 
285fa534816SMatthew G. Knepley /*@
286fa534816SMatthew G. Knepley   DMPlexNaturalToGlobalBegin - Rearranges a Vector in the natural order to the Global order.
287fa534816SMatthew G. Knepley 
288fa534816SMatthew G. Knepley   Collective on dm
289fa534816SMatthew G. Knepley 
290fa534816SMatthew G. Knepley   Input Parameters:
291fa534816SMatthew G. Knepley + dm - The distributed DMPlex
292fa534816SMatthew G. Knepley - nv - The natural Vec
293fa534816SMatthew G. Knepley 
294fa534816SMatthew G. Knepley   Output Parameters:
295fa534816SMatthew G. Knepley . gv - The global Vec
296fa534816SMatthew G. Knepley 
29742ea106eSTristan Konolige   Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute().
298fa534816SMatthew G. Knepley 
299fa534816SMatthew G. Knepley   Level: intermediate
300fa534816SMatthew G. Knepley 
301c2e3fba1SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalEnd()`
302fa534816SMatthew G. Knepley @*/
303*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexNaturalToGlobalBegin(DM dm, Vec nv, Vec gv)
304*d71ae5a4SJacob Faibussowitsch {
305fa534816SMatthew G. Knepley   const PetscScalar *inarray;
306fa534816SMatthew G. Knepley   PetscScalar       *outarray;
307a6a55facSMatthew G. Knepley   PetscMPIInt        size;
308fa534816SMatthew G. Knepley 
309fa534816SMatthew G. Knepley   PetscFunctionBegin;
3109566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_NaturalToGlobalBegin, dm, 0, 0, 0));
3119566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
312fa534816SMatthew G. Knepley   if (dm->sfNatural) {
313a5b23f4aSJose E. Roman     /* We only have access to the SF that goes from Global to Natural.
314b64f75a9SBlaise Bourdin        Instead of inverting dm->sfNatural, we can call PetscSFReduceBegin/End with MPI_Op MPI_SUM.
315b64f75a9SBlaise Bourdin        Here the SUM really does nothing since sfNatural is one to one, as long as gV is set to zero first. */
3169566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(gv));
3179566063dSJacob Faibussowitsch     PetscCall(VecGetArray(gv, &outarray));
3189566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(nv, &inarray));
3199566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_SUM));
3209566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(nv, &inarray));
3219566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(gv, &outarray));
322a6a55facSMatthew G. Knepley   } else if (size == 1) {
3239566063dSJacob Faibussowitsch     PetscCall(VecCopy(nv, gv));
324f7d195e4SLawrence Mitchell   } else {
325f7d195e4SLawrence Mitchell     PetscCheck(!dm->useNatural, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "DM global to natural SF not present.\nIf DMPlexDistribute() was called and a section was defined, report to petsc-maint@mcs.anl.gov.");
326f7d195e4SLawrence Mitchell     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().");
327f7d195e4SLawrence Mitchell   }
3289566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_NaturalToGlobalBegin, dm, 0, 0, 0));
329fa534816SMatthew G. Knepley   PetscFunctionReturn(0);
330fa534816SMatthew G. Knepley }
331fa534816SMatthew G. Knepley 
332fa534816SMatthew G. Knepley /*@
333fa534816SMatthew G. Knepley   DMPlexNaturalToGlobalEnd - Rearranges a Vector in the natural order to the Global order.
334fa534816SMatthew G. Knepley 
335fa534816SMatthew G. Knepley   Collective on dm
336fa534816SMatthew G. Knepley 
337fa534816SMatthew G. Knepley   Input Parameters:
338fa534816SMatthew G. Knepley + dm - The distributed DMPlex
339fa534816SMatthew G. Knepley - nv - The natural Vec
340fa534816SMatthew G. Knepley 
341fa534816SMatthew G. Knepley   Output Parameters:
342fa534816SMatthew G. Knepley . gv - The global Vec
343fa534816SMatthew G. Knepley 
34442ea106eSTristan Konolige   Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute().
345fa534816SMatthew G. Knepley 
346fa534816SMatthew G. Knepley   Level: intermediate
347fa534816SMatthew G. Knepley 
348db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()`
349fa534816SMatthew G. Knepley  @*/
350*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexNaturalToGlobalEnd(DM dm, Vec nv, Vec gv)
351*d71ae5a4SJacob Faibussowitsch {
352fa534816SMatthew G. Knepley   const PetscScalar *inarray;
353fa534816SMatthew G. Knepley   PetscScalar       *outarray;
354a6a55facSMatthew G. Knepley   PetscMPIInt        size;
355fa534816SMatthew G. Knepley 
356fa534816SMatthew G. Knepley   PetscFunctionBegin;
3579566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_NaturalToGlobalEnd, dm, 0, 0, 0));
3589566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
359fa534816SMatthew G. Knepley   if (dm->sfNatural) {
3609566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(nv, &inarray));
3619566063dSJacob Faibussowitsch     PetscCall(VecGetArray(gv, &outarray));
3629566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_SUM));
3639566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(nv, &inarray));
3649566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(gv, &outarray));
365a6a55facSMatthew G. Knepley   } else if (size == 1) {
366f7d195e4SLawrence Mitchell   } else {
367f7d195e4SLawrence Mitchell     PetscCheck(!dm->useNatural, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "DM global to natural SF not present.\nIf DMPlexDistribute() was called and a section was defined, report to petsc-maint@mcs.anl.gov.");
368f7d195e4SLawrence Mitchell     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().");
369f7d195e4SLawrence Mitchell   }
3709566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_NaturalToGlobalEnd, dm, 0, 0, 0));
371fa534816SMatthew G. Knepley   PetscFunctionReturn(0);
372fa534816SMatthew G. Knepley }
37309cf685aSAlexis Marboeuf 
37409cf685aSAlexis Marboeuf /*@
3755cb80ecdSBarry Smith   DMPlexCreateNaturalVector - Provide a `Vec` capable of holding the natural ordering and distribution.
37609cf685aSAlexis Marboeuf 
37709cf685aSAlexis Marboeuf   Collective on dm
37809cf685aSAlexis Marboeuf 
3795cb80ecdSBarry Smith   Input Parameter:
3805cb80ecdSBarry Smith . dm - The distributed `DMPlex`
38109cf685aSAlexis Marboeuf 
3825cb80ecdSBarry Smith   Output Parameter:
3835cb80ecdSBarry Smith . nv - The natural `Vec`
38409cf685aSAlexis Marboeuf 
3855cb80ecdSBarry Smith   Note: The user must call `DMSetUseNatural`(dm, PETSC_TRUE) before `DMPlexDistribute()`.
38609cf685aSAlexis Marboeuf 
38709cf685aSAlexis Marboeuf   Level: intermediate
38809cf685aSAlexis Marboeuf 
38909cf685aSAlexis Marboeuf .seealso: `DMPlexDistribute()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()`
39009cf685aSAlexis Marboeuf  @*/
391*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreateNaturalVector(DM dm, Vec *nv)
392*d71ae5a4SJacob Faibussowitsch {
39309cf685aSAlexis Marboeuf   PetscMPIInt size;
39409cf685aSAlexis Marboeuf 
39509cf685aSAlexis Marboeuf   PetscFunctionBegin;
39609cf685aSAlexis Marboeuf   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
39709cf685aSAlexis Marboeuf   if (dm->sfNatural) {
3982e8d78feSBlaise Bourdin     PetscInt nleaves, bs;
3992e8d78feSBlaise Bourdin     Vec      v;
4002e8d78feSBlaise Bourdin     PetscCall(DMGetLocalVector(dm, &v));
4012e8d78feSBlaise Bourdin     PetscCall(VecGetBlockSize(v, &bs));
4022e8d78feSBlaise Bourdin     PetscCall(DMRestoreLocalVector(dm, &v));
4032e8d78feSBlaise Bourdin 
40409cf685aSAlexis Marboeuf     PetscCall(PetscSFGetGraph(dm->sfNatural, NULL, &nleaves, NULL, NULL));
40509cf685aSAlexis Marboeuf     PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), nv));
40609cf685aSAlexis Marboeuf     PetscCall(VecSetSizes(*nv, nleaves, PETSC_DETERMINE));
40709cf685aSAlexis Marboeuf     PetscCall(VecSetBlockSize(*nv, bs));
40809cf685aSAlexis Marboeuf     PetscCall(VecSetType(*nv, dm->vectype));
40909cf685aSAlexis Marboeuf     PetscCall(VecSetDM(*nv, dm));
41009cf685aSAlexis Marboeuf   } else if (size == 1) {
4112e8d78feSBlaise Bourdin     PetscCall(DMCreateLocalVector(dm, nv));
41209cf685aSAlexis Marboeuf   } else {
41309cf685aSAlexis Marboeuf     PetscCheck(!dm->useNatural, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "DM global to natural SF not present.\nIf DMPlexDistribute() was called and a section was defined, report to petsc-maint@mcs.anl.gov.");
41409cf685aSAlexis Marboeuf     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().");
41509cf685aSAlexis Marboeuf   }
41609cf685aSAlexis Marboeuf   PetscFunctionReturn(0);
41709cf685aSAlexis Marboeuf }
418