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