xref: /petsc/src/dm/impls/plex/plexnatural.c (revision fa5348168dac2721715c89160270e4baf3e66010)
1*fa534816SMatthew G. Knepley #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2*fa534816SMatthew G. Knepley 
3*fa534816SMatthew G. Knepley #undef __FUNCT__
4*fa534816SMatthew G. Knepley #define __FUNCT__ "DMPlexCreateGlobalToNaturalSF"
5*fa534816SMatthew G. Knepley /*@
6*fa534816SMatthew G. Knepley   DMPlexCreateGlobalToNaturalSF - Creates the SF for mapping Global Vec to the Natural Vec
7*fa534816SMatthew G. Knepley 
8*fa534816SMatthew G. Knepley   Input Parameters:
9*fa534816SMatthew G. Knepley + dm          - The DM
10*fa534816SMatthew G. Knepley . section     - The PetscSection before the mesh was distributed
11*fa534816SMatthew G. Knepley - sfMigration - The PetscSF used to distribute the mesh
12*fa534816SMatthew G. Knepley 
13*fa534816SMatthew G. Knepley   Output Parameters:
14*fa534816SMatthew G. Knepley . sfNatural - PetscSF for mapping the Vec in PETSc ordering to the canonical ordering
15*fa534816SMatthew G. Knepley 
16*fa534816SMatthew G. Knepley   Level: intermediate
17*fa534816SMatthew G. Knepley 
18*fa534816SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField()
19*fa534816SMatthew G. Knepley  @*/
20*fa534816SMatthew G. Knepley PetscErrorCode DMPlexCreateGlobalToNaturalSF(DM dm, PetscSection section, PetscSF sfMigration, PetscSF *sfNatural)
21*fa534816SMatthew G. Knepley {
22*fa534816SMatthew G. Knepley   MPI_Comm       comm;
23*fa534816SMatthew G. Knepley   Vec            gv;
24*fa534816SMatthew G. Knepley   PetscSF        sfEmbed, sfSeqToNatural, sfField, sfFieldInv;
25*fa534816SMatthew G. Knepley   PetscSection   gSection, newSection;
26*fa534816SMatthew G. Knepley   PetscInt      *spoints, *remoteOffsets;
27*fa534816SMatthew G. Knepley   PetscInt       ssize, pStart, pEnd, p;
28*fa534816SMatthew G. Knepley   PetscErrorCode ierr;
29*fa534816SMatthew G. Knepley 
30*fa534816SMatthew G. Knepley   PetscFunctionBegin;
31*fa534816SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
32*fa534816SMatthew G. Knepley   /* Get a pruned version of migration SF */
33*fa534816SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dm, &gSection);CHKERRQ(ierr);
34*fa534816SMatthew G. Knepley   ierr = PetscSectionGetChart(gSection, &pStart, &pEnd);CHKERRQ(ierr);
35*fa534816SMatthew G. Knepley   for (p = pStart, ssize = 0; p < pEnd; ++p) {
36*fa534816SMatthew G. Knepley     PetscInt dof, off;
37*fa534816SMatthew G. Knepley 
38*fa534816SMatthew G. Knepley     ierr = PetscSectionGetDof(gSection, p, &dof);CHKERRQ(ierr);
39*fa534816SMatthew G. Knepley     ierr = PetscSectionGetOffset(gSection, p, &off);CHKERRQ(ierr);
40*fa534816SMatthew G. Knepley     if ((dof > 0) && (off >= 0)) ++ssize;
41*fa534816SMatthew G. Knepley   }
42*fa534816SMatthew G. Knepley   ierr = PetscMalloc1(ssize, &spoints);CHKERRQ(ierr);
43*fa534816SMatthew G. Knepley   for (p = pStart, ssize = 0; p < pEnd; ++p) {
44*fa534816SMatthew G. Knepley     PetscInt dof, off;
45*fa534816SMatthew G. Knepley 
46*fa534816SMatthew G. Knepley     ierr = PetscSectionGetDof(gSection, p, &dof);CHKERRQ(ierr);
47*fa534816SMatthew G. Knepley     ierr = PetscSectionGetOffset(gSection, p, &off);CHKERRQ(ierr);
48*fa534816SMatthew G. Knepley     if ((dof > 0) && (off >= 0)) spoints[ssize++] = p;
49*fa534816SMatthew G. Knepley   }
50*fa534816SMatthew G. Knepley   ierr = PetscSFCreateEmbeddedLeafSF(sfMigration, ssize, spoints, &sfEmbed);CHKERRQ(ierr);
51*fa534816SMatthew G. Knepley   /* Create a new section from distributing the original section */
52*fa534816SMatthew G. Knepley   ierr = PetscSectionCreate(comm, &newSection);CHKERRQ(ierr);
53*fa534816SMatthew G. Knepley   ierr = PetscSFDistributeSection(sfEmbed, section, &remoteOffsets, newSection);CHKERRQ(ierr);
54*fa534816SMatthew G. Knepley   ierr = PetscSFDestroy(&sfEmbed);CHKERRQ(ierr);
55*fa534816SMatthew G. Knepley   /* Create the SF for seq to natural */
56*fa534816SMatthew G. Knepley   ierr = DMGetGlobalVector(dm, &gv);CHKERRQ(ierr);
57*fa534816SMatthew G. Knepley   ierr = PetscSFCreateFromZero(comm, gv, &sfSeqToNatural);CHKERRQ(ierr);
58*fa534816SMatthew G. Knepley   ierr = DMRestoreGlobalVector(dm, &gv);CHKERRQ(ierr);
59*fa534816SMatthew G. Knepley   /* Create the SF associated with this section */
60*fa534816SMatthew G. Knepley   ierr = PetscSFCreateSectionSF(sfEmbed, section, remoteOffsets, newSection, &sfField);CHKERRQ(ierr);
61*fa534816SMatthew G. Knepley   ierr = PetscSFDestroy(&sfEmbed);CHKERRQ(ierr);
62*fa534816SMatthew G. Knepley   /* Invert the field SF so it's now from distributed to sequential */
63*fa534816SMatthew G. Knepley   ierr = PetscSFCreateInverseSF(sfField, &sfFieldInv);CHKERRQ(ierr);
64*fa534816SMatthew G. Knepley   ierr = PetscSFDestroy(&sfField);CHKERRQ(ierr);
65*fa534816SMatthew G. Knepley   /* Multiply the sfFieldInv with the */
66*fa534816SMatthew G. Knepley   ierr = PetscSFCompose(sfFieldInv, sfSeqToNatural, sfNatural);CHKERRQ(ierr);
67*fa534816SMatthew G. Knepley   /* Clean up */
68*fa534816SMatthew G. Knepley   ierr = PetscSectionDestroy(&newSection);CHKERRQ(ierr);
69*fa534816SMatthew G. Knepley   ierr = PetscSFDestroy(&sfFieldInv);CHKERRQ(ierr);
70*fa534816SMatthew G. Knepley   ierr = PetscSFDestroy(&sfSeqToNatural);CHKERRQ(ierr);
71*fa534816SMatthew G. Knepley   PetscFunctionReturn(0);
72*fa534816SMatthew G. Knepley }
73*fa534816SMatthew G. Knepley 
74*fa534816SMatthew G. Knepley #undef __FUNCT__
75*fa534816SMatthew G. Knepley #define __FUNCT__ "DMPlexGlobalToNaturalBegin"
76*fa534816SMatthew G. Knepley /*@
77*fa534816SMatthew G. Knepley   DMPlexGlobalToNaturalBegin - Rearranges a global Vector in the natural order.
78*fa534816SMatthew G. Knepley 
79*fa534816SMatthew G. Knepley   Collective on dm
80*fa534816SMatthew G. Knepley 
81*fa534816SMatthew G. Knepley   Input Parameters:
82*fa534816SMatthew G. Knepley + dm - The distributed DMPlex
83*fa534816SMatthew G. Knepley - gv - The global Vec
84*fa534816SMatthew G. Knepley 
85*fa534816SMatthew G. Knepley   Output Parameters:
86*fa534816SMatthew G. Knepley . nv - Vec in the canonical ordering distributed over all processors associated with gv
87*fa534816SMatthew G. Knepley 
88*fa534816SMatthew G. Knepley   Note: The user must call DMPlexSetUseNaturalSF(dm, PETSC_TRUE) before DMPlexDistribute().
89*fa534816SMatthew G. Knepley 
90*fa534816SMatthew G. Knepley   Level: intermediate
91*fa534816SMatthew G. Knepley 
92*fa534816SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexNaturalToGlobalBegin(), DMPlexGlobalToNaturalEnd()
93*fa534816SMatthew G. Knepley @*/
94*fa534816SMatthew G. Knepley PetscErrorCode DMPlexGlobalToNaturalBegin(DM dm, Vec gv, Vec nv)
95*fa534816SMatthew G. Knepley {
96*fa534816SMatthew G. Knepley   const PetscScalar *inarray;
97*fa534816SMatthew G. Knepley   PetscScalar       *outarray;
98*fa534816SMatthew G. Knepley   PetscErrorCode     ierr;
99*fa534816SMatthew G. Knepley 
100*fa534816SMatthew G. Knepley   PetscFunctionBegin;
101*fa534816SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_GlobalToNaturalBegin,dm,0,0,0);CHKERRQ(ierr);
102*fa534816SMatthew G. Knepley   if (dm->sfNatural) {
103*fa534816SMatthew G. Knepley     ierr = VecGetArray(nv, &outarray);CHKERRQ(ierr);
104*fa534816SMatthew G. Knepley     ierr = VecGetArrayRead(gv, &inarray);CHKERRQ(ierr);
105*fa534816SMatthew G. Knepley     ierr = PetscSFBcastBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray);CHKERRQ(ierr);
106*fa534816SMatthew G. Knepley     ierr = VecRestoreArrayRead(gv, &inarray);CHKERRQ(ierr);
107*fa534816SMatthew G. Knepley     ierr = VecRestoreArray(nv, &outarray);CHKERRQ(ierr);
108*fa534816SMatthew 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");
109*fa534816SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_GlobalToNaturalBegin,dm,0,0,0);CHKERRQ(ierr);
110*fa534816SMatthew G. Knepley   PetscFunctionReturn(0);
111*fa534816SMatthew G. Knepley }
112*fa534816SMatthew G. Knepley 
113*fa534816SMatthew G. Knepley #undef __FUNCT__
114*fa534816SMatthew G. Knepley #define __FUNCT__ "DMPlexGlobalToNaturalEnd"
115*fa534816SMatthew G. Knepley /*@
116*fa534816SMatthew G. Knepley   DMPlexGlobalToNaturalEnd - Rearranges a global Vector in the natural order.
117*fa534816SMatthew G. Knepley 
118*fa534816SMatthew G. Knepley   Collective on dm
119*fa534816SMatthew G. Knepley 
120*fa534816SMatthew G. Knepley   Input Parameters:
121*fa534816SMatthew G. Knepley + dm - The distributed DMPlex
122*fa534816SMatthew G. Knepley - gv - The global Vec
123*fa534816SMatthew G. Knepley 
124*fa534816SMatthew G. Knepley   Output Parameters:
125*fa534816SMatthew G. Knepley . nv - The natural Vec
126*fa534816SMatthew G. Knepley 
127*fa534816SMatthew G. Knepley   Note: The user must call DMPlexSetUseNaturalSF(dm, PETSC_TRUE) before DMPlexDistribute().
128*fa534816SMatthew G. Knepley 
129*fa534816SMatthew G. Knepley   Level: intermediate
130*fa534816SMatthew G. Knepley 
131*fa534816SMatthew G. Knepley  .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexNaturalToGlobalBegin(), DMPlexGlobalToNaturalBegin()
132*fa534816SMatthew G. Knepley  @*/
133*fa534816SMatthew G. Knepley PetscErrorCode DMPlexGlobalToNaturalEnd(DM dm, Vec gv, Vec nv)
134*fa534816SMatthew G. Knepley {
135*fa534816SMatthew G. Knepley   const PetscScalar *inarray;
136*fa534816SMatthew G. Knepley   PetscScalar       *outarray;
137*fa534816SMatthew G. Knepley   PetscErrorCode     ierr;
138*fa534816SMatthew G. Knepley 
139*fa534816SMatthew G. Knepley   PetscFunctionBegin;
140*fa534816SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_GlobalToNaturalEnd,dm,0,0,0);CHKERRQ(ierr);
141*fa534816SMatthew G. Knepley   if (dm->sfNatural) {
142*fa534816SMatthew G. Knepley     ierr = VecGetArrayRead(gv, &inarray);CHKERRQ(ierr);
143*fa534816SMatthew G. Knepley     ierr = VecGetArray(nv, &outarray);CHKERRQ(ierr);CHKERRQ(ierr);
144*fa534816SMatthew G. Knepley     ierr = PetscSFBcastEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray);CHKERRQ(ierr);
145*fa534816SMatthew G. Knepley     ierr = VecRestoreArrayRead(gv, &inarray);CHKERRQ(ierr);
146*fa534816SMatthew G. Knepley     ierr = VecRestoreArray(nv, &outarray);CHKERRQ(ierr);
147*fa534816SMatthew G. Knepley   }
148*fa534816SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_GlobalToNaturalEnd,dm,0,0,0);CHKERRQ(ierr);
149*fa534816SMatthew G. Knepley   PetscFunctionReturn(0);
150*fa534816SMatthew G. Knepley }
151*fa534816SMatthew G. Knepley 
152*fa534816SMatthew G. Knepley #undef __FUNCT__
153*fa534816SMatthew G. Knepley #define __FUNCT__ "DMPlexNaturalToGlobalBegin"
154*fa534816SMatthew G. Knepley /*@
155*fa534816SMatthew G. Knepley   DMPlexNaturalToGlobalBegin - Rearranges a Vector in the natural order to the Global order.
156*fa534816SMatthew G. Knepley 
157*fa534816SMatthew G. Knepley   Collective on dm
158*fa534816SMatthew G. Knepley 
159*fa534816SMatthew G. Knepley   Input Parameters:
160*fa534816SMatthew G. Knepley + dm - The distributed DMPlex
161*fa534816SMatthew G. Knepley - nv - The natural Vec
162*fa534816SMatthew G. Knepley 
163*fa534816SMatthew G. Knepley   Output Parameters:
164*fa534816SMatthew G. Knepley . gv - The global Vec
165*fa534816SMatthew G. Knepley 
166*fa534816SMatthew G. Knepley   Note: The user must call DMPlexSetUseNaturalSF(dm, PETSC_TRUE) before DMPlexDistribute().
167*fa534816SMatthew G. Knepley 
168*fa534816SMatthew G. Knepley   Level: intermediate
169*fa534816SMatthew G. Knepley 
170*fa534816SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexNaturalToGlobalBegin(),DMPlexGlobalToNaturalEnd()
171*fa534816SMatthew G. Knepley @*/
172*fa534816SMatthew G. Knepley PetscErrorCode DMPlexNaturalToGlobalBegin(DM dm, Vec nv, Vec gv)
173*fa534816SMatthew G. Knepley {
174*fa534816SMatthew G. Knepley   const PetscScalar *inarray;
175*fa534816SMatthew G. Knepley   PetscScalar       *outarray;
176*fa534816SMatthew G. Knepley   PetscErrorCode     ierr;
177*fa534816SMatthew G. Knepley 
178*fa534816SMatthew G. Knepley   PetscFunctionBegin;
179*fa534816SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_NaturalToGlobalBegin,dm,0,0,0);CHKERRQ(ierr);
180*fa534816SMatthew G. Knepley   if (dm->sfNatural) {
181*fa534816SMatthew G. Knepley     ierr = VecGetArray(gv, &outarray);CHKERRQ(ierr);
182*fa534816SMatthew G. Knepley     ierr = VecGetArrayRead(nv, &inarray);CHKERRQ(ierr);
183*fa534816SMatthew G. Knepley     ierr = PetscSFReduceBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray, MPI_SUM);CHKERRQ(ierr);
184*fa534816SMatthew G. Knepley     ierr = VecRestoreArrayRead(nv, &inarray);CHKERRQ(ierr);
185*fa534816SMatthew G. Knepley     ierr = VecRestoreArray(gv, &outarray);CHKERRQ(ierr);
186*fa534816SMatthew 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");
187*fa534816SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_NaturalToGlobalBegin,dm,0,0,0);CHKERRQ(ierr);
188*fa534816SMatthew G. Knepley   PetscFunctionReturn(0);
189*fa534816SMatthew G. Knepley }
190*fa534816SMatthew G. Knepley 
191*fa534816SMatthew G. Knepley #undef __FUNCT__
192*fa534816SMatthew G. Knepley #define __FUNCT__ "DMPlexNaturalToGlobalEnd"
193*fa534816SMatthew G. Knepley /*@
194*fa534816SMatthew G. Knepley   DMPlexNaturalToGlobalEnd - Rearranges a Vector in the natural order to the Global order.
195*fa534816SMatthew G. Knepley 
196*fa534816SMatthew G. Knepley   Collective on dm
197*fa534816SMatthew G. Knepley 
198*fa534816SMatthew G. Knepley   Input Parameters:
199*fa534816SMatthew G. Knepley + dm - The distributed DMPlex
200*fa534816SMatthew G. Knepley - nv - The natural Vec
201*fa534816SMatthew G. Knepley 
202*fa534816SMatthew G. Knepley   Output Parameters:
203*fa534816SMatthew G. Knepley . gv - The global Vec
204*fa534816SMatthew G. Knepley 
205*fa534816SMatthew G. Knepley   Note: The user must call DMPlexSetUseNaturalSF(dm, PETSC_TRUE) before DMPlexDistribute().
206*fa534816SMatthew G. Knepley 
207*fa534816SMatthew G. Knepley   Level: intermediate
208*fa534816SMatthew G. Knepley 
209*fa534816SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexNaturalToGlobalBegin(), DMPlexGlobalToNaturalBegin()
210*fa534816SMatthew G. Knepley  @*/
211*fa534816SMatthew G. Knepley PetscErrorCode DMPlexNaturalToGlobalEnd(DM dm, Vec nv, Vec gv)
212*fa534816SMatthew G. Knepley {
213*fa534816SMatthew G. Knepley   const PetscScalar *inarray;
214*fa534816SMatthew G. Knepley   PetscScalar       *outarray;
215*fa534816SMatthew G. Knepley   PetscErrorCode     ierr;
216*fa534816SMatthew G. Knepley 
217*fa534816SMatthew G. Knepley   PetscFunctionBegin;
218*fa534816SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_NaturalToGlobalEnd,dm,0,0,0);CHKERRQ(ierr);
219*fa534816SMatthew G. Knepley   if (dm->sfNatural) {
220*fa534816SMatthew G. Knepley     ierr = VecGetArrayRead(nv, &inarray);CHKERRQ(ierr);
221*fa534816SMatthew G. Knepley     ierr = VecGetArray(gv, &outarray);CHKERRQ(ierr);CHKERRQ(ierr);
222*fa534816SMatthew G. Knepley     ierr = PetscSFReduceEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray, MPI_SUM);CHKERRQ(ierr);
223*fa534816SMatthew G. Knepley     ierr = VecRestoreArrayRead(nv, &inarray);CHKERRQ(ierr);
224*fa534816SMatthew G. Knepley     ierr = VecRestoreArray(gv, &outarray);CHKERRQ(ierr);
225*fa534816SMatthew G. Knepley   }
226*fa534816SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_NaturalToGlobalEnd,dm,0,0,0);CHKERRQ(ierr);
227*fa534816SMatthew G. Knepley   PetscFunctionReturn(0);
228*fa534816SMatthew G. Knepley }
229