xref: /petsc/src/dm/impls/plex/plexnatural.c (revision 736995cdfdc1197aa072eb11621cf7b72dee6669)
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 
6f94b4a02SBlaise Bourdin + dm          - The DM
7f94b4a02SBlaise Bourdin . naturalSF   - The PetscSF
8f94b4a02SBlaise Bourdin   Level: intermediate
9f94b4a02SBlaise Bourdin 
10f94b4a02SBlaise Bourdin .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexCreateMigrationSF(), DMPlexGetMigrationSF()
11f94b4a02SBlaise Bourdin @*/
12f94b4a02SBlaise Bourdin PetscErrorCode DMPlexSetMigrationSF(DM dm, PetscSF migrationSF)
13f94b4a02SBlaise Bourdin {
14*736995cdSBlaise Bourdin   PetscErrorCode ierr;
15f94b4a02SBlaise Bourdin   PetscFunctionBegin;
16f94b4a02SBlaise Bourdin   dm->sfMigration = migrationSF;
17*736995cdSBlaise Bourdin   ierr = PetscObjectReference((PetscObject) migrationSF);CHKERRQ(ierr);
18f94b4a02SBlaise Bourdin   PetscFunctionReturn(0);
19f94b4a02SBlaise Bourdin }
20f94b4a02SBlaise Bourdin 
21f94b4a02SBlaise Bourdin /*@
22f94b4a02SBlaise Bourdin   DMPlexGetMigrationSF - Gets the SF for migrating from a parent DM into this DM
23f94b4a02SBlaise Bourdin 
24f94b4a02SBlaise Bourdin + dm          - The DM
25f94b4a02SBlaise Bourdin . *migrationSF   - The PetscSF
26f94b4a02SBlaise Bourdin   Level: intermediate
27f94b4a02SBlaise Bourdin 
28f94b4a02SBlaise Bourdin .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexCreateMigrationSF(), DMPlexSetMigrationSF
29f94b4a02SBlaise Bourdin @*/
30f94b4a02SBlaise Bourdin PetscErrorCode DMPlexGetMigrationSF(DM dm, PetscSF *migrationSF)
31f94b4a02SBlaise Bourdin {
32f94b4a02SBlaise Bourdin   PetscFunctionBegin;
33f94b4a02SBlaise Bourdin   *migrationSF = dm->sfMigration;
34f94b4a02SBlaise Bourdin   PetscFunctionReturn(0);
35f94b4a02SBlaise Bourdin }
36f94b4a02SBlaise Bourdin 
37f94b4a02SBlaise Bourdin /*@
38f94b4a02SBlaise Bourdin   DMPlexSetGlobalToNaturalSF - Sets the SF for mapping Global Vec to the Natural Vec
39f94b4a02SBlaise Bourdin 
40f94b4a02SBlaise Bourdin + dm          - The DM
41f94b4a02SBlaise Bourdin . naturalSF   - The PetscSF
42f94b4a02SBlaise Bourdin   Level: intermediate
43f94b4a02SBlaise Bourdin 
44f94b4a02SBlaise Bourdin .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexCreateGlobalToNaturalSF(), DMPlexGetGlobaltoNaturalSF()
45f94b4a02SBlaise Bourdin @*/
46f94b4a02SBlaise Bourdin PetscErrorCode DMPlexSetGlobalToNaturalSF(DM dm, PetscSF naturalSF)
47f94b4a02SBlaise Bourdin {
48*736995cdSBlaise Bourdin   PetscErrorCode ierr;
49f94b4a02SBlaise Bourdin   PetscFunctionBegin;
50f94b4a02SBlaise Bourdin   dm->sfNatural = naturalSF;
51*736995cdSBlaise Bourdin   ierr = PetscObjectReference((PetscObject) naturalSF);CHKERRQ(ierr);
52f94b4a02SBlaise Bourdin   dm->useNatural = PETSC_TRUE;
53f94b4a02SBlaise Bourdin   PetscFunctionReturn(0);
54f94b4a02SBlaise Bourdin }
55f94b4a02SBlaise Bourdin 
56f94b4a02SBlaise Bourdin /*@
57f94b4a02SBlaise Bourdin   DMPlexGetGlobalToNaturalSF - Gets the SF for mapping Global Vec to the Natural Vec
58f94b4a02SBlaise Bourdin 
59f94b4a02SBlaise Bourdin + dm          - The DM
60f94b4a02SBlaise Bourdin . *naturalSF   - The PetscSF
61f94b4a02SBlaise Bourdin   Level: intermediate
62f94b4a02SBlaise Bourdin 
63f94b4a02SBlaise Bourdin .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexCreateGlobalToNaturalSF(), DMPlexSetGlobaltoNaturalSF
64f94b4a02SBlaise Bourdin @*/
65f94b4a02SBlaise Bourdin PetscErrorCode DMPlexGetGlobalToNaturalSF(DM dm, PetscSF *naturalSF)
66f94b4a02SBlaise Bourdin {
67f94b4a02SBlaise Bourdin   PetscFunctionBegin;
68f94b4a02SBlaise Bourdin   *naturalSF = dm->sfNatural;
69f94b4a02SBlaise Bourdin   PetscFunctionReturn(0);
70f94b4a02SBlaise Bourdin }
71f94b4a02SBlaise Bourdin 
72f94b4a02SBlaise Bourdin /*@
73fa534816SMatthew G. Knepley   DMPlexCreateGlobalToNaturalSF - Creates the SF for mapping Global Vec to the Natural Vec
74fa534816SMatthew G. Knepley 
75fa534816SMatthew G. Knepley   Input Parameters:
76fa534816SMatthew G. Knepley + dm          - The DM
77fa534816SMatthew G. Knepley . section     - The PetscSection before the mesh was distributed
78fa534816SMatthew G. Knepley - sfMigration - The PetscSF used to distribute the mesh
79fa534816SMatthew G. Knepley 
80fa534816SMatthew G. Knepley   Output Parameters:
81fa534816SMatthew G. Knepley . sfNatural - PetscSF for mapping the Vec in PETSc ordering to the canonical ordering
82fa534816SMatthew G. Knepley 
83fa534816SMatthew G. Knepley   Level: intermediate
84fa534816SMatthew G. Knepley 
85fa534816SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField()
86fa534816SMatthew G. Knepley  @*/
87fa534816SMatthew G. Knepley PetscErrorCode DMPlexCreateGlobalToNaturalSF(DM dm, PetscSection section, PetscSF sfMigration, PetscSF *sfNatural)
88fa534816SMatthew G. Knepley {
89fa534816SMatthew G. Knepley   MPI_Comm       comm;
90fa534816SMatthew G. Knepley   Vec            gv;
91e5b44f4fSMatthew G. Knepley   PetscSF        sf, sfEmbed, sfSeqToNatural, sfField, sfFieldInv;
92e5b44f4fSMatthew G. Knepley   PetscSection   gSection, sectionDist, gLocSection;
93fa534816SMatthew G. Knepley   PetscInt      *spoints, *remoteOffsets;
94fa534816SMatthew G. Knepley   PetscInt       ssize, pStart, pEnd, p;
95fa534816SMatthew G. Knepley   PetscErrorCode ierr;
96fa534816SMatthew G. Knepley 
97fa534816SMatthew G. Knepley   PetscFunctionBegin;
98fa534816SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
99e5b44f4fSMatthew G. Knepley   /* ierr = PetscPrintf(comm, "Point migration SF\n");CHKERRQ(ierr);
100e5b44f4fSMatthew G. Knepley    ierr = PetscSFView(sfMigration, 0);CHKERRQ(ierr); */
101e5b44f4fSMatthew G. Knepley   /* Create a new section from distributing the original section */
102e5b44f4fSMatthew G. Knepley   ierr = PetscSectionCreate(comm, &sectionDist);CHKERRQ(ierr);
103e5b44f4fSMatthew G. Knepley   ierr = PetscSFDistributeSection(sfMigration, section, &remoteOffsets, sectionDist);CHKERRQ(ierr);
104e5b44f4fSMatthew G. Knepley   /* ierr = PetscPrintf(comm, "Distributed Section\n");CHKERRQ(ierr);
105e5b44f4fSMatthew G. Knepley    ierr = PetscSectionView(sectionDist, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */
106e5b44f4fSMatthew G. Knepley   ierr = DMSetDefaultSection(dm, sectionDist);CHKERRQ(ierr);
107fa534816SMatthew G. Knepley   /* Get a pruned version of migration SF */
108fa534816SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dm, &gSection);CHKERRQ(ierr);
109fa534816SMatthew G. Knepley   ierr = PetscSectionGetChart(gSection, &pStart, &pEnd);CHKERRQ(ierr);
110fa534816SMatthew G. Knepley   for (p = pStart, ssize = 0; p < pEnd; ++p) {
111fa534816SMatthew G. Knepley     PetscInt dof, off;
112fa534816SMatthew G. Knepley 
113fa534816SMatthew G. Knepley     ierr = PetscSectionGetDof(gSection, p, &dof);CHKERRQ(ierr);
114fa534816SMatthew G. Knepley     ierr = PetscSectionGetOffset(gSection, p, &off);CHKERRQ(ierr);
115fa534816SMatthew G. Knepley     if ((dof > 0) && (off >= 0)) ++ssize;
116fa534816SMatthew G. Knepley   }
117fa534816SMatthew G. Knepley   ierr = PetscMalloc1(ssize, &spoints);CHKERRQ(ierr);
118fa534816SMatthew G. Knepley   for (p = pStart, ssize = 0; p < pEnd; ++p) {
119fa534816SMatthew G. Knepley     PetscInt dof, off;
120fa534816SMatthew G. Knepley 
121fa534816SMatthew G. Knepley     ierr = PetscSectionGetDof(gSection, p, &dof);CHKERRQ(ierr);
122fa534816SMatthew G. Knepley     ierr = PetscSectionGetOffset(gSection, p, &off);CHKERRQ(ierr);
123fa534816SMatthew G. Knepley     if ((dof > 0) && (off >= 0)) spoints[ssize++] = p;
124fa534816SMatthew G. Knepley   }
125fa534816SMatthew G. Knepley   ierr = PetscSFCreateEmbeddedLeafSF(sfMigration, ssize, spoints, &sfEmbed);CHKERRQ(ierr);
126e5b44f4fSMatthew G. Knepley   ierr = PetscFree(spoints);CHKERRQ(ierr);
127e5b44f4fSMatthew G. Knepley   /* ierr = PetscPrintf(comm, "Embedded SF\n");CHKERRQ(ierr);
128e5b44f4fSMatthew G. Knepley    ierr = PetscSFView(sfEmbed, 0);CHKERRQ(ierr); */
129fa534816SMatthew G. Knepley   /* Create the SF for seq to natural */
130fa534816SMatthew G. Knepley   ierr = DMGetGlobalVector(dm, &gv);CHKERRQ(ierr);
131fa534816SMatthew G. Knepley   ierr = PetscSFCreateFromZero(comm, gv, &sfSeqToNatural);CHKERRQ(ierr);
132fa534816SMatthew G. Knepley   ierr = DMRestoreGlobalVector(dm, &gv);CHKERRQ(ierr);
133e5b44f4fSMatthew G. Knepley   /* ierr = PetscPrintf(comm, "Seq-to-Natural SF\n");CHKERRQ(ierr);
134e5b44f4fSMatthew G. Knepley    ierr = PetscSFView(sfSeqToNatural, 0);CHKERRQ(ierr); */
135fa534816SMatthew G. Knepley   /* Create the SF associated with this section */
136e5b44f4fSMatthew G. Knepley   ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
137e5b44f4fSMatthew G. Knepley   ierr = PetscSectionCreateGlobalSection(sectionDist, sf, PETSC_FALSE, PETSC_TRUE, &gLocSection);CHKERRQ(ierr);
138e5b44f4fSMatthew G. Knepley   ierr = PetscSFCreateSectionSF(sfEmbed, section, remoteOffsets, gLocSection, &sfField);CHKERRQ(ierr);
1390c374c54SMatthew G. Knepley   ierr = PetscFree(remoteOffsets);CHKERRQ(ierr);
140fa534816SMatthew G. Knepley   ierr = PetscSFDestroy(&sfEmbed);CHKERRQ(ierr);
141e5b44f4fSMatthew G. Knepley   ierr = PetscSectionDestroy(&gLocSection);CHKERRQ(ierr);
142e5b44f4fSMatthew G. Knepley   ierr = PetscSectionDestroy(&sectionDist);CHKERRQ(ierr);
143e5b44f4fSMatthew G. Knepley   /* ierr = PetscPrintf(comm, "Field SF\n");CHKERRQ(ierr);
144e5b44f4fSMatthew G. Knepley    ierr = PetscSFView(sfField, 0);CHKERRQ(ierr); */
145fa534816SMatthew G. Knepley   /* Invert the field SF so it's now from distributed to sequential */
146fa534816SMatthew G. Knepley   ierr = PetscSFCreateInverseSF(sfField, &sfFieldInv);CHKERRQ(ierr);
147fa534816SMatthew G. Knepley   ierr = PetscSFDestroy(&sfField);CHKERRQ(ierr);
148e5b44f4fSMatthew G. Knepley   /* ierr = PetscPrintf(comm, "Inverse Field SF\n");CHKERRQ(ierr);
149e5b44f4fSMatthew G. Knepley    ierr = PetscSFView(sfFieldInv, 0);CHKERRQ(ierr); */
150fa534816SMatthew G. Knepley   /* Multiply the sfFieldInv with the */
151fa534816SMatthew G. Knepley   ierr = PetscSFCompose(sfFieldInv, sfSeqToNatural, sfNatural);CHKERRQ(ierr);
152e5b44f4fSMatthew G. Knepley   ierr = PetscObjectViewFromOptions((PetscObject) *sfNatural, NULL, "-globaltonatural_sf_view");CHKERRQ(ierr);
153fa534816SMatthew G. Knepley   /* Clean up */
154fa534816SMatthew G. Knepley   ierr = PetscSFDestroy(&sfFieldInv);CHKERRQ(ierr);
155fa534816SMatthew G. Knepley   ierr = PetscSFDestroy(&sfSeqToNatural);CHKERRQ(ierr);
156fa534816SMatthew G. Knepley   PetscFunctionReturn(0);
157fa534816SMatthew G. Knepley }
158fa534816SMatthew G. Knepley 
159fa534816SMatthew G. Knepley /*@
160fa534816SMatthew G. Knepley   DMPlexGlobalToNaturalBegin - Rearranges a global Vector in the natural order.
161fa534816SMatthew G. Knepley 
162fa534816SMatthew G. Knepley   Collective on dm
163fa534816SMatthew G. Knepley 
164fa534816SMatthew G. Knepley   Input Parameters:
165fa534816SMatthew G. Knepley + dm - The distributed DMPlex
166fa534816SMatthew G. Knepley - gv - The global Vec
167fa534816SMatthew G. Knepley 
168fa534816SMatthew G. Knepley   Output Parameters:
169fa534816SMatthew G. Knepley . nv - Vec in the canonical ordering distributed over all processors associated with gv
170fa534816SMatthew G. Knepley 
171fa534816SMatthew G. Knepley   Note: The user must call DMPlexSetUseNaturalSF(dm, PETSC_TRUE) before DMPlexDistribute().
172fa534816SMatthew G. Knepley 
173fa534816SMatthew G. Knepley   Level: intermediate
174fa534816SMatthew G. Knepley 
175fa534816SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexNaturalToGlobalBegin(), DMPlexGlobalToNaturalEnd()
176fa534816SMatthew G. Knepley @*/
177fa534816SMatthew G. Knepley PetscErrorCode DMPlexGlobalToNaturalBegin(DM dm, Vec gv, Vec nv)
178fa534816SMatthew G. Knepley {
179fa534816SMatthew G. Knepley   const PetscScalar *inarray;
180fa534816SMatthew G. Knepley   PetscScalar       *outarray;
181fa534816SMatthew G. Knepley   PetscErrorCode     ierr;
182fa534816SMatthew G. Knepley 
183fa534816SMatthew G. Knepley   PetscFunctionBegin;
184fa534816SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_GlobalToNaturalBegin,dm,0,0,0);CHKERRQ(ierr);
185fa534816SMatthew G. Knepley   if (dm->sfNatural) {
186fa534816SMatthew G. Knepley     ierr = VecGetArray(nv, &outarray);CHKERRQ(ierr);
187fa534816SMatthew G. Knepley     ierr = VecGetArrayRead(gv, &inarray);CHKERRQ(ierr);
188fa534816SMatthew G. Knepley     ierr = PetscSFBcastBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray);CHKERRQ(ierr);
189fa534816SMatthew G. Knepley     ierr = VecRestoreArrayRead(gv, &inarray);CHKERRQ(ierr);
190fa534816SMatthew G. Knepley     ierr = VecRestoreArray(nv, &outarray);CHKERRQ(ierr);
191e5b44f4fSMatthew G. Knepley   } else SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().\n");
192fa534816SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_GlobalToNaturalBegin,dm,0,0,0);CHKERRQ(ierr);
193fa534816SMatthew G. Knepley   PetscFunctionReturn(0);
194fa534816SMatthew G. Knepley }
195fa534816SMatthew G. Knepley 
196fa534816SMatthew G. Knepley /*@
197fa534816SMatthew G. Knepley   DMPlexGlobalToNaturalEnd - Rearranges a global Vector in the natural order.
198fa534816SMatthew G. Knepley 
199fa534816SMatthew G. Knepley   Collective on dm
200fa534816SMatthew G. Knepley 
201fa534816SMatthew G. Knepley   Input Parameters:
202fa534816SMatthew G. Knepley + dm - The distributed DMPlex
203fa534816SMatthew G. Knepley - gv - The global Vec
204fa534816SMatthew G. Knepley 
205fa534816SMatthew G. Knepley   Output Parameters:
206fa534816SMatthew G. Knepley . nv - The natural Vec
207fa534816SMatthew G. Knepley 
208fa534816SMatthew G. Knepley   Note: The user must call DMPlexSetUseNaturalSF(dm, PETSC_TRUE) before DMPlexDistribute().
209fa534816SMatthew G. Knepley 
210fa534816SMatthew G. Knepley   Level: intermediate
211fa534816SMatthew G. Knepley 
212fa534816SMatthew G. Knepley  .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexNaturalToGlobalBegin(), DMPlexGlobalToNaturalBegin()
213fa534816SMatthew G. Knepley  @*/
214fa534816SMatthew G. Knepley PetscErrorCode DMPlexGlobalToNaturalEnd(DM dm, Vec gv, Vec nv)
215fa534816SMatthew G. Knepley {
216fa534816SMatthew G. Knepley   const PetscScalar *inarray;
217fa534816SMatthew G. Knepley   PetscScalar       *outarray;
218fa534816SMatthew G. Knepley   PetscErrorCode     ierr;
219fa534816SMatthew G. Knepley 
220fa534816SMatthew G. Knepley   PetscFunctionBegin;
221fa534816SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_GlobalToNaturalEnd,dm,0,0,0);CHKERRQ(ierr);
222fa534816SMatthew G. Knepley   if (dm->sfNatural) {
223fa534816SMatthew G. Knepley     ierr = VecGetArrayRead(gv, &inarray);CHKERRQ(ierr);
224fa534816SMatthew G. Knepley     ierr = VecGetArray(nv, &outarray);CHKERRQ(ierr);CHKERRQ(ierr);
225fa534816SMatthew G. Knepley     ierr = PetscSFBcastEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray);CHKERRQ(ierr);
226fa534816SMatthew G. Knepley     ierr = VecRestoreArrayRead(gv, &inarray);CHKERRQ(ierr);
227fa534816SMatthew G. Knepley     ierr = VecRestoreArray(nv, &outarray);CHKERRQ(ierr);
228fa534816SMatthew G. Knepley   }
229fa534816SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_GlobalToNaturalEnd,dm,0,0,0);CHKERRQ(ierr);
230fa534816SMatthew G. Knepley   PetscFunctionReturn(0);
231fa534816SMatthew G. Knepley }
232fa534816SMatthew G. Knepley 
233fa534816SMatthew G. Knepley /*@
234fa534816SMatthew G. Knepley   DMPlexNaturalToGlobalBegin - Rearranges a Vector in the natural order to the Global order.
235fa534816SMatthew G. Knepley 
236fa534816SMatthew G. Knepley   Collective on dm
237fa534816SMatthew G. Knepley 
238fa534816SMatthew G. Knepley   Input Parameters:
239fa534816SMatthew G. Knepley + dm - The distributed DMPlex
240fa534816SMatthew G. Knepley - nv - The natural Vec
241fa534816SMatthew G. Knepley 
242fa534816SMatthew G. Knepley   Output Parameters:
243fa534816SMatthew G. Knepley . gv - The global Vec
244fa534816SMatthew G. Knepley 
245fa534816SMatthew G. Knepley   Note: The user must call DMPlexSetUseNaturalSF(dm, PETSC_TRUE) before DMPlexDistribute().
246fa534816SMatthew G. Knepley 
247fa534816SMatthew G. Knepley   Level: intermediate
248fa534816SMatthew G. Knepley 
249fa534816SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexNaturalToGlobalBegin(),DMPlexGlobalToNaturalEnd()
250fa534816SMatthew G. Knepley @*/
251fa534816SMatthew G. Knepley PetscErrorCode DMPlexNaturalToGlobalBegin(DM dm, Vec nv, Vec gv)
252fa534816SMatthew G. Knepley {
253fa534816SMatthew G. Knepley   const PetscScalar *inarray;
254fa534816SMatthew G. Knepley   PetscScalar       *outarray;
255fa534816SMatthew G. Knepley   PetscErrorCode     ierr;
256fa534816SMatthew G. Knepley 
257fa534816SMatthew G. Knepley   PetscFunctionBegin;
258fa534816SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_NaturalToGlobalBegin,dm,0,0,0);CHKERRQ(ierr);
259fa534816SMatthew G. Knepley   if (dm->sfNatural) {
260b64f75a9SBlaise Bourdin     /* We only have acces to the SF that goes from Global to Natural.
261b64f75a9SBlaise Bourdin        Instead of inverting dm->sfNatural, we can call PetscSFReduceBegin/End with MPI_Op MPI_SUM.
262b64f75a9SBlaise Bourdin        Here the SUM really does nothing since sfNatural is one to one, as long as gV is set to zero first. */
263b64f75a9SBlaise Bourdin     ierr = VecZeroEntries(gv);CHKERRQ(ierr);
264fa534816SMatthew G. Knepley     ierr = VecGetArray(gv, &outarray);CHKERRQ(ierr);
265fa534816SMatthew G. Knepley     ierr = VecGetArrayRead(nv, &inarray);CHKERRQ(ierr);
266fa534816SMatthew G. Knepley     ierr = PetscSFReduceBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray, MPI_SUM);CHKERRQ(ierr);
267fa534816SMatthew G. Knepley     ierr = VecRestoreArrayRead(nv, &inarray);CHKERRQ(ierr);
268fa534816SMatthew G. Knepley     ierr = VecRestoreArray(gv, &outarray);CHKERRQ(ierr);
269fa534816SMatthew G. Knepley   } else SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMPlexSetUseNaturalSF() before DMPlexDistribute().\n");
270fa534816SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_NaturalToGlobalBegin,dm,0,0,0);CHKERRQ(ierr);
271fa534816SMatthew G. Knepley   PetscFunctionReturn(0);
272fa534816SMatthew G. Knepley }
273fa534816SMatthew G. Knepley 
274fa534816SMatthew G. Knepley /*@
275fa534816SMatthew G. Knepley   DMPlexNaturalToGlobalEnd - Rearranges a Vector in the natural order to the Global order.
276fa534816SMatthew G. Knepley 
277fa534816SMatthew G. Knepley   Collective on dm
278fa534816SMatthew G. Knepley 
279fa534816SMatthew G. Knepley   Input Parameters:
280fa534816SMatthew G. Knepley + dm - The distributed DMPlex
281fa534816SMatthew G. Knepley - nv - The natural Vec
282fa534816SMatthew G. Knepley 
283fa534816SMatthew G. Knepley   Output Parameters:
284fa534816SMatthew G. Knepley . gv - The global Vec
285fa534816SMatthew G. Knepley 
286fa534816SMatthew G. Knepley   Note: The user must call DMPlexSetUseNaturalSF(dm, PETSC_TRUE) before DMPlexDistribute().
287fa534816SMatthew G. Knepley 
288fa534816SMatthew G. Knepley   Level: intermediate
289fa534816SMatthew G. Knepley 
290fa534816SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexNaturalToGlobalBegin(), DMPlexGlobalToNaturalBegin()
291fa534816SMatthew G. Knepley  @*/
292fa534816SMatthew G. Knepley PetscErrorCode DMPlexNaturalToGlobalEnd(DM dm, Vec nv, Vec gv)
293fa534816SMatthew G. Knepley {
294fa534816SMatthew G. Knepley   const PetscScalar *inarray;
295fa534816SMatthew G. Knepley   PetscScalar       *outarray;
296fa534816SMatthew G. Knepley   PetscErrorCode     ierr;
297fa534816SMatthew G. Knepley 
298fa534816SMatthew G. Knepley   PetscFunctionBegin;
299fa534816SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_NaturalToGlobalEnd,dm,0,0,0);CHKERRQ(ierr);
300fa534816SMatthew G. Knepley   if (dm->sfNatural) {
301fa534816SMatthew G. Knepley     ierr = VecGetArrayRead(nv, &inarray);CHKERRQ(ierr);
302fa534816SMatthew G. Knepley     ierr = VecGetArray(gv, &outarray);CHKERRQ(ierr);CHKERRQ(ierr);
303fa534816SMatthew G. Knepley     ierr = PetscSFReduceEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray, MPI_SUM);CHKERRQ(ierr);
304fa534816SMatthew G. Knepley     ierr = VecRestoreArrayRead(nv, &inarray);CHKERRQ(ierr);
305fa534816SMatthew G. Knepley     ierr = VecRestoreArray(gv, &outarray);CHKERRQ(ierr);
306fa534816SMatthew G. Knepley   }
307fa534816SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_NaturalToGlobalEnd,dm,0,0,0);CHKERRQ(ierr);
308fa534816SMatthew G. Knepley   PetscFunctionReturn(0);
309fa534816SMatthew G. Knepley }
310