xref: /petsc/src/dm/impls/plex/tests/ex47.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1231de6a2SMatthew G. Knepley static char help[] = "The main goal of this code is to retrieve the original element numbers as found in the "
2231de6a2SMatthew G. Knepley                      "initial partitions (sInitialPartition)... but after the call to DMPlexDistribute";
3231de6a2SMatthew G. Knepley 
4231de6a2SMatthew G. Knepley #include <petsc.h>
5231de6a2SMatthew G. Knepley 
6231de6a2SMatthew G. Knepley PetscReal sCoords2x5Mesh[18][2] = {
7231de6a2SMatthew G. Knepley   {0.00000000000000000e+00, 0.00000000000000000e+00},
8231de6a2SMatthew G. Knepley   {2.00000000000000000e+00, 0.00000000000000000e+00},
9231de6a2SMatthew G. Knepley   {0.00000000000000000e+00, 1.00000000000000000e+00},
10231de6a2SMatthew G. Knepley   {2.00000000000000000e+00, 1.00000000000000000e+00},
11231de6a2SMatthew G. Knepley   {9.99999999997387978e-01, 0.00000000000000000e+00},
12231de6a2SMatthew G. Knepley   {9.99999999997387978e-01, 1.00000000000000000e+00},
13231de6a2SMatthew G. Knepley   {0.00000000000000000e+00, 2.00000000000000011e-01},
14231de6a2SMatthew G. Knepley   {0.00000000000000000e+00, 4.00000000000000022e-01},
15231de6a2SMatthew G. Knepley   {0.00000000000000000e+00, 5.99999999999999978e-01},
16231de6a2SMatthew G. Knepley   {0.00000000000000000e+00, 8.00000000000000044e-01},
17231de6a2SMatthew G. Knepley   {2.00000000000000000e+00, 2.00000000000000011e-01},
18231de6a2SMatthew G. Knepley   {2.00000000000000000e+00, 4.00000000000000022e-01},
19231de6a2SMatthew G. Knepley   {2.00000000000000000e+00, 5.99999999999999978e-01},
20231de6a2SMatthew G. Knepley   {2.00000000000000000e+00, 8.00000000000000044e-01},
21231de6a2SMatthew G. Knepley   {9.99999999997387756e-01, 2.00000000000000011e-01},
22231de6a2SMatthew G. Knepley   {9.99999999997387978e-01, 4.00000000000000022e-01},
23231de6a2SMatthew G. Knepley   {9.99999999997387978e-01, 6.00000000000000089e-01},
249371c9d4SSatish Balay   {9.99999999997388089e-01, 8.00000000000000044e-01}
259371c9d4SSatish Balay };
26231de6a2SMatthew G. Knepley 
27231de6a2SMatthew G. Knepley //Connectivity of a 2x5 rectangular mesh of quads :
28231de6a2SMatthew G. Knepley const PetscInt sConnectivity2x5Mesh[10][4] = {
29231de6a2SMatthew G. Knepley   {0,  4,  14, 6 },
30231de6a2SMatthew G. Knepley   {6,  14, 15, 7 },
31231de6a2SMatthew G. Knepley   {7,  15, 16, 8 },
32231de6a2SMatthew G. Knepley   {8,  16, 17, 9 },
33231de6a2SMatthew G. Knepley   {9,  17, 5,  2 },
34231de6a2SMatthew G. Knepley   {4,  1,  10, 14},
35231de6a2SMatthew G. Knepley   {14, 10, 11, 15},
36231de6a2SMatthew G. Knepley   {15, 11, 12, 16},
37231de6a2SMatthew G. Knepley   {16, 12, 13, 17},
389371c9d4SSatish Balay   {17, 13, 3,  5 }
399371c9d4SSatish Balay };
40231de6a2SMatthew G. Knepley 
41231de6a2SMatthew G. Knepley const PetscInt sInitialPartition2x5Mesh[2][5] = {
42231de6a2SMatthew G. Knepley   {0, 2, 4, 6, 8},
43231de6a2SMatthew G. Knepley   {1, 3, 5, 7, 9}
44231de6a2SMatthew G. Knepley };
45231de6a2SMatthew G. Knepley 
46231de6a2SMatthew G. Knepley const PetscInt sNLoclCells2x5Mesh = 5;
47231de6a2SMatthew G. Knepley const PetscInt sNGlobVerts2x5Mesh = 18;
48231de6a2SMatthew G. Knepley 
49*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
50*d71ae5a4SJacob Faibussowitsch {
51231de6a2SMatthew G. Knepley   const PetscInt  Nc                 = sNLoclCells2x5Mesh; //Same on each rank for this example...
52231de6a2SMatthew G. Knepley   const PetscInt  Nv                 = sNGlobVerts2x5Mesh;
539371c9d4SSatish Balay   const PetscInt *InitPartForRank[2] = {&sInitialPartition2x5Mesh[0][0], &sInitialPartition2x5Mesh[1][0]};
54231de6a2SMatthew G. Knepley   const PetscInt(*Conn)[4]           = sConnectivity2x5Mesh;
55231de6a2SMatthew G. Knepley 
56231de6a2SMatthew G. Knepley   const PetscInt   Ncor = 4;
57231de6a2SMatthew G. Knepley   const PetscInt   dim  = 2;
58231de6a2SMatthew G. Knepley   DM               dm, idm, ddm;
59231de6a2SMatthew G. Knepley   PetscSF          sfVert, sfMig, sfPart;
60231de6a2SMatthew G. Knepley   PetscPartitioner part;
61231de6a2SMatthew G. Knepley   PetscSection     s;
62231de6a2SMatthew G. Knepley   PetscInt        *cells, c;
63231de6a2SMatthew G. Knepley   PetscMPIInt      size, rank;
64231de6a2SMatthew G. Knepley   PetscBool        box = PETSC_FALSE, field = PETSC_FALSE;
65231de6a2SMatthew G. Knepley 
66327415f7SBarry Smith   PetscFunctionBeginUser;
67231de6a2SMatthew G. Knepley   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
68231de6a2SMatthew G. Knepley   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
69231de6a2SMatthew G. Knepley   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
70231de6a2SMatthew G. Knepley   PetscCheck(size == 2, PETSC_COMM_WORLD, PETSC_ERR_SUP, "This is a 2 processors example only");
71231de6a2SMatthew G. Knepley   PetscCall(PetscOptionsGetBool(NULL, NULL, "-box", &box, NULL));
72231de6a2SMatthew G. Knepley   PetscCall(PetscOptionsGetBool(NULL, NULL, "-field", &field, NULL));
73231de6a2SMatthew G. Knepley 
74231de6a2SMatthew G. Knepley   PetscCall(DMPlexCreate(PETSC_COMM_WORLD, &dm));
75231de6a2SMatthew G. Knepley   if (box) {
76231de6a2SMatthew G. Knepley     PetscCall(DMSetType(dm, DMPLEX));
77231de6a2SMatthew G. Knepley     PetscCall(DMSetFromOptions(dm));
78231de6a2SMatthew G. Knepley   } else {
79231de6a2SMatthew G. Knepley     PetscCall(PetscMalloc1(Nc * Ncor, &cells));
80231de6a2SMatthew G. Knepley     for (c = 0; c < Nc; ++c) {
81231de6a2SMatthew G. Knepley       PetscInt cell = (InitPartForRank[rank])[c], cor;
82231de6a2SMatthew G. Knepley 
83ad540459SPierre Jolivet       for (cor = 0; cor < Ncor; ++cor) cells[c * Ncor + cor] = Conn[cell][cor];
84231de6a2SMatthew G. Knepley     }
85231de6a2SMatthew G. Knepley     PetscCall(DMSetDimension(dm, dim));
86231de6a2SMatthew G. Knepley     PetscCall(DMPlexBuildFromCellListParallel(dm, Nc, PETSC_DECIDE, Nv, Ncor, cells, &sfVert, NULL));
87231de6a2SMatthew G. Knepley     PetscCall(PetscSFDestroy(&sfVert));
88231de6a2SMatthew G. Knepley     PetscCall(PetscFree(cells));
89231de6a2SMatthew G. Knepley     PetscCall(DMPlexInterpolate(dm, &idm));
90231de6a2SMatthew G. Knepley     PetscCall(DMDestroy(&dm));
91231de6a2SMatthew G. Knepley     dm = idm;
92231de6a2SMatthew G. Knepley   }
93231de6a2SMatthew G. Knepley   PetscCall(DMSetUseNatural(dm, PETSC_TRUE));
94231de6a2SMatthew G. Knepley   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
95231de6a2SMatthew G. Knepley 
96231de6a2SMatthew G. Knepley   if (field) {
97231de6a2SMatthew G. Knepley     const PetscInt Nf         = 1;
98231de6a2SMatthew G. Knepley     const PetscInt numComp[1] = {1};
99231de6a2SMatthew G. Knepley     const PetscInt numDof[3]  = {0, 0, 1};
100231de6a2SMatthew G. Knepley     const PetscInt numBC      = 0;
101231de6a2SMatthew G. Knepley 
102231de6a2SMatthew G. Knepley     PetscCall(DMSetNumFields(dm, Nf));
103231de6a2SMatthew G. Knepley     PetscCall(DMPlexCreateSection(dm, NULL, numComp, numDof, numBC, NULL, NULL, NULL, NULL, &s));
104231de6a2SMatthew G. Knepley     PetscCall(DMSetLocalSection(dm, s));
105231de6a2SMatthew G. Knepley     PetscCall(PetscSectionView(s, PETSC_VIEWER_STDOUT_WORLD));
106231de6a2SMatthew G. Knepley     PetscCall(PetscSectionDestroy(&s));
107231de6a2SMatthew G. Knepley   }
108231de6a2SMatthew G. Knepley 
109231de6a2SMatthew G. Knepley   PetscCall(DMPlexGetPartitioner(dm, &part));
110231de6a2SMatthew G. Knepley   PetscCall(PetscPartitionerSetFromOptions(part));
111231de6a2SMatthew G. Knepley 
112231de6a2SMatthew G. Knepley   PetscCall(DMPlexDistribute(dm, 0, &sfMig, &ddm));
113231de6a2SMatthew G. Knepley   PetscCall(PetscSFView(sfMig, PETSC_VIEWER_STDOUT_WORLD));
114231de6a2SMatthew G. Knepley   PetscCall(PetscSFCreateInverseSF(sfMig, &sfPart));
115231de6a2SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)sfPart, "Inverse Migration SF"));
116231de6a2SMatthew G. Knepley   PetscCall(PetscSFView(sfPart, PETSC_VIEWER_STDOUT_WORLD));
117231de6a2SMatthew G. Knepley 
118231de6a2SMatthew G. Knepley   Vec          lGlobalVec, lNatVec;
119231de6a2SMatthew G. Knepley   PetscScalar *lNatVecArray;
120231de6a2SMatthew G. Knepley 
121231de6a2SMatthew G. Knepley   {
122231de6a2SMatthew G. Knepley     PetscSection s;
123231de6a2SMatthew G. Knepley 
124231de6a2SMatthew G. Knepley     PetscCall(DMGetGlobalSection(dm, &s));
125231de6a2SMatthew G. Knepley     PetscCall(PetscSectionView(s, PETSC_VIEWER_STDOUT_WORLD));
126231de6a2SMatthew G. Knepley   }
127231de6a2SMatthew G. Knepley   PetscCall(DMGetGlobalVector(dm, &lNatVec));
128231de6a2SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)lNatVec, "Natural Vector (initial partition)"));
129231de6a2SMatthew G. Knepley 
130231de6a2SMatthew G. Knepley   //Copying the initial partition into the "natural" vector:
131231de6a2SMatthew G. Knepley   PetscCall(VecGetArray(lNatVec, &lNatVecArray));
132231de6a2SMatthew G. Knepley   for (c = 0; c < Nc; ++c) lNatVecArray[c] = (InitPartForRank[rank])[c];
133231de6a2SMatthew G. Knepley   PetscCall(VecRestoreArray(lNatVec, &lNatVecArray));
134231de6a2SMatthew G. Knepley 
135231de6a2SMatthew G. Knepley   PetscCall(DMGetGlobalVector(ddm, &lGlobalVec));
136231de6a2SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)lGlobalVec, "Global Vector (reordered element numbers in the petsc distributed order)"));
137231de6a2SMatthew G. Knepley   PetscCall(VecZeroEntries(lGlobalVec));
138231de6a2SMatthew G. Knepley 
139231de6a2SMatthew G. Knepley   // The call to DMPlexNaturalToGlobalBegin/End does not produce our expected result...
140231de6a2SMatthew G. Knepley   // In lGlobalVec, we expect to have:
141231de6a2SMatthew G. Knepley   /*
142231de6a2SMatthew G. Knepley    * Process [0]
143231de6a2SMatthew G. Knepley    * 2.
144231de6a2SMatthew G. Knepley    * 4.
145231de6a2SMatthew G. Knepley    * 8.
146231de6a2SMatthew G. Knepley    * 3.
147231de6a2SMatthew G. Knepley    * 9.
148231de6a2SMatthew G. Knepley    * Process [1]
149231de6a2SMatthew G. Knepley    * 1.
150231de6a2SMatthew G. Knepley    * 5.
151231de6a2SMatthew G. Knepley    * 7.
152231de6a2SMatthew G. Knepley    * 0.
153231de6a2SMatthew G. Knepley    * 6.
154231de6a2SMatthew G. Knepley    *
155231de6a2SMatthew G. Knepley    * but we obtained:
156231de6a2SMatthew G. Knepley    *
157231de6a2SMatthew G. Knepley    * Process [0]
158231de6a2SMatthew G. Knepley    * 2.
159231de6a2SMatthew G. Knepley    * 4.
160231de6a2SMatthew G. Knepley    * 8.
161231de6a2SMatthew G. Knepley    * 0.
162231de6a2SMatthew G. Knepley    * 0.
163231de6a2SMatthew G. Knepley    * Process [1]
164231de6a2SMatthew G. Knepley    * 0.
165231de6a2SMatthew G. Knepley    * 0.
166231de6a2SMatthew G. Knepley    * 0.
167231de6a2SMatthew G. Knepley    * 0.
168231de6a2SMatthew G. Knepley    * 0.
169231de6a2SMatthew G. Knepley    */
170231de6a2SMatthew G. Knepley 
171231de6a2SMatthew G. Knepley   {
172231de6a2SMatthew G. Knepley     PetscSF nsf;
173231de6a2SMatthew G. Knepley 
174231de6a2SMatthew G. Knepley     PetscCall(DMPlexGetGlobalToNaturalSF(ddm, &nsf));
175231de6a2SMatthew G. Knepley     PetscCall(PetscSFView(nsf, NULL));
176231de6a2SMatthew G. Knepley   }
177231de6a2SMatthew G. Knepley   PetscCall(DMPlexNaturalToGlobalBegin(ddm, lNatVec, lGlobalVec));
178231de6a2SMatthew G. Knepley   PetscCall(DMPlexNaturalToGlobalEnd(ddm, lNatVec, lGlobalVec));
179231de6a2SMatthew G. Knepley 
180231de6a2SMatthew G. Knepley   PetscCall(VecView(lNatVec, PETSC_VIEWER_STDOUT_WORLD));
181231de6a2SMatthew G. Knepley   PetscCall(VecView(lGlobalVec, PETSC_VIEWER_STDOUT_WORLD));
182231de6a2SMatthew G. Knepley 
183231de6a2SMatthew G. Knepley   PetscCall(DMRestoreGlobalVector(dm, &lNatVec));
184231de6a2SMatthew G. Knepley   PetscCall(DMRestoreGlobalVector(ddm, &lGlobalVec));
185231de6a2SMatthew G. Knepley 
186231de6a2SMatthew G. Knepley   PetscCall(PetscSFDestroy(&sfMig));
187231de6a2SMatthew G. Knepley   PetscCall(PetscSFDestroy(&sfPart));
188231de6a2SMatthew G. Knepley   PetscCall(DMDestroy(&dm));
189231de6a2SMatthew G. Knepley   PetscCall(DMDestroy(&ddm));
190231de6a2SMatthew G. Knepley   PetscCall(PetscFinalize());
191231de6a2SMatthew G. Knepley   return 0;
192231de6a2SMatthew G. Knepley }
193231de6a2SMatthew G. Knepley 
194231de6a2SMatthew G. Knepley /*TEST
195231de6a2SMatthew G. Knepley 
196231de6a2SMatthew G. Knepley   testset:
197231de6a2SMatthew G. Knepley     args: -field -petscpartitioner_type simple
198231de6a2SMatthew G. Knepley     nsize: 2
199231de6a2SMatthew G. Knepley 
200231de6a2SMatthew G. Knepley     test:
201231de6a2SMatthew G. Knepley       suffix: 0
202231de6a2SMatthew G. Knepley       args:
203231de6a2SMatthew G. Knepley 
204231de6a2SMatthew G. Knepley     test:
205231de6a2SMatthew G. Knepley       suffix: 1
206231de6a2SMatthew G. Knepley       args: -box -dm_plex_simplex 0 -dm_plex_box_faces 2,5 -dm_distribute
207231de6a2SMatthew G. Knepley 
208231de6a2SMatthew G. Knepley TEST*/
209