xref: /petsc/src/dm/impls/plex/plexnatural.c (revision e78587018951896e4d4c45b50c69c7c99ccd87c0)
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 @*/
16f94b4a02SBlaise Bourdin PetscErrorCode DMPlexSetMigrationSF(DM dm, PetscSF migrationSF)
17f94b4a02SBlaise Bourdin {
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 @*/
37f94b4a02SBlaise Bourdin PetscErrorCode DMPlexGetMigrationSF(DM dm, PetscSF *migrationSF)
38f94b4a02SBlaise Bourdin {
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 @*/
55f94b4a02SBlaise Bourdin PetscErrorCode DMPlexSetGlobalToNaturalSF(DM dm, PetscSF naturalSF)
56f94b4a02SBlaise Bourdin {
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 @*/
77f94b4a02SBlaise Bourdin PetscErrorCode DMPlexGetGlobalToNaturalSF(DM dm, PetscSF *naturalSF)
78f94b4a02SBlaise Bourdin {
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:
88*e7858701SMatthew G. Knepley + dm          - The redistributed DM
89*e7858701SMatthew 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  @*/
101fa534816SMatthew G. Knepley PetscErrorCode DMPlexCreateGlobalToNaturalSF(DM dm, PetscSection section, PetscSF sfMigration, PetscSF *sfNatural)
102fa534816SMatthew G. Knepley {
103fa534816SMatthew G. Knepley   MPI_Comm     comm;
104*e7858701SMatthew G. Knepley   Vec          tmpVec;
105*e7858701SMatthew G. Knepley   PetscSF      sf, sfEmbed, sfField;
106e5b44f4fSMatthew G. Knepley   PetscSection gSection, sectionDist, gLocSection;
107fa534816SMatthew G. Knepley   PetscInt    *spoints, *remoteOffsets;
1085c3f5608SAlexisMarb   PetscInt     ssize, pStart, pEnd, p, globalSize;
109*e7858701SMatthew G. Knepley   PetscBool    destroyFlag = PETSC_FALSE, debug = PETSC_FALSE;
110fa534816SMatthew G. Knepley 
111fa534816SMatthew G. Knepley   PetscFunctionBegin;
1129566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
1135c3f5608SAlexisMarb   if (!sfMigration) {
114*e7858701SMatthew G. Knepley     /* If sfMigration is missing, sfNatural cannot be computed and is set to NULL */
1155c3f5608SAlexisMarb     *sfNatural = NULL;
1165c3f5608SAlexisMarb     PetscFunctionReturn(0);
1175c3f5608SAlexisMarb   } else if (!section) {
118*e7858701SMatthew G. Knepley     /* If the sequential section is not provided (NULL), it is reconstructed from the parallel section */
1195c3f5608SAlexisMarb     PetscSF sfMigrationInv;
1205c3f5608SAlexisMarb     PetscSection localSection;
1215c3f5608SAlexisMarb 
1229566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dm, &localSection));
1239566063dSJacob Faibussowitsch     PetscCall(PetscSFCreateInverseSF(sfMigration, &sfMigrationInv));
1249566063dSJacob Faibussowitsch     PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject) dm), &section));
1259566063dSJacob Faibussowitsch     PetscCall(PetscSFDistributeSection(sfMigrationInv, localSection, NULL, section));
1269566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&sfMigrationInv));
1275c3f5608SAlexisMarb     destroyFlag = PETSC_TRUE;
1285c3f5608SAlexisMarb   }
129*e7858701SMatthew G. Knepley   if (debug) PetscCall(PetscSFView(sfMigration, NULL));
130e5b44f4fSMatthew G. Knepley   /* Create a new section from distributing the original section */
1319566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(comm, &sectionDist));
1329566063dSJacob Faibussowitsch   PetscCall(PetscSFDistributeSection(sfMigration, section, &remoteOffsets, sectionDist));
133*e7858701SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject) sectionDist, "Migrated Section"));
134*e7858701SMatthew G. Knepley   if (debug) PetscCall(PetscSectionView(sectionDist, NULL));
1359566063dSJacob Faibussowitsch   PetscCall(DMSetLocalSection(dm, sectionDist));
136*e7858701SMatthew G. Knepley   /* If a sequential section is provided but no dof is affected, sfNatural cannot be computed and is set to NULL */
1379566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &tmpVec));
1389566063dSJacob Faibussowitsch   PetscCall(VecGetSize(tmpVec, &globalSize));
1399566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dm, &tmpVec));
14014744f87SBlaise Bourdin   if (globalSize) {
141*e7858701SMatthew G. Knepley     const PetscInt *leaves;
142*e7858701SMatthew G. Knepley     PetscInt       *sortleaves, *indices;
143*e7858701SMatthew G. Knepley     PetscInt        Nl;
144*e7858701SMatthew G. Knepley 
145fa534816SMatthew G. Knepley     /* Get a pruned version of migration SF */
1469566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalSection(dm, &gSection));
147*e7858701SMatthew G. Knepley     if (debug) PetscCall(PetscSectionView(gSection, NULL));
148*e7858701SMatthew G. Knepley     PetscCall(PetscSFGetGraph(sfMigration, NULL, &Nl, &leaves, NULL));
1499566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetChart(gSection, &pStart, &pEnd));
150fa534816SMatthew G. Knepley     for (p = pStart, ssize = 0; p < pEnd; ++p) {
151fa534816SMatthew G. Knepley       PetscInt dof, off;
152fa534816SMatthew G. Knepley 
1539566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(gSection, p, &dof));
1549566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(gSection, p, &off));
155fa534816SMatthew G. Knepley       if ((dof > 0) && (off >= 0)) ++ssize;
156fa534816SMatthew G. Knepley     }
157*e7858701SMatthew G. Knepley     PetscCall(PetscMalloc3(ssize, &spoints, Nl, &sortleaves, Nl, &indices));
158*e7858701SMatthew G. Knepley     for (p = 0; p < Nl; ++p) {sortleaves[p] = leaves ? leaves[p] : p; indices[p] = p;}
159*e7858701SMatthew G. Knepley     PetscCall(PetscSortIntWithArray(Nl, sortleaves, indices));
160fa534816SMatthew G. Knepley     for (p = pStart, ssize = 0; p < pEnd; ++p) {
161*e7858701SMatthew G. Knepley       PetscInt dof, off, loc;
162fa534816SMatthew G. Knepley 
1639566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(gSection, p, &dof));
1649566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(gSection, p, &off));
165*e7858701SMatthew G. Knepley       if ((dof > 0) && (off >= 0)) {
166*e7858701SMatthew G. Knepley         PetscCall(PetscFindInt(p, Nl, sortleaves, &loc));
167*e7858701SMatthew 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);
168*e7858701SMatthew G. Knepley         spoints[ssize++] = indices[loc];
169*e7858701SMatthew G. Knepley       }
170fa534816SMatthew G. Knepley     }
1719566063dSJacob Faibussowitsch     PetscCall(PetscSFCreateEmbeddedLeafSF(sfMigration, ssize, spoints, &sfEmbed));
172*e7858701SMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject) sfEmbed, "Embedded SF"));
173*e7858701SMatthew G. Knepley     PetscCall(PetscFree3(spoints, sortleaves, indices));
174*e7858701SMatthew G. Knepley     if (debug) PetscCall(PetscSFView(sfEmbed, NULL));
175*e7858701SMatthew G. Knepley     /* Create the SF associated with this section
176*e7858701SMatthew G. Knepley          Roots are natural dofs, leaves are global dofs */
1779566063dSJacob Faibussowitsch     PetscCall(DMGetPointSF(dm, &sf));
1789566063dSJacob Faibussowitsch     PetscCall(PetscSectionCreateGlobalSection(sectionDist, sf, PETSC_FALSE, PETSC_TRUE, &gLocSection));
1799566063dSJacob Faibussowitsch     PetscCall(PetscSFCreateSectionSF(sfEmbed, section, remoteOffsets, gLocSection, &sfField));
1809566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&sfEmbed));
1819566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&gLocSection));
182*e7858701SMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject) sfField, "Natural-to-Global SF"));
183*e7858701SMatthew G. Knepley     if (debug) PetscCall(PetscSFView(sfField, NULL));
184*e7858701SMatthew G. Knepley     /* Invert the field SF
185*e7858701SMatthew G. Knepley          Roots are global dofs, leaves are natural dofs */
186*e7858701SMatthew G. Knepley     PetscCall(PetscSFCreateInverseSF(sfField, sfNatural));
187*e7858701SMatthew 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 */
190*e7858701SMatthew 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 @*/
218fa534816SMatthew G. Knepley PetscErrorCode DMPlexGlobalToNaturalBegin(DM dm, Vec gv, Vec nv)
219fa534816SMatthew G. Knepley {
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));
23528b400f6SJacob Faibussowitsch   } else 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.");
236546078acSJacob Faibussowitsch   else SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().");
2379566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_GlobalToNaturalBegin,dm,0,0,0));
238fa534816SMatthew G. Knepley   PetscFunctionReturn(0);
239fa534816SMatthew G. Knepley }
240fa534816SMatthew G. Knepley 
241fa534816SMatthew G. Knepley /*@
242fa534816SMatthew G. Knepley   DMPlexGlobalToNaturalEnd - Rearranges a global Vector in the natural order.
243fa534816SMatthew G. Knepley 
244fa534816SMatthew G. Knepley   Collective on dm
245fa534816SMatthew G. Knepley 
246fa534816SMatthew G. Knepley   Input Parameters:
247fa534816SMatthew G. Knepley + dm - The distributed DMPlex
248fa534816SMatthew G. Knepley - gv - The global Vec
249fa534816SMatthew G. Knepley 
25097bb3fdcSJose E. Roman   Output Parameter:
251fa534816SMatthew G. Knepley . nv - The natural Vec
252fa534816SMatthew G. Knepley 
25342ea106eSTristan Konolige   Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute().
254fa534816SMatthew G. Knepley 
255fa534816SMatthew G. Knepley   Level: intermediate
256fa534816SMatthew G. Knepley 
257db781477SPatrick Sanan  .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()`
258fa534816SMatthew G. Knepley  @*/
259fa534816SMatthew G. Knepley PetscErrorCode DMPlexGlobalToNaturalEnd(DM dm, Vec gv, Vec nv)
260fa534816SMatthew G. Knepley {
261fa534816SMatthew G. Knepley   const PetscScalar *inarray;
262fa534816SMatthew G. Knepley   PetscScalar       *outarray;
263a6a55facSMatthew G. Knepley   PetscMPIInt        size;
264fa534816SMatthew G. Knepley 
265fa534816SMatthew G. Knepley   PetscFunctionBegin;
2669566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_GlobalToNaturalEnd,dm,0,0,0));
2679566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size));
268fa534816SMatthew G. Knepley   if (dm->sfNatural) {
2699566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(gv, &inarray));
2709566063dSJacob Faibussowitsch     PetscCall(VecGetArray(nv, &outarray));
2719566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray,MPI_REPLACE));
2729566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(gv, &inarray));
2739566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(nv, &outarray));
274a6a55facSMatthew G. Knepley   } else if (size == 1) {
27528b400f6SJacob Faibussowitsch   } else 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.");
276546078acSJacob Faibussowitsch   else SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().");
2779566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_GlobalToNaturalEnd,dm,0,0,0));
278fa534816SMatthew G. Knepley   PetscFunctionReturn(0);
279fa534816SMatthew G. Knepley }
280fa534816SMatthew G. Knepley 
281fa534816SMatthew G. Knepley /*@
282fa534816SMatthew G. Knepley   DMPlexNaturalToGlobalBegin - Rearranges a Vector in the natural order to the Global order.
283fa534816SMatthew G. Knepley 
284fa534816SMatthew G. Knepley   Collective on dm
285fa534816SMatthew G. Knepley 
286fa534816SMatthew G. Knepley   Input Parameters:
287fa534816SMatthew G. Knepley + dm - The distributed DMPlex
288fa534816SMatthew G. Knepley - nv - The natural Vec
289fa534816SMatthew G. Knepley 
290fa534816SMatthew G. Knepley   Output Parameters:
291fa534816SMatthew G. Knepley . gv - The global Vec
292fa534816SMatthew G. Knepley 
29342ea106eSTristan Konolige   Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute().
294fa534816SMatthew G. Knepley 
295fa534816SMatthew G. Knepley   Level: intermediate
296fa534816SMatthew G. Knepley 
297c2e3fba1SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalEnd()`
298fa534816SMatthew G. Knepley @*/
299fa534816SMatthew G. Knepley PetscErrorCode DMPlexNaturalToGlobalBegin(DM dm, Vec nv, Vec gv)
300fa534816SMatthew G. Knepley {
301fa534816SMatthew G. Knepley   const PetscScalar *inarray;
302fa534816SMatthew G. Knepley   PetscScalar       *outarray;
303a6a55facSMatthew G. Knepley   PetscMPIInt        size;
304fa534816SMatthew G. Knepley 
305fa534816SMatthew G. Knepley   PetscFunctionBegin;
3069566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_NaturalToGlobalBegin,dm,0,0,0));
3079566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size));
308fa534816SMatthew G. Knepley   if (dm->sfNatural) {
309a5b23f4aSJose E. Roman     /* We only have access to the SF that goes from Global to Natural.
310b64f75a9SBlaise Bourdin        Instead of inverting dm->sfNatural, we can call PetscSFReduceBegin/End with MPI_Op MPI_SUM.
311b64f75a9SBlaise Bourdin        Here the SUM really does nothing since sfNatural is one to one, as long as gV is set to zero first. */
3129566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(gv));
3139566063dSJacob Faibussowitsch     PetscCall(VecGetArray(gv, &outarray));
3149566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(nv, &inarray));
3159566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray, MPI_SUM));
3169566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(nv, &inarray));
3179566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(gv, &outarray));
318a6a55facSMatthew G. Knepley   } else if (size == 1) {
3199566063dSJacob Faibussowitsch     PetscCall(VecCopy(nv, gv));
32028b400f6SJacob Faibussowitsch   } else 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.");
321546078acSJacob Faibussowitsch   else SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().");
3229566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_NaturalToGlobalBegin,dm,0,0,0));
323fa534816SMatthew G. Knepley   PetscFunctionReturn(0);
324fa534816SMatthew G. Knepley }
325fa534816SMatthew G. Knepley 
326fa534816SMatthew G. Knepley /*@
327fa534816SMatthew G. Knepley   DMPlexNaturalToGlobalEnd - Rearranges a Vector in the natural order to the Global order.
328fa534816SMatthew G. Knepley 
329fa534816SMatthew G. Knepley   Collective on dm
330fa534816SMatthew G. Knepley 
331fa534816SMatthew G. Knepley   Input Parameters:
332fa534816SMatthew G. Knepley + dm - The distributed DMPlex
333fa534816SMatthew G. Knepley - nv - The natural Vec
334fa534816SMatthew G. Knepley 
335fa534816SMatthew G. Knepley   Output Parameters:
336fa534816SMatthew G. Knepley . gv - The global Vec
337fa534816SMatthew G. Knepley 
33842ea106eSTristan Konolige   Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute().
339fa534816SMatthew G. Knepley 
340fa534816SMatthew G. Knepley   Level: intermediate
341fa534816SMatthew G. Knepley 
342db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()`
343fa534816SMatthew G. Knepley  @*/
344fa534816SMatthew G. Knepley PetscErrorCode DMPlexNaturalToGlobalEnd(DM dm, Vec nv, Vec gv)
345fa534816SMatthew G. Knepley {
346fa534816SMatthew G. Knepley   const PetscScalar *inarray;
347fa534816SMatthew G. Knepley   PetscScalar       *outarray;
348a6a55facSMatthew G. Knepley   PetscMPIInt        size;
349fa534816SMatthew G. Knepley 
350fa534816SMatthew G. Knepley   PetscFunctionBegin;
3519566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_NaturalToGlobalEnd,dm,0,0,0));
3529566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size));
353fa534816SMatthew G. Knepley   if (dm->sfNatural) {
3549566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(nv, &inarray));
3559566063dSJacob Faibussowitsch     PetscCall(VecGetArray(gv, &outarray));
3569566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray, MPI_SUM));
3579566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(nv, &inarray));
3589566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(gv, &outarray));
359a6a55facSMatthew G. Knepley   } else if (size == 1) {
36028b400f6SJacob Faibussowitsch   } else 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.");
361546078acSJacob Faibussowitsch   else SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().");
3629566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_NaturalToGlobalEnd,dm,0,0,0));
363fa534816SMatthew G. Knepley   PetscFunctionReturn(0);
364fa534816SMatthew G. Knepley }
365