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