xref: /petsc/src/dm/impls/plex/plexnatural.c (revision 5c3f560880f646f22ca964863138187637ebcb2e)
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 
14f94b4a02SBlaise Bourdin .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexCreateMigrationSF(), DMPlexGetMigrationSF()
15f94b4a02SBlaise Bourdin @*/
16f94b4a02SBlaise Bourdin PetscErrorCode DMPlexSetMigrationSF(DM dm, PetscSF migrationSF)
17f94b4a02SBlaise Bourdin {
18736995cdSBlaise Bourdin   PetscErrorCode ierr;
19f94b4a02SBlaise Bourdin   PetscFunctionBegin;
20f94b4a02SBlaise Bourdin   dm->sfMigration = migrationSF;
21736995cdSBlaise Bourdin   ierr = PetscObjectReference((PetscObject) migrationSF);CHKERRQ(ierr);
22f94b4a02SBlaise Bourdin   PetscFunctionReturn(0);
23f94b4a02SBlaise Bourdin }
24f94b4a02SBlaise Bourdin 
25f94b4a02SBlaise Bourdin /*@
26f94b4a02SBlaise Bourdin   DMPlexGetMigrationSF - Gets the SF for migrating from a parent DM into this DM
27f94b4a02SBlaise Bourdin 
285d3b26e6SMatthew G. Knepley   Input Parameter:
295d3b26e6SMatthew G. Knepley . dm          - The DM
305d3b26e6SMatthew G. Knepley 
315d3b26e6SMatthew G. Knepley   Output Parameter:
325d3b26e6SMatthew G. Knepley . migrationSF - The PetscSF
335d3b26e6SMatthew G. Knepley 
34f94b4a02SBlaise Bourdin   Level: intermediate
35f94b4a02SBlaise Bourdin 
36f94b4a02SBlaise Bourdin .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexCreateMigrationSF(), DMPlexSetMigrationSF
37f94b4a02SBlaise Bourdin @*/
38f94b4a02SBlaise Bourdin PetscErrorCode DMPlexGetMigrationSF(DM dm, PetscSF *migrationSF)
39f94b4a02SBlaise Bourdin {
40f94b4a02SBlaise Bourdin   PetscFunctionBegin;
41f94b4a02SBlaise Bourdin   *migrationSF = dm->sfMigration;
42f94b4a02SBlaise Bourdin   PetscFunctionReturn(0);
43f94b4a02SBlaise Bourdin }
44f94b4a02SBlaise Bourdin 
45f94b4a02SBlaise Bourdin /*@
46f94b4a02SBlaise Bourdin   DMPlexSetGlobalToNaturalSF - Sets the SF for mapping Global Vec to the Natural Vec
47f94b4a02SBlaise Bourdin 
485d3b26e6SMatthew G. Knepley   Input Parameters:
49f94b4a02SBlaise Bourdin + dm          - The DM
505d3b26e6SMatthew G. Knepley - naturalSF   - The PetscSF
515d3b26e6SMatthew G. Knepley 
52f94b4a02SBlaise Bourdin   Level: intermediate
53f94b4a02SBlaise Bourdin 
54f94b4a02SBlaise Bourdin .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexCreateGlobalToNaturalSF(), DMPlexGetGlobaltoNaturalSF()
55f94b4a02SBlaise Bourdin @*/
56f94b4a02SBlaise Bourdin PetscErrorCode DMPlexSetGlobalToNaturalSF(DM dm, PetscSF naturalSF)
57f94b4a02SBlaise Bourdin {
58736995cdSBlaise Bourdin   PetscErrorCode ierr;
59f94b4a02SBlaise Bourdin   PetscFunctionBegin;
60f94b4a02SBlaise Bourdin   dm->sfNatural = naturalSF;
61736995cdSBlaise Bourdin   ierr = PetscObjectReference((PetscObject) naturalSF);CHKERRQ(ierr);
62f94b4a02SBlaise Bourdin   dm->useNatural = PETSC_TRUE;
63f94b4a02SBlaise Bourdin   PetscFunctionReturn(0);
64f94b4a02SBlaise Bourdin }
65f94b4a02SBlaise Bourdin 
66f94b4a02SBlaise Bourdin /*@
67f94b4a02SBlaise Bourdin   DMPlexGetGlobalToNaturalSF - Gets the SF for mapping Global Vec to the Natural Vec
68f94b4a02SBlaise Bourdin 
695d3b26e6SMatthew G. Knepley   Input Parameter:
705d3b26e6SMatthew G. Knepley . dm          - The DM
715d3b26e6SMatthew G. Knepley 
725d3b26e6SMatthew G. Knepley   Output Parameter:
735d3b26e6SMatthew G. Knepley . naturalSF   - The PetscSF
745d3b26e6SMatthew G. Knepley 
75f94b4a02SBlaise Bourdin   Level: intermediate
76f94b4a02SBlaise Bourdin 
77f94b4a02SBlaise Bourdin .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexCreateGlobalToNaturalSF(), DMPlexSetGlobaltoNaturalSF
78f94b4a02SBlaise Bourdin @*/
79f94b4a02SBlaise Bourdin PetscErrorCode DMPlexGetGlobalToNaturalSF(DM dm, PetscSF *naturalSF)
80f94b4a02SBlaise Bourdin {
81f94b4a02SBlaise Bourdin   PetscFunctionBegin;
82f94b4a02SBlaise Bourdin   *naturalSF = dm->sfNatural;
83f94b4a02SBlaise Bourdin   PetscFunctionReturn(0);
84f94b4a02SBlaise Bourdin }
85f94b4a02SBlaise Bourdin 
86f94b4a02SBlaise Bourdin /*@
87fa534816SMatthew G. Knepley   DMPlexCreateGlobalToNaturalSF - Creates the SF for mapping Global Vec to the Natural Vec
88fa534816SMatthew G. Knepley 
89fa534816SMatthew G. Knepley   Input Parameters:
90fa534816SMatthew G. Knepley + dm          - The DM
91*5c3f5608SAlexisMarb . section     - The PetscSection describing the Vec before the mesh was distributed,
92*5c3f5608SAlexisMarb                 or NULL if not available
93*5c3f5608SAlexisMarb - sfMigration - The PetscSF used to distribute the mesh, or NULL if it cannot be computed
94fa534816SMatthew G. Knepley 
955d3b26e6SMatthew G. Knepley   Output Parameter:
96fa534816SMatthew G. Knepley . sfNatural   - PetscSF for mapping the Vec in PETSc ordering to the canonical ordering
97fa534816SMatthew G. Knepley 
985d3b26e6SMatthew G. Knepley   Note: This is not typically called by the user.
995d3b26e6SMatthew G. Knepley 
100fa534816SMatthew G. Knepley   Level: intermediate
101fa534816SMatthew G. Knepley 
102fa534816SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField()
103fa534816SMatthew G. Knepley  @*/
104fa534816SMatthew G. Knepley PetscErrorCode DMPlexCreateGlobalToNaturalSF(DM dm, PetscSection section, PetscSF sfMigration, PetscSF *sfNatural)
105fa534816SMatthew G. Knepley {
106fa534816SMatthew G. Knepley   MPI_Comm       comm;
107*5c3f5608SAlexisMarb   Vec            gv, tmpVec;
108e5b44f4fSMatthew G. Knepley   PetscSF        sf, sfEmbed, sfSeqToNatural, sfField, sfFieldInv;
109e5b44f4fSMatthew G. Knepley   PetscSection   gSection, sectionDist, gLocSection;
110fa534816SMatthew G. Knepley   PetscInt      *spoints, *remoteOffsets;
111*5c3f5608SAlexisMarb   PetscInt       ssize, pStart, pEnd, p, globalSize;
1123a350576SJunchao Zhang   PetscLayout    map;
113*5c3f5608SAlexisMarb   PetscBool      destroyFlag = PETSC_FALSE;
114fa534816SMatthew G. Knepley   PetscErrorCode ierr;
115fa534816SMatthew G. Knepley 
116fa534816SMatthew G. Knepley   PetscFunctionBegin;
117fa534816SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
118*5c3f5608SAlexisMarb   if (!sfMigration) {
119*5c3f5608SAlexisMarb     /* If sfMigration is missing,
120*5c3f5608SAlexisMarb     sfNatural cannot be computed and is set to NULL */
121*5c3f5608SAlexisMarb     *sfNatural = NULL;
122*5c3f5608SAlexisMarb     PetscFunctionReturn(0);
123*5c3f5608SAlexisMarb   } else if (!section) {
124*5c3f5608SAlexisMarb     /* If the sequential section is not provided (NULL),
125*5c3f5608SAlexisMarb     it is reconstructed from the parallel section */
126*5c3f5608SAlexisMarb     PetscSF sfMigrationInv;
127*5c3f5608SAlexisMarb     PetscSection localSection;
128*5c3f5608SAlexisMarb 
129*5c3f5608SAlexisMarb     ierr = DMGetLocalSection(dm, &localSection);CHKERRQ(ierr);
130*5c3f5608SAlexisMarb     ierr = PetscSFCreateInverseSF(sfMigration, &sfMigrationInv);CHKERRQ(ierr);
131*5c3f5608SAlexisMarb     ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &section);CHKERRQ(ierr);
132*5c3f5608SAlexisMarb     ierr = PetscSFDistributeSection(sfMigrationInv, localSection, NULL, section);CHKERRQ(ierr);
133*5c3f5608SAlexisMarb     ierr = PetscSFDestroy(&sfMigrationInv);CHKERRQ(ierr);
134*5c3f5608SAlexisMarb     destroyFlag = PETSC_TRUE;
135*5c3f5608SAlexisMarb   }
136e5b44f4fSMatthew G. Knepley   /* ierr = PetscPrintf(comm, "Point migration SF\n");CHKERRQ(ierr);
137e5b44f4fSMatthew G. Knepley    ierr = PetscSFView(sfMigration, 0);CHKERRQ(ierr); */
138e5b44f4fSMatthew G. Knepley   /* Create a new section from distributing the original section */
139e5b44f4fSMatthew G. Knepley   ierr = PetscSectionCreate(comm, &sectionDist);CHKERRQ(ierr);
140e5b44f4fSMatthew G. Knepley   ierr = PetscSFDistributeSection(sfMigration, section, &remoteOffsets, sectionDist);CHKERRQ(ierr);
141e5b44f4fSMatthew G. Knepley   /* ierr = PetscPrintf(comm, "Distributed Section\n");CHKERRQ(ierr);
142e5b44f4fSMatthew G. Knepley    ierr = PetscSectionView(sectionDist, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */
14392fd8e1eSJed Brown   ierr = DMSetLocalSection(dm, sectionDist);CHKERRQ(ierr);
144*5c3f5608SAlexisMarb   /* If a sequential section is provided but no dof is affected,
145*5c3f5608SAlexisMarb   sfNatural cannot be computed and is set to NULL */
146*5c3f5608SAlexisMarb   ierr = DMCreateGlobalVector(dm, &tmpVec);CHKERRQ(ierr);
147*5c3f5608SAlexisMarb   ierr = VecGetSize(tmpVec, &globalSize);CHKERRQ(ierr);
148*5c3f5608SAlexisMarb   ierr = DMRestoreGlobalVector(dm, &tmpVec);CHKERRQ(ierr);
149*5c3f5608SAlexisMarb   if (!globalSize) {
150*5c3f5608SAlexisMarb     *sfNatural = NULL;
151*5c3f5608SAlexisMarb     if (destroyFlag) {ierr = PetscSectionDestroy(&section);CHKERRQ(ierr);}
152*5c3f5608SAlexisMarb     PetscFunctionReturn(0);
153*5c3f5608SAlexisMarb   }
154fa534816SMatthew G. Knepley   /* Get a pruned version of migration SF */
155e87a4003SBarry Smith   ierr = DMGetGlobalSection(dm, &gSection);CHKERRQ(ierr);
156fa534816SMatthew G. Knepley   ierr = PetscSectionGetChart(gSection, &pStart, &pEnd);CHKERRQ(ierr);
157fa534816SMatthew G. Knepley   for (p = pStart, ssize = 0; p < pEnd; ++p) {
158fa534816SMatthew G. Knepley     PetscInt dof, off;
159fa534816SMatthew G. Knepley 
160fa534816SMatthew G. Knepley     ierr = PetscSectionGetDof(gSection, p, &dof);CHKERRQ(ierr);
161fa534816SMatthew G. Knepley     ierr = PetscSectionGetOffset(gSection, p, &off);CHKERRQ(ierr);
162fa534816SMatthew G. Knepley     if ((dof > 0) && (off >= 0)) ++ssize;
163fa534816SMatthew G. Knepley   }
164fa534816SMatthew G. Knepley   ierr = PetscMalloc1(ssize, &spoints);CHKERRQ(ierr);
165fa534816SMatthew G. Knepley   for (p = pStart, ssize = 0; p < pEnd; ++p) {
166fa534816SMatthew G. Knepley     PetscInt dof, off;
167fa534816SMatthew G. Knepley 
168fa534816SMatthew G. Knepley     ierr = PetscSectionGetDof(gSection, p, &dof);CHKERRQ(ierr);
169fa534816SMatthew G. Knepley     ierr = PetscSectionGetOffset(gSection, p, &off);CHKERRQ(ierr);
170fa534816SMatthew G. Knepley     if ((dof > 0) && (off >= 0)) spoints[ssize++] = p;
171fa534816SMatthew G. Knepley   }
172fa534816SMatthew G. Knepley   ierr = PetscSFCreateEmbeddedLeafSF(sfMigration, ssize, spoints, &sfEmbed);CHKERRQ(ierr);
173e5b44f4fSMatthew G. Knepley   ierr = PetscFree(spoints);CHKERRQ(ierr);
174e5b44f4fSMatthew G. Knepley   /* ierr = PetscPrintf(comm, "Embedded SF\n");CHKERRQ(ierr);
175e5b44f4fSMatthew G. Knepley    ierr = PetscSFView(sfEmbed, 0);CHKERRQ(ierr); */
176fa534816SMatthew G. Knepley   /* Create the SF for seq to natural */
177fa534816SMatthew G. Knepley   ierr = DMGetGlobalVector(dm, &gv);CHKERRQ(ierr);
1783a350576SJunchao Zhang   ierr = VecGetLayout(gv,&map);CHKERRQ(ierr);
1793a350576SJunchao Zhang   /* Note that entries of gv are leaves in sfSeqToNatural, entries of the seq vec are roots */
1803a350576SJunchao Zhang   ierr = PetscSFCreate(comm, &sfSeqToNatural);CHKERRQ(ierr);
1813a350576SJunchao Zhang   ierr = PetscSFSetGraphWithPattern(sfSeqToNatural, map, PETSCSF_PATTERN_GATHER);CHKERRQ(ierr);
182fa534816SMatthew G. Knepley   ierr = DMRestoreGlobalVector(dm, &gv);CHKERRQ(ierr);
183e5b44f4fSMatthew G. Knepley   /* ierr = PetscPrintf(comm, "Seq-to-Natural SF\n");CHKERRQ(ierr);
184e5b44f4fSMatthew G. Knepley    ierr = PetscSFView(sfSeqToNatural, 0);CHKERRQ(ierr); */
185fa534816SMatthew G. Knepley   /* Create the SF associated with this section */
186e5b44f4fSMatthew G. Knepley   ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
187e5b44f4fSMatthew G. Knepley   ierr = PetscSectionCreateGlobalSection(sectionDist, sf, PETSC_FALSE, PETSC_TRUE, &gLocSection);CHKERRQ(ierr);
188e5b44f4fSMatthew G. Knepley   ierr = PetscSFCreateSectionSF(sfEmbed, section, remoteOffsets, gLocSection, &sfField);CHKERRQ(ierr);
1890c374c54SMatthew G. Knepley   ierr = PetscFree(remoteOffsets);CHKERRQ(ierr);
190fa534816SMatthew G. Knepley   ierr = PetscSFDestroy(&sfEmbed);CHKERRQ(ierr);
191e5b44f4fSMatthew G. Knepley   ierr = PetscSectionDestroy(&gLocSection);CHKERRQ(ierr);
192e5b44f4fSMatthew G. Knepley   ierr = PetscSectionDestroy(&sectionDist);CHKERRQ(ierr);
193e5b44f4fSMatthew G. Knepley   /* ierr = PetscPrintf(comm, "Field SF\n");CHKERRQ(ierr);
194e5b44f4fSMatthew G. Knepley    ierr = PetscSFView(sfField, 0);CHKERRQ(ierr); */
195fa534816SMatthew G. Knepley   /* Invert the field SF so it's now from distributed to sequential */
196fa534816SMatthew G. Knepley   ierr = PetscSFCreateInverseSF(sfField, &sfFieldInv);CHKERRQ(ierr);
197fa534816SMatthew G. Knepley   ierr = PetscSFDestroy(&sfField);CHKERRQ(ierr);
198e5b44f4fSMatthew G. Knepley   /* ierr = PetscPrintf(comm, "Inverse Field SF\n");CHKERRQ(ierr);
199e5b44f4fSMatthew G. Knepley    ierr = PetscSFView(sfFieldInv, 0);CHKERRQ(ierr); */
200fa534816SMatthew G. Knepley   /* Multiply the sfFieldInv with the */
2013a350576SJunchao Zhang   ierr = PetscSFComposeInverse(sfFieldInv, sfSeqToNatural, sfNatural);CHKERRQ(ierr);
202e5b44f4fSMatthew G. Knepley   ierr = PetscObjectViewFromOptions((PetscObject) *sfNatural, NULL, "-globaltonatural_sf_view");CHKERRQ(ierr);
203fa534816SMatthew G. Knepley   /* Clean up */
204fa534816SMatthew G. Knepley   ierr = PetscSFDestroy(&sfFieldInv);CHKERRQ(ierr);
205fa534816SMatthew G. Knepley   ierr = PetscSFDestroy(&sfSeqToNatural);CHKERRQ(ierr);
206*5c3f5608SAlexisMarb   if (destroyFlag) {ierr = PetscSectionDestroy(&section);CHKERRQ(ierr);}
207fa534816SMatthew G. Knepley   PetscFunctionReturn(0);
208fa534816SMatthew G. Knepley }
209fa534816SMatthew G. Knepley 
210fa534816SMatthew G. Knepley /*@
211fa534816SMatthew G. Knepley   DMPlexGlobalToNaturalBegin - Rearranges a global Vector in the natural order.
212fa534816SMatthew G. Knepley 
213fa534816SMatthew G. Knepley   Collective on dm
214fa534816SMatthew G. Knepley 
215fa534816SMatthew G. Knepley   Input Parameters:
216fa534816SMatthew G. Knepley + dm - The distributed DMPlex
217fa534816SMatthew G. Knepley - gv - The global Vec
218fa534816SMatthew G. Knepley 
219fa534816SMatthew G. Knepley   Output Parameters:
220fa534816SMatthew G. Knepley . nv - Vec in the canonical ordering distributed over all processors associated with gv
221fa534816SMatthew G. Knepley 
22242ea106eSTristan Konolige   Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute().
223fa534816SMatthew G. Knepley 
224fa534816SMatthew G. Knepley   Level: intermediate
225fa534816SMatthew G. Knepley 
226fa534816SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexNaturalToGlobalBegin(), DMPlexGlobalToNaturalEnd()
227fa534816SMatthew G. Knepley @*/
228fa534816SMatthew G. Knepley PetscErrorCode DMPlexGlobalToNaturalBegin(DM dm, Vec gv, Vec nv)
229fa534816SMatthew G. Knepley {
230fa534816SMatthew G. Knepley   const PetscScalar *inarray;
231fa534816SMatthew G. Knepley   PetscScalar       *outarray;
232a6a55facSMatthew G. Knepley   PetscMPIInt        size;
233fa534816SMatthew G. Knepley   PetscErrorCode     ierr;
234fa534816SMatthew G. Knepley 
235fa534816SMatthew G. Knepley   PetscFunctionBegin;
236fa534816SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_GlobalToNaturalBegin,dm,0,0,0);CHKERRQ(ierr);
237ffc4695bSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);CHKERRMPI(ierr);
238fa534816SMatthew G. Knepley   if (dm->sfNatural) {
239fa534816SMatthew G. Knepley     ierr = VecGetArray(nv, &outarray);CHKERRQ(ierr);
240fa534816SMatthew G. Knepley     ierr = VecGetArrayRead(gv, &inarray);CHKERRQ(ierr);
241fa534816SMatthew G. Knepley     ierr = PetscSFBcastBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray);CHKERRQ(ierr);
242fa534816SMatthew G. Knepley     ierr = VecRestoreArrayRead(gv, &inarray);CHKERRQ(ierr);
243fa534816SMatthew G. Knepley     ierr = VecRestoreArray(nv, &outarray);CHKERRQ(ierr);
244a6a55facSMatthew G. Knepley   } else if (size == 1) {
245cafebca8SMatthew G. Knepley     ierr = VecCopy(gv, nv);CHKERRQ(ierr);
246*5c3f5608SAlexisMarb   } else if (dm->useNatural) SETERRQ(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.\n");
247a6a55facSMatthew 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");
248fa534816SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_GlobalToNaturalBegin,dm,0,0,0);CHKERRQ(ierr);
249fa534816SMatthew G. Knepley   PetscFunctionReturn(0);
250fa534816SMatthew G. Knepley }
251fa534816SMatthew G. Knepley 
252fa534816SMatthew G. Knepley /*@
253fa534816SMatthew G. Knepley   DMPlexGlobalToNaturalEnd - Rearranges a global Vector in the natural order.
254fa534816SMatthew G. Knepley 
255fa534816SMatthew G. Knepley   Collective on dm
256fa534816SMatthew G. Knepley 
257fa534816SMatthew G. Knepley   Input Parameters:
258fa534816SMatthew G. Knepley + dm - The distributed DMPlex
259fa534816SMatthew G. Knepley - gv - The global Vec
260fa534816SMatthew G. Knepley 
261fa534816SMatthew G. Knepley   Output Parameters:
262fa534816SMatthew G. Knepley . nv - The natural Vec
263fa534816SMatthew G. Knepley 
26442ea106eSTristan Konolige   Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute().
265fa534816SMatthew G. Knepley 
266fa534816SMatthew G. Knepley   Level: intermediate
267fa534816SMatthew G. Knepley 
268fa534816SMatthew G. Knepley  .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexNaturalToGlobalBegin(), DMPlexGlobalToNaturalBegin()
269fa534816SMatthew G. Knepley  @*/
270fa534816SMatthew G. Knepley PetscErrorCode DMPlexGlobalToNaturalEnd(DM dm, Vec gv, Vec nv)
271fa534816SMatthew G. Knepley {
272fa534816SMatthew G. Knepley   const PetscScalar *inarray;
273fa534816SMatthew G. Knepley   PetscScalar       *outarray;
274a6a55facSMatthew G. Knepley   PetscMPIInt        size;
275fa534816SMatthew G. Knepley   PetscErrorCode     ierr;
276fa534816SMatthew G. Knepley 
277fa534816SMatthew G. Knepley   PetscFunctionBegin;
278fa534816SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_GlobalToNaturalEnd,dm,0,0,0);CHKERRQ(ierr);
279ffc4695bSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);CHKERRMPI(ierr);
280fa534816SMatthew G. Knepley   if (dm->sfNatural) {
281fa534816SMatthew G. Knepley     ierr = VecGetArrayRead(gv, &inarray);CHKERRQ(ierr);
282c3b366b1Sprj-     ierr = VecGetArray(nv, &outarray);CHKERRQ(ierr);
283fa534816SMatthew G. Knepley     ierr = PetscSFBcastEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray);CHKERRQ(ierr);
284fa534816SMatthew G. Knepley     ierr = VecRestoreArrayRead(gv, &inarray);CHKERRQ(ierr);
285fa534816SMatthew G. Knepley     ierr = VecRestoreArray(nv, &outarray);CHKERRQ(ierr);
286a6a55facSMatthew G. Knepley   } else if (size == 1) {
287*5c3f5608SAlexisMarb   } else if (dm->useNatural) SETERRQ(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.\n");
288a6a55facSMatthew 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");
289fa534816SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_GlobalToNaturalEnd,dm,0,0,0);CHKERRQ(ierr);
290fa534816SMatthew G. Knepley   PetscFunctionReturn(0);
291fa534816SMatthew G. Knepley }
292fa534816SMatthew G. Knepley 
293fa534816SMatthew G. Knepley /*@
294fa534816SMatthew G. Knepley   DMPlexNaturalToGlobalBegin - Rearranges a Vector in the natural order to the Global order.
295fa534816SMatthew G. Knepley 
296fa534816SMatthew G. Knepley   Collective on dm
297fa534816SMatthew G. Knepley 
298fa534816SMatthew G. Knepley   Input Parameters:
299fa534816SMatthew G. Knepley + dm - The distributed DMPlex
300fa534816SMatthew G. Knepley - nv - The natural Vec
301fa534816SMatthew G. Knepley 
302fa534816SMatthew G. Knepley   Output Parameters:
303fa534816SMatthew G. Knepley . gv - The global Vec
304fa534816SMatthew G. Knepley 
30542ea106eSTristan Konolige   Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute().
306fa534816SMatthew G. Knepley 
307fa534816SMatthew G. Knepley   Level: intermediate
308fa534816SMatthew G. Knepley 
309fa534816SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexNaturalToGlobalBegin(),DMPlexGlobalToNaturalEnd()
310fa534816SMatthew G. Knepley @*/
311fa534816SMatthew G. Knepley PetscErrorCode DMPlexNaturalToGlobalBegin(DM dm, Vec nv, Vec gv)
312fa534816SMatthew G. Knepley {
313fa534816SMatthew G. Knepley   const PetscScalar *inarray;
314fa534816SMatthew G. Knepley   PetscScalar       *outarray;
315a6a55facSMatthew G. Knepley   PetscMPIInt        size;
316fa534816SMatthew G. Knepley   PetscErrorCode     ierr;
317fa534816SMatthew G. Knepley 
318fa534816SMatthew G. Knepley   PetscFunctionBegin;
319fa534816SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_NaturalToGlobalBegin,dm,0,0,0);CHKERRQ(ierr);
320ffc4695bSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);CHKERRMPI(ierr);
321fa534816SMatthew G. Knepley   if (dm->sfNatural) {
322b64f75a9SBlaise Bourdin     /* We only have acces to the SF that goes from Global to Natural.
323b64f75a9SBlaise Bourdin        Instead of inverting dm->sfNatural, we can call PetscSFReduceBegin/End with MPI_Op MPI_SUM.
324b64f75a9SBlaise Bourdin        Here the SUM really does nothing since sfNatural is one to one, as long as gV is set to zero first. */
325b64f75a9SBlaise Bourdin     ierr = VecZeroEntries(gv);CHKERRQ(ierr);
326fa534816SMatthew G. Knepley     ierr = VecGetArray(gv, &outarray);CHKERRQ(ierr);
327fa534816SMatthew G. Knepley     ierr = VecGetArrayRead(nv, &inarray);CHKERRQ(ierr);
328fa534816SMatthew G. Knepley     ierr = PetscSFReduceBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray, MPI_SUM);CHKERRQ(ierr);
329fa534816SMatthew G. Knepley     ierr = VecRestoreArrayRead(nv, &inarray);CHKERRQ(ierr);
330fa534816SMatthew G. Knepley     ierr = VecRestoreArray(gv, &outarray);CHKERRQ(ierr);
331a6a55facSMatthew G. Knepley   } else if (size == 1) {
332a6a55facSMatthew G. Knepley     ierr = VecCopy(nv, gv);CHKERRQ(ierr);
333*5c3f5608SAlexisMarb   } else if (dm->useNatural) SETERRQ(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.\n");
334a6a55facSMatthew 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");
335fa534816SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_NaturalToGlobalBegin,dm,0,0,0);CHKERRQ(ierr);
336fa534816SMatthew G. Knepley   PetscFunctionReturn(0);
337fa534816SMatthew G. Knepley }
338fa534816SMatthew G. Knepley 
339fa534816SMatthew G. Knepley /*@
340fa534816SMatthew G. Knepley   DMPlexNaturalToGlobalEnd - Rearranges a Vector in the natural order to the Global order.
341fa534816SMatthew G. Knepley 
342fa534816SMatthew G. Knepley   Collective on dm
343fa534816SMatthew G. Knepley 
344fa534816SMatthew G. Knepley   Input Parameters:
345fa534816SMatthew G. Knepley + dm - The distributed DMPlex
346fa534816SMatthew G. Knepley - nv - The natural Vec
347fa534816SMatthew G. Knepley 
348fa534816SMatthew G. Knepley   Output Parameters:
349fa534816SMatthew G. Knepley . gv - The global Vec
350fa534816SMatthew G. Knepley 
35142ea106eSTristan Konolige   Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute().
352fa534816SMatthew G. Knepley 
353fa534816SMatthew G. Knepley   Level: intermediate
354fa534816SMatthew G. Knepley 
355fa534816SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexNaturalToGlobalBegin(), DMPlexGlobalToNaturalBegin()
356fa534816SMatthew G. Knepley  @*/
357fa534816SMatthew G. Knepley PetscErrorCode DMPlexNaturalToGlobalEnd(DM dm, Vec nv, Vec gv)
358fa534816SMatthew G. Knepley {
359fa534816SMatthew G. Knepley   const PetscScalar *inarray;
360fa534816SMatthew G. Knepley   PetscScalar       *outarray;
361a6a55facSMatthew G. Knepley   PetscMPIInt        size;
362fa534816SMatthew G. Knepley   PetscErrorCode     ierr;
363fa534816SMatthew G. Knepley 
364fa534816SMatthew G. Knepley   PetscFunctionBegin;
365fa534816SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_NaturalToGlobalEnd,dm,0,0,0);CHKERRQ(ierr);
366ffc4695bSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);CHKERRMPI(ierr);
367fa534816SMatthew G. Knepley   if (dm->sfNatural) {
368fa534816SMatthew G. Knepley     ierr = VecGetArrayRead(nv, &inarray);CHKERRQ(ierr);
369c3b366b1Sprj-     ierr = VecGetArray(gv, &outarray);CHKERRQ(ierr);
370fa534816SMatthew G. Knepley     ierr = PetscSFReduceEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray, MPI_SUM);CHKERRQ(ierr);
371fa534816SMatthew G. Knepley     ierr = VecRestoreArrayRead(nv, &inarray);CHKERRQ(ierr);
372fa534816SMatthew G. Knepley     ierr = VecRestoreArray(gv, &outarray);CHKERRQ(ierr);
373a6a55facSMatthew G. Knepley   } else if (size == 1) {
374*5c3f5608SAlexisMarb   } else if (dm->useNatural) SETERRQ(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.\n");
375a6a55facSMatthew 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");
376fa534816SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_NaturalToGlobalEnd,dm,0,0,0);CHKERRQ(ierr);
377fa534816SMatthew G. Knepley   PetscFunctionReturn(0);
378fa534816SMatthew G. Knepley }
379