xref: /petsc/src/dm/impls/plex/tests/ex47.c (revision 231de6a21270c79e12cf3957f353a9789180e905)
1*231de6a2SMatthew G. Knepley static char help[] = "The main goal of this code is to retrieve the original element numbers as found in the "
2*231de6a2SMatthew G. Knepley                      "initial partitions (sInitialPartition)... but after the call to DMPlexDistribute";
3*231de6a2SMatthew G. Knepley 
4*231de6a2SMatthew G. Knepley #include <petsc.h>
5*231de6a2SMatthew G. Knepley 
6*231de6a2SMatthew G. Knepley PetscReal sCoords2x5Mesh[18][2] = {
7*231de6a2SMatthew G. Knepley  {0.00000000000000000e+00, 0.00000000000000000e+00},
8*231de6a2SMatthew G. Knepley  {2.00000000000000000e+00, 0.00000000000000000e+00},
9*231de6a2SMatthew G. Knepley  {0.00000000000000000e+00, 1.00000000000000000e+00},
10*231de6a2SMatthew G. Knepley  {2.00000000000000000e+00, 1.00000000000000000e+00},
11*231de6a2SMatthew G. Knepley  {9.99999999997387978e-01, 0.00000000000000000e+00},
12*231de6a2SMatthew G. Knepley  {9.99999999997387978e-01, 1.00000000000000000e+00},
13*231de6a2SMatthew G. Knepley  {0.00000000000000000e+00, 2.00000000000000011e-01},
14*231de6a2SMatthew G. Knepley  {0.00000000000000000e+00, 4.00000000000000022e-01},
15*231de6a2SMatthew G. Knepley  {0.00000000000000000e+00, 5.99999999999999978e-01},
16*231de6a2SMatthew G. Knepley  {0.00000000000000000e+00, 8.00000000000000044e-01},
17*231de6a2SMatthew G. Knepley  {2.00000000000000000e+00, 2.00000000000000011e-01},
18*231de6a2SMatthew G. Knepley  {2.00000000000000000e+00, 4.00000000000000022e-01},
19*231de6a2SMatthew G. Knepley  {2.00000000000000000e+00, 5.99999999999999978e-01},
20*231de6a2SMatthew G. Knepley  {2.00000000000000000e+00, 8.00000000000000044e-01},
21*231de6a2SMatthew G. Knepley  {9.99999999997387756e-01, 2.00000000000000011e-01},
22*231de6a2SMatthew G. Knepley  {9.99999999997387978e-01, 4.00000000000000022e-01},
23*231de6a2SMatthew G. Knepley  {9.99999999997387978e-01, 6.00000000000000089e-01},
24*231de6a2SMatthew G. Knepley  {9.99999999997388089e-01, 8.00000000000000044e-01}};
25*231de6a2SMatthew G. Knepley 
26*231de6a2SMatthew G. Knepley //Connectivity of a 2x5 rectangular mesh of quads :
27*231de6a2SMatthew G. Knepley const PetscInt sConnectivity2x5Mesh[10][4] = {
28*231de6a2SMatthew G. Knepley   {0,4,14,6},
29*231de6a2SMatthew G. Knepley   {6,14,15,7},
30*231de6a2SMatthew G. Knepley   {7,15,16,8},
31*231de6a2SMatthew G. Knepley   {8,16,17,9},
32*231de6a2SMatthew G. Knepley   {9,17,5,2},
33*231de6a2SMatthew G. Knepley   {4,1,10,14},
34*231de6a2SMatthew G. Knepley   {14,10,11,15},
35*231de6a2SMatthew G. Knepley   {15,11,12,16},
36*231de6a2SMatthew G. Knepley   {16,12,13,17},
37*231de6a2SMatthew G. Knepley   {17,13,3,5}};
38*231de6a2SMatthew G. Knepley 
39*231de6a2SMatthew G. Knepley const PetscInt sInitialPartition2x5Mesh[2][5] = {
40*231de6a2SMatthew G. Knepley   {0,2,4,6,8},
41*231de6a2SMatthew G. Knepley   {1,3,5,7,9}
42*231de6a2SMatthew G. Knepley };
43*231de6a2SMatthew G. Knepley 
44*231de6a2SMatthew G. Knepley const PetscInt sNLoclCells2x5Mesh = 5;
45*231de6a2SMatthew G. Knepley const PetscInt sNGlobVerts2x5Mesh = 18;
46*231de6a2SMatthew G. Knepley 
47*231de6a2SMatthew G. Knepley int main(int argc, char **argv)
48*231de6a2SMatthew G. Knepley {
49*231de6a2SMatthew G. Knepley   const PetscInt   Nc                 = sNLoclCells2x5Mesh; //Same on each rank for this example...
50*231de6a2SMatthew G. Knepley   const PetscInt   Nv                 = sNGlobVerts2x5Mesh;
51*231de6a2SMatthew G. Knepley   const PetscInt*  InitPartForRank[2] = {&sInitialPartition2x5Mesh[0][0],
52*231de6a2SMatthew G. Knepley                                          &sInitialPartition2x5Mesh[1][0]};
53*231de6a2SMatthew G. Knepley   const PetscInt (*Conn)[4]           = sConnectivity2x5Mesh;
54*231de6a2SMatthew G. Knepley 
55*231de6a2SMatthew G. Knepley   const PetscInt   Ncor = 4;
56*231de6a2SMatthew G. Knepley   const PetscInt   dim  = 2;
57*231de6a2SMatthew G. Knepley   DM               dm, idm, ddm;
58*231de6a2SMatthew G. Knepley   PetscSF          sfVert, sfMig, sfPart;
59*231de6a2SMatthew G. Knepley   PetscPartitioner part;
60*231de6a2SMatthew G. Knepley   PetscSection     s;
61*231de6a2SMatthew G. Knepley   PetscInt        *cells, c;
62*231de6a2SMatthew G. Knepley   PetscMPIInt      size, rank;
63*231de6a2SMatthew G. Knepley   PetscBool        box = PETSC_FALSE, field = PETSC_FALSE;
64*231de6a2SMatthew G. Knepley 
65*231de6a2SMatthew G. Knepley   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
66*231de6a2SMatthew G. Knepley   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
67*231de6a2SMatthew G. Knepley   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
68*231de6a2SMatthew G. Knepley   PetscCheck(size == 2, PETSC_COMM_WORLD, PETSC_ERR_SUP, "This is a 2 processors example only");
69*231de6a2SMatthew G. Knepley   PetscCall(PetscOptionsGetBool(NULL, NULL, "-box", &box, NULL));
70*231de6a2SMatthew G. Knepley   PetscCall(PetscOptionsGetBool(NULL, NULL, "-field", &field, NULL));
71*231de6a2SMatthew G. Knepley 
72*231de6a2SMatthew G. Knepley   PetscCall(DMPlexCreate(PETSC_COMM_WORLD, &dm));
73*231de6a2SMatthew G. Knepley   if (box) {
74*231de6a2SMatthew G. Knepley     PetscCall(DMSetType(dm, DMPLEX));
75*231de6a2SMatthew G. Knepley     PetscCall(DMSetFromOptions(dm));
76*231de6a2SMatthew G. Knepley   } else {
77*231de6a2SMatthew G. Knepley     PetscCall(PetscMalloc1(Nc * Ncor, &cells));
78*231de6a2SMatthew G. Knepley     for (c = 0; c < Nc; ++c) {
79*231de6a2SMatthew G. Knepley       PetscInt cell = (InitPartForRank[rank])[c], cor;
80*231de6a2SMatthew G. Knepley 
81*231de6a2SMatthew G. Knepley       for (cor = 0; cor < Ncor; ++cor) {
82*231de6a2SMatthew G. Knepley         cells[c*Ncor + cor] = Conn[cell][cor];
83*231de6a2SMatthew G. Knepley       }
84*231de6a2SMatthew G. Knepley     }
85*231de6a2SMatthew G. Knepley     PetscCall(DMSetDimension(dm, dim));
86*231de6a2SMatthew G. Knepley     PetscCall(DMPlexBuildFromCellListParallel(dm, Nc, PETSC_DECIDE, Nv, Ncor, cells, &sfVert, NULL));
87*231de6a2SMatthew G. Knepley     PetscCall(PetscSFDestroy(&sfVert));
88*231de6a2SMatthew G. Knepley     PetscCall(PetscFree(cells));
89*231de6a2SMatthew G. Knepley     PetscCall(DMPlexInterpolate(dm, &idm));
90*231de6a2SMatthew G. Knepley     PetscCall(DMDestroy(&dm));
91*231de6a2SMatthew G. Knepley     dm = idm;
92*231de6a2SMatthew G. Knepley   }
93*231de6a2SMatthew G. Knepley   PetscCall(DMSetUseNatural(dm, PETSC_TRUE));
94*231de6a2SMatthew G. Knepley   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
95*231de6a2SMatthew G. Knepley 
96*231de6a2SMatthew G. Knepley   if (field) {
97*231de6a2SMatthew G. Knepley    const PetscInt Nf         = 1;
98*231de6a2SMatthew G. Knepley    const PetscInt numComp[1] = {1};
99*231de6a2SMatthew G. Knepley    const PetscInt numDof[3]  = {0, 0, 1};
100*231de6a2SMatthew G. Knepley    const PetscInt numBC      = 0;
101*231de6a2SMatthew G. Knepley 
102*231de6a2SMatthew G. Knepley    PetscCall(DMSetNumFields(dm, Nf));
103*231de6a2SMatthew G. Knepley    PetscCall(DMPlexCreateSection(dm, NULL, numComp, numDof, numBC, NULL, NULL, NULL, NULL, &s));
104*231de6a2SMatthew G. Knepley    PetscCall(DMSetLocalSection(dm, s));
105*231de6a2SMatthew G. Knepley    PetscCall(PetscSectionView(s, PETSC_VIEWER_STDOUT_WORLD));
106*231de6a2SMatthew G. Knepley    PetscCall(PetscSectionDestroy(&s));
107*231de6a2SMatthew G. Knepley   }
108*231de6a2SMatthew G. Knepley 
109*231de6a2SMatthew G. Knepley   PetscCall(DMPlexGetPartitioner(dm, &part));
110*231de6a2SMatthew G. Knepley   PetscCall(PetscPartitionerSetFromOptions(part));
111*231de6a2SMatthew G. Knepley 
112*231de6a2SMatthew G. Knepley   PetscCall(DMPlexDistribute(dm, 0, &sfMig, &ddm));
113*231de6a2SMatthew G. Knepley   PetscCall(PetscSFView(sfMig, PETSC_VIEWER_STDOUT_WORLD));
114*231de6a2SMatthew G. Knepley   PetscCall(PetscSFCreateInverseSF(sfMig, &sfPart));
115*231de6a2SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject) sfPart, "Inverse Migration SF"));
116*231de6a2SMatthew G. Knepley   PetscCall(PetscSFView(sfPart, PETSC_VIEWER_STDOUT_WORLD));
117*231de6a2SMatthew G. Knepley 
118*231de6a2SMatthew G. Knepley   Vec          lGlobalVec, lNatVec;
119*231de6a2SMatthew G. Knepley   PetscScalar *lNatVecArray;
120*231de6a2SMatthew G. Knepley 
121*231de6a2SMatthew G. Knepley   {
122*231de6a2SMatthew G. Knepley     PetscSection s;
123*231de6a2SMatthew G. Knepley 
124*231de6a2SMatthew G. Knepley     PetscCall(DMGetGlobalSection(dm, &s));
125*231de6a2SMatthew G. Knepley     PetscCall(PetscSectionView(s, PETSC_VIEWER_STDOUT_WORLD));
126*231de6a2SMatthew G. Knepley   }
127*231de6a2SMatthew G. Knepley   PetscCall(DMGetGlobalVector(dm, &lNatVec));
128*231de6a2SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject) lNatVec, "Natural Vector (initial partition)"));
129*231de6a2SMatthew G. Knepley 
130*231de6a2SMatthew G. Knepley   //Copying the initial partition into the "natural" vector:
131*231de6a2SMatthew G. Knepley   PetscCall(VecGetArray(lNatVec, &lNatVecArray));
132*231de6a2SMatthew G. Knepley   for (c = 0; c < Nc; ++c) lNatVecArray[c] = (InitPartForRank[rank])[c];
133*231de6a2SMatthew G. Knepley   PetscCall(VecRestoreArray(lNatVec, &lNatVecArray));
134*231de6a2SMatthew G. Knepley 
135*231de6a2SMatthew G. Knepley   PetscCall(DMGetGlobalVector(ddm,&lGlobalVec));
136*231de6a2SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject) lGlobalVec, "Global Vector (reordered element numbers in the petsc distributed order)"));
137*231de6a2SMatthew G. Knepley   PetscCall(VecZeroEntries(lGlobalVec));
138*231de6a2SMatthew G. Knepley 
139*231de6a2SMatthew G. Knepley   // The call to DMPlexNaturalToGlobalBegin/End does not produce our expected result...
140*231de6a2SMatthew G. Knepley   // In lGlobalVec, we expect to have:
141*231de6a2SMatthew G. Knepley   /*
142*231de6a2SMatthew G. Knepley    * Process [0]
143*231de6a2SMatthew G. Knepley    * 2.
144*231de6a2SMatthew G. Knepley    * 4.
145*231de6a2SMatthew G. Knepley    * 8.
146*231de6a2SMatthew G. Knepley    * 3.
147*231de6a2SMatthew G. Knepley    * 9.
148*231de6a2SMatthew G. Knepley    * Process [1]
149*231de6a2SMatthew G. Knepley    * 1.
150*231de6a2SMatthew G. Knepley    * 5.
151*231de6a2SMatthew G. Knepley    * 7.
152*231de6a2SMatthew G. Knepley    * 0.
153*231de6a2SMatthew G. Knepley    * 6.
154*231de6a2SMatthew G. Knepley    *
155*231de6a2SMatthew G. Knepley    * but we obtained:
156*231de6a2SMatthew G. Knepley    *
157*231de6a2SMatthew G. Knepley    * Process [0]
158*231de6a2SMatthew G. Knepley    * 2.
159*231de6a2SMatthew G. Knepley    * 4.
160*231de6a2SMatthew G. Knepley    * 8.
161*231de6a2SMatthew G. Knepley    * 0.
162*231de6a2SMatthew G. Knepley    * 0.
163*231de6a2SMatthew G. Knepley    * Process [1]
164*231de6a2SMatthew G. Knepley    * 0.
165*231de6a2SMatthew G. Knepley    * 0.
166*231de6a2SMatthew G. Knepley    * 0.
167*231de6a2SMatthew G. Knepley    * 0.
168*231de6a2SMatthew G. Knepley    * 0.
169*231de6a2SMatthew G. Knepley    */
170*231de6a2SMatthew G. Knepley 
171*231de6a2SMatthew G. Knepley    {
172*231de6a2SMatthew G. Knepley      PetscSF nsf;
173*231de6a2SMatthew G. Knepley 
174*231de6a2SMatthew G. Knepley      PetscCall(DMPlexGetGlobalToNaturalSF(ddm, &nsf));
175*231de6a2SMatthew G. Knepley      PetscCall(PetscSFView(nsf, NULL));
176*231de6a2SMatthew G. Knepley    }
177*231de6a2SMatthew G. Knepley   PetscCall(DMPlexNaturalToGlobalBegin(ddm, lNatVec, lGlobalVec));
178*231de6a2SMatthew G. Knepley   PetscCall(DMPlexNaturalToGlobalEnd  (ddm, lNatVec, lGlobalVec));
179*231de6a2SMatthew G. Knepley 
180*231de6a2SMatthew G. Knepley   PetscCall(VecView(lNatVec, PETSC_VIEWER_STDOUT_WORLD));
181*231de6a2SMatthew G. Knepley   PetscCall(VecView(lGlobalVec, PETSC_VIEWER_STDOUT_WORLD));
182*231de6a2SMatthew G. Knepley 
183*231de6a2SMatthew G. Knepley   PetscCall(DMRestoreGlobalVector(dm,&lNatVec));
184*231de6a2SMatthew G. Knepley   PetscCall(DMRestoreGlobalVector(ddm,&lGlobalVec));
185*231de6a2SMatthew G. Knepley 
186*231de6a2SMatthew G. Knepley   PetscCall(PetscSFDestroy(&sfMig));
187*231de6a2SMatthew G. Knepley   PetscCall(PetscSFDestroy(&sfPart));
188*231de6a2SMatthew G. Knepley   PetscCall(DMDestroy(&dm));
189*231de6a2SMatthew G. Knepley   PetscCall(DMDestroy(&ddm));
190*231de6a2SMatthew G. Knepley   PetscCall(PetscFinalize());
191*231de6a2SMatthew G. Knepley   return 0;
192*231de6a2SMatthew G. Knepley }
193*231de6a2SMatthew G. Knepley 
194*231de6a2SMatthew G. Knepley /*TEST
195*231de6a2SMatthew G. Knepley 
196*231de6a2SMatthew G. Knepley   testset:
197*231de6a2SMatthew G. Knepley     args: -field -petscpartitioner_type simple
198*231de6a2SMatthew G. Knepley     nsize: 2
199*231de6a2SMatthew G. Knepley 
200*231de6a2SMatthew G. Knepley     test:
201*231de6a2SMatthew G. Knepley       suffix: 0
202*231de6a2SMatthew G. Knepley       args:
203*231de6a2SMatthew G. Knepley 
204*231de6a2SMatthew G. Knepley     test:
205*231de6a2SMatthew G. Knepley       suffix: 1
206*231de6a2SMatthew G. Knepley       args: -box -dm_plex_simplex 0 -dm_plex_box_faces 2,5 -dm_distribute
207*231de6a2SMatthew G. Knepley 
208*231de6a2SMatthew G. Knepley TEST*/
209