Tpetra parallel linear algebra Version of the Day
Loading...
Searching...
No Matches
Tpetra_Distribution2D.hpp
1// @HEADER
2// *****************************************************************************
3// Tpetra: Templated Linear Algebra Services Package
4//
5// Copyright 2008 NTESS and the Tpetra contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10// 2D matrix distribution
11// Assumes square matrix
12// Karen Devine, SNL
13//
14
15#ifndef __TPETRA_DISTRIBUTION2D_HPP
16#define __TPETRA_DISTRIBUTION2D_HPP
17
18namespace Tpetra
19{
20
22template <typename gno_t, typename scalar_t>
23class Distribution2D : public Distribution<gno_t,scalar_t> {
24// Processors will be laid out logically first down columns then
25// across rows. For example, assume np = 8, npRows=2, npCols=4.
26// Then the processors will be laid out in 2D as
27// 0 2 4 6
28// 1 3 5 7
29//
30// The matrix will be distributed using np=8 row blocks:
31// 0 2 4 6
32// 1 3 5 7
33// 0 2 4 6
34// 1 3 5 7
35// 0 2 4 6
36// 1 3 5 7
37// 0 2 4 6
38// 1 3 5 7
39//
40// The vector will be distributed linearly or randomly. The row and
41// column maps will be built to allow only row- or column-based
42// communication in the matvec.
43
44public:
45 using Distribution<gno_t,scalar_t>::me;
46 using Distribution<gno_t,scalar_t>::np;
47 using Distribution<gno_t,scalar_t>::comm;
48 using Distribution<gno_t,scalar_t>::nrows;
49 using Distribution<gno_t,scalar_t>::Mine;
50
51 Distribution2D(size_t nrows_,
52 const Teuchos::RCP<const Teuchos::Comm<int> > &comm_,
53 const Teuchos::ParameterList &params) :
54 Distribution<gno_t,scalar_t>(nrows_, comm_, params),
55 npRows(-1), npCols(-1)
56 {
57
58 {
59 const Teuchos::ParameterEntry *pe = params.getEntryPtr("nProcessorRows");
60 if (pe != NULL)
61 npRows = pe->getValue<int>(&npRows);
62 }
63
64 {
65 const Teuchos::ParameterEntry *pe = params.getEntryPtr("nProcessorCols");
66 if (pe != NULL)
67 npCols = pe->getValue<int>(&npCols);
68 }
69
70 // Compute the processor configuration npRows * npCols
71
72 if (npRows == -1 && npCols == -1) { // Compute both npRows and npCols
73 // First guess
74 npRows = (int)(sqrt(np));
75 npCols = np / npRows;
76 // Adjust npRows so that npRows * npCols == np
77 while (npRows * npCols != np) {
78 npRows++;
79 npCols = np / npRows;
80 }
81 }
82 else { // User specified either npRows or npCols
83 if (npRows == -1) // npCols specified; compute npRows
84 npRows = np / npCols;
85 else if (npCols == -1) // npRows specified; compute npCols
86 npCols = np / npRows;
87
88 if (npCols * npRows != np) {
89 TEUCHOS_TEST_FOR_EXCEPTION(npRows * npCols != np, std::logic_error,
90 "nProcessorCols " << npCols <<
91 " * nProcessorRows " << npRows <<
92 " = " << npCols * npRows <<
93 " must equal nProcessors " << np <<
94 " for 2D distribution");
95 }
96 }
97 if (me == 0)
98 std::cout << "\n2D Distribution: "
99 << "\n npRows = " << npRows
100 << "\n npCols = " << npCols
101 << "\n np = " << np << std::endl;
102
103 mypCol = this->TWODPCOL(me);
104 mypRow = this->TWODPROW(me);
105 }
106
107 virtual ~Distribution2D() {};
108
109protected:
110
111 // Return the processor row for rank p
112 inline int TWODPROW(int p) {return (p % npRows);}
113
114 // Return the processor column for rank p
115 inline int TWODPCOL(int p) {return (p / npRows);}
116
117 // Return the rank for processor row i and processor column j
118 inline int TWODPRANK(int i, int j) {return (j * npRows + (j+i) % npRows);}
119
120 int npRows; // Number of processor rows
121 int npCols; // Number of processor columns
122 int mypRow; // This rank's processor row
123 int mypCol; // This rank's processor column
124};
125
127template <typename gno_t, typename scalar_t>
128class Distribution2DLinear : public Distribution2D<gno_t,scalar_t> {
129
130public:
131 using Distribution<gno_t,scalar_t>::me;
132 using Distribution<gno_t,scalar_t>::np;
133 using Distribution<gno_t,scalar_t>::nrows;
134 using Distribution2D<gno_t,scalar_t>::npRows;
135 using Distribution2D<gno_t,scalar_t>::npCols;
136 using Distribution2D<gno_t,scalar_t>::mypRow;
137 using Distribution2D<gno_t,scalar_t>::mypCol;
138
139 Distribution2DLinear(size_t nrows_,
140 const Teuchos::RCP<const Teuchos::Comm<int> > &comm_,
141 const Teuchos::ParameterList &params) :
142 Distribution2D<gno_t,scalar_t>(nrows_, comm_, params)
143 {
144 // Build vector describing distribution of vector entries across ranks
145 entries.assign(np+1, nrows / np);
146 int nExtraEntries = nrows % np;
147
148 // Distribute the extra entries evenly among processors.
149 // To evenly distribute them extra entries among processor rows and
150 // columns, we distribute them along diagonals of the matrix distribution.
151 // For our example, assume we have seven extra values (the max possible
152 // with np=8). Then we give one extra entry each to ranks
153 // [0, 3, 4, 7, 1, 2, 5]. For fewer extra entries, we follow the same
154 // order of assignment, and just stop early.
155
156 for (int cnt = 0, i = 0; (cnt < nExtraEntries) && (i < npRows); i++) {
157 for (int j = 0; (cnt < nExtraEntries) && (j < npCols); cnt++, j++) {
158 int rankForExtra = Distribution2D<gno_t,scalar_t>::TWODPRANK(i, j);
159 entries[rankForExtra+1]++; // Store in rankForExtra+1 to simplify
160 // prefix sum.
161 }
162 }
163
164 // Perform prefix sum of entries.
165 entries[0] = 0;
166 for (int i = 1; i <= np; i++)
167 entries[i] = entries[i-1] + entries[i];
168 // Now entries contains the first vector entry for each rank.
169
170 // Column map: Easy; consecutive entries for all ranks in column.
171 int firstRank = mypCol * npRows; // First rank in my column
172 myFirstCol = entries[firstRank];
173
174 gno_t nMyCols = 0;
175 for (int i = firstRank; i < firstRank + npRows; i++)
176 nMyCols += entries[i+1] - entries[i];
177 myLastCol = myFirstCol + nMyCols - 1;
178 }
179
180 inline enum DistributionType DistType() { return TwoDLinear; }
181
182 bool Mine(gno_t i, gno_t j) {
183 int idx = int(float(i) * float(np) / float(nrows));
184 while (i < entries[idx]) idx--;
185 while (i >= entries[idx+1]) idx++;
186 return ((mypRow == Distribution2D<gno_t,scalar_t>::TWODPROW(idx)) // RowMine
187 && (j >= myFirstCol && j <= myLastCol)); // ColMine
188 }
189 inline bool Mine(gno_t i, gno_t j, int p) {return Mine(i,j);}
190
191 inline bool VecMine(gno_t i) {
192 return(i >= entries[me] && i < entries[me+1]);
193 }
194
195
196private:
197 std::vector<gno_t> entries; // Describes vector entries' distribution to ranks
198 // Organized like vtxdist
199 gno_t myFirstCol; // First column owned by this rank
200 gno_t myLastCol; // Last column owned by this rank
201};
202
204template <typename gno_t, typename scalar_t>
205class Distribution2DRandom : public Distribution2D<gno_t,scalar_t> {
206
207public:
208 using Distribution<gno_t,scalar_t>::me;
209 using Distribution2D<gno_t,scalar_t>::mypRow;
210 using Distribution2D<gno_t,scalar_t>::mypCol;
211
212 Distribution2DRandom(size_t nrows_,
213 const Teuchos::RCP<const Teuchos::Comm<int> > &comm_,
214 const Teuchos::ParameterList &params) :
215 Distribution2D<gno_t,scalar_t>(nrows_, comm_, params)
216 { if (me == 0) std::cout << " randomize = true" << std::endl; }
217
218 inline enum DistributionType DistType() { return TwoDRandom; }
219
220 inline bool Mine(gno_t i, gno_t j) {
221 return ((mypRow == this->TWODPROW(this->HashToProc(i))) && // RowMine
222 (mypCol == this->TWODPCOL(this->HashToProc(j)))); // ColMine
223 }
224 inline bool Mine(gno_t i, gno_t j, int p) {return Mine(i,j);}
225
226 inline bool VecMine(gno_t i) { return (me == this->HashToProc(i)); }
227
228};
229
231
232template <typename gno_t, typename scalar_t>
233class Distribution2DVec : public Distribution2D<gno_t,scalar_t>
234{
235// Distribute non-zeros in a 2D manner based on the vector distribution
236// and the nprows x npcols configuration;
237// rows are assigned to same process owning the vector entry.
238public:
239 using Distribution<gno_t,scalar_t>::me;
240 using Distribution<gno_t,scalar_t>::np;
241 using Distribution<gno_t,scalar_t>::comm;
242 using Distribution<gno_t,scalar_t>::nrows;
243 using Distribution2D<gno_t,scalar_t>::npRows;
244 using Distribution2D<gno_t,scalar_t>::npCols;
245
246 Distribution2DVec(size_t nrows_,
247 const Teuchos::RCP<const Teuchos::Comm<int> > &comm_,
248 const Teuchos::ParameterList &params,
249 std::string &distributionfile) :
250 Distribution2D<gno_t,scalar_t>(nrows_, comm_, params)
251 {
252 if (me == 0) std::cout << "\n 2DVec Distribution: "
253 << "\n np = " << np << std::endl;
254 std::ifstream fpin;
255 if (me == 0) {
256 fpin.open(distributionfile.c_str(), std::ios::in);
257 if (!fpin.is_open()) {
258 std::cout << "Error: distributionfile " << distributionfile
259 << " not found" << std::endl;
260 exit(-1);
261 }
262 }
263
264 // Read the vector part assignment and broadcast it to all processes.
265 // Broadcast in chunks of bcastsize values.
266 // TODO: Make the vector part assignment more scalable instead of
267 // TODO: storing entire vector on every process.
268
269 vecpart = new int[nrows];
270
271 const int bcastsize = 1000000;
272
273 gno_t start = 0;
274 int cnt = 0;
275 for (size_t i = 0; i < nrows; i++) {
276 if (me == 0) fpin >> vecpart[i];
277 cnt++;
278 if (cnt == bcastsize || i == nrows-1) {
279 Teuchos::broadcast(*comm, 0, cnt, &(vecpart[start]));
280 start += cnt;
281 cnt = 0;
282 }
283 }
284
285 if (me == 0) fpin.close();
286 }
287
288 ~Distribution2DVec() {delete [] vecpart;}
289
290 inline enum DistributionType DistType() { return TwoDVec; }
291
292 bool Mine(gno_t i, gno_t j) {
293 return (me == (vecpart[i] % npRows + (vecpart[j] / npRows) * npRows));
294 }
295 inline bool Mine(gno_t i, gno_t j, int p) {return Mine(i,j);}
296
297 inline bool VecMine(gno_t i) { return(vecpart[i] == me); }
298
299private:
300 int *vecpart;
301
302};
303
304}
305#endif
void start()
Start the deep_copy counter.
Namespace Tpetra contains the class and methods constituting the Tpetra library.