PetaVision  Alpha
BorderExchange.cpp
1 /*
2  * BorderExchange.cpp
3  */
4 
5 #include "BorderExchange.hpp"
6 #include "include/pv_common.h"
7 #include "utils/PVAssert.hpp"
8 #include "utils/conversions.h"
9 
10 namespace PV {
11 
12 BorderExchange::BorderExchange(MPIBlock const &mpiBlock, PVLayerLoc const &loc) {
13  mMPIBlock = &mpiBlock; // TODO: copy the block instead of storing a reference.
14  mLayerLoc = loc;
15  newDatatypes();
16  initNeighbors();
17 }
18 
19 BorderExchange::~BorderExchange() { freeDatatypes(); }
20 
21 void BorderExchange::newDatatypes() {
22 #ifdef PV_USE_MPI
23  int count, blocklength, stride;
24  mDatatypes.resize(NUM_NEIGHBORHOOD);
25 
26  int const leftBorder = mLayerLoc.halo.lt;
27  int const rightBorder = mLayerLoc.halo.rt;
28  int const bottomBorder = mLayerLoc.halo.dn;
29  int const topBorder = mLayerLoc.halo.up;
30 
31  int const nf = mLayerLoc.nf;
32 
33  count = mLayerLoc.ny;
34  blocklength = nf * mLayerLoc.nx;
35  stride = nf * (mLayerLoc.nx + leftBorder + rightBorder);
36 
37  /* local interior */
38  MPI_Type_vector(count, blocklength, stride, MPI_FLOAT, &mDatatypes[LOCAL]);
39  MPI_Type_commit(&mDatatypes[LOCAL]);
40 
41  count = topBorder;
42 
43  /* northwest */
44  blocklength = nf * leftBorder;
45  MPI_Type_vector(count, blocklength, stride, MPI_FLOAT, &mDatatypes[NORTHWEST]);
46  MPI_Type_commit(&mDatatypes[NORTHWEST]);
47 
48  /* north */
49  blocklength = nf * mLayerLoc.nx;
50  MPI_Type_vector(count, blocklength, stride, MPI_FLOAT, &mDatatypes[NORTH]);
51  MPI_Type_commit(&mDatatypes[NORTH]);
52 
53  /* northeast */
54  blocklength = nf * rightBorder;
55  MPI_Type_vector(count, blocklength, stride, MPI_FLOAT, &mDatatypes[NORTHEAST]);
56  MPI_Type_commit(&mDatatypes[NORTHEAST]);
57 
58  count = mLayerLoc.ny;
59 
60  /* west */
61  blocklength = nf * leftBorder;
62  MPI_Type_vector(count, blocklength, stride, MPI_FLOAT, &mDatatypes[WEST]);
63  MPI_Type_commit(&mDatatypes[WEST]);
64 
65  /* east */
66  blocklength = nf * rightBorder;
67  MPI_Type_vector(count, blocklength, stride, MPI_FLOAT, &mDatatypes[EAST]);
68  MPI_Type_commit(&mDatatypes[EAST]);
69 
70  count = bottomBorder;
71 
72  /* southwest */
73  blocklength = nf * leftBorder;
74  MPI_Type_vector(count, blocklength, stride, MPI_FLOAT, &mDatatypes[SOUTHWEST]);
75  MPI_Type_commit(&mDatatypes[SOUTHWEST]);
76 
77  /* south */
78  blocklength = nf * mLayerLoc.nx;
79  MPI_Type_vector(count, blocklength, stride, MPI_FLOAT, &mDatatypes[SOUTH]);
80  MPI_Type_commit(&mDatatypes[SOUTH]);
81 
82  /* southeast */
83  blocklength = nf * rightBorder;
84  MPI_Type_vector(count, blocklength, stride, MPI_FLOAT, &mDatatypes[SOUTHEAST]);
85  MPI_Type_commit(&mDatatypes[SOUTHEAST]);
86 #else // PV_USE_MPI
87  mDatatypes.clear();
88 #endif // PV_USE_MPI
89 }
90 
91 void BorderExchange::freeDatatypes() {
92 #ifdef PV_USE_MPI
93  for (auto &d : mDatatypes) {
94  MPI_Type_free(&d);
95  }
96  mDatatypes.clear();
97 #endif // PV_USE_MPI
98 }
99 
100 void BorderExchange::initNeighbors() {
101  neighbors.resize(NUM_NEIGHBORHOOD, -1);
102  mNumNeighbors = 0U;
103 
104  // initialize neighbor and border lists
105  // (local borders and remote neighbors form the complete neighborhood)
106 
107  for (int i = 0; i < NUM_NEIGHBORHOOD; i++) {
108  int n = neighborIndex(mMPIBlock->getRank(), i);
109  if (n >= 0) {
110  neighbors[i] = n;
111  mNumNeighbors++;
112  }
113  else {
114  neighbors[i] = mMPIBlock->getRank();
115  }
116  }
117 }
118 
119 void BorderExchange::exchange(float *data, std::vector<MPI_Request> &req) {
120 #ifdef PV_USE_MPI
121  PVHalo const &halo = mLayerLoc.halo;
122  if (halo.lt == 0 && halo.rt == 0 && halo.dn == 0 && halo.up == 0) {
123  return;
124  }
125 
126  req.clear();
127  // Start at n=1 because n=0 is the interior
128  for (int n = 1; n < NUM_NEIGHBORHOOD; n++) {
129  if (neighbors[n] == mMPIBlock->getRank())
130  continue; // don't send interior/self
131  float *recvBuf = data + recvOffset(n);
132  auto sz = req.size();
133  req.resize(sz + 1);
134  int revDir = reverseDirection(mMPIBlock->getRank(), n);
135  MPI_Irecv(
136  recvBuf,
137  1,
138  mDatatypes[n],
139  neighbors[n],
140  exchangeCounter * 16 + mTags[revDir],
141  mMPIBlock->getComm(),
142  &(req.data())[sz]);
143  }
144 
145  for (int n = 1; n < NUM_NEIGHBORHOOD; n++) {
146  if (neighbors[n] == mMPIBlock->getRank())
147  continue; // don't send interior/self
148  float *sendBuf = data + sendOffset(n);
149  auto sz = req.size();
150  req.resize(sz + 1);
151  MPI_Isend(
152  sendBuf,
153  1,
154  mDatatypes[n],
155  neighbors[n],
156  exchangeCounter * 16 + mTags[n],
157  mMPIBlock->getComm(),
158  &(req.data())[sz]);
159  }
160 
161  exchangeCounter = (exchangeCounter == 2047) ? 1024 : exchangeCounter + 1;
162 
163 // don't recv interior
164 #endif // PV_USE_MPI
165 }
166 
167 int BorderExchange::wait(std::vector<MPI_Request> &req) {
168  int status = MPI_Waitall(req.size(), req.data(), MPI_STATUSES_IGNORE);
169  req.clear();
170  return status;
171 }
172 
177 int BorderExchange::neighborIndex(int commId, int direction) {
178  int batchDimension = mMPIBlock->getBatchDimension();
179  int numRows = mMPIBlock->getNumRows();
180  int numColumns = mMPIBlock->getNumColumns();
181  int rankRowColumn = commId % (numRows * numColumns);
182  int row = rowFromRank(rankRowColumn, numRows, numColumns);
183  int column = columnFromRank(rankRowColumn, numRows, numColumns);
184  int neighborRank;
185  switch (direction) {
186  case LOCAL: return commId;
187  case NORTHWEST: neighborRank = northwest(row, column, numRows, numColumns); break;
188  case NORTH: neighborRank = north(row, column, numRows, numColumns); break;
189  case NORTHEAST: neighborRank = northeast(row, column, numRows, numColumns); break;
190  case WEST: neighborRank = west(row, column, numRows, numColumns); break;
191  case EAST: neighborRank = east(row, column, numRows, numColumns); break;
192  case SOUTHWEST: neighborRank = southwest(row, column, numRows, numColumns); break;
193  case SOUTH: neighborRank = south(row, column, numRows, numColumns); break;
194  case SOUTHEAST: neighborRank = southeast(row, column, numRows, numColumns); break;
195  default: pvAssert(0); break;
196  }
197  if (neighborRank >= 0) {
198  int rankBatchStart = commId - rankRowColumn;
199  neighborRank += rankBatchStart;
200  }
201  return neighborRank;
202 }
203 
204 int BorderExchange::northwest(int row, int column, int numRows, int numColumns) {
205  if (hasNorthwesternNeighbor(row, column, numRows, numColumns)) {
206  int const neighborRow = row - (row > 0);
207  int const neighborColumn = column - (column > 0);
208  return rankFromRowAndColumn(neighborRow, neighborColumn, numRows, numColumns);
209  }
210  else {
211  return -1;
212  }
213 }
214 
215 int BorderExchange::north(int row, int column, int numRows, int numColumns) {
216  if (hasNorthernNeighbor(row, column, numRows, numColumns)) {
217  int const neighborRow = row - 1;
218  int const neighborColumn = column;
219  return rankFromRowAndColumn(neighborRow, neighborColumn, numRows, numColumns);
220  }
221  else {
222  return -1;
223  }
224 }
225 
226 int BorderExchange::northeast(int row, int column, int numRows, int numColumns) {
227  if (hasNortheasternNeighbor(row, column, numRows, numColumns)) {
228  int const neighborRow = row - (row > 0);
229  int const neighborColumn = column + (column < numColumns - 1);
230  return rankFromRowAndColumn(neighborRow, neighborColumn, numRows, numColumns);
231  }
232  else {
233  return -1;
234  }
235 }
236 
237 int BorderExchange::west(int row, int column, int numRows, int numColumns) {
238  if (hasWesternNeighbor(row, column, numRows, numColumns)) {
239  int const neighborRow = row;
240  int const neighborColumn = column - 1;
241  return rankFromRowAndColumn(neighborRow, neighborColumn, numRows, numColumns);
242  }
243  else {
244  return -1;
245  }
246 }
247 
248 int BorderExchange::east(int row, int column, int numRows, int numColumns) {
249  if (hasEasternNeighbor(row, column, numRows, numColumns)) {
250  int const neighborRow = row;
251  int const neighborColumn = column + 1;
252  return rankFromRowAndColumn(neighborRow, neighborColumn, numRows, numColumns);
253  }
254  else {
255  return -1;
256  }
257 }
258 
259 int BorderExchange::southwest(int row, int column, int numRows, int numColumns) {
260  if (hasSouthwesternNeighbor(row, column, numRows, numColumns)) {
261  int const neighborRow = row + (row < numRows - 1);
262  int const neighborColumn = column - (column > 0);
263  return rankFromRowAndColumn(neighborRow, neighborColumn, numRows, numColumns);
264  }
265  else {
266  return -1;
267  }
268 }
269 
270 int BorderExchange::south(int row, int column, int numRows, int numColumns) {
271  if (hasSouthernNeighbor(row, column, numRows, numColumns)) {
272  int const neighborRow = row + 1;
273  int const neighborColumn = column;
274  return rankFromRowAndColumn(neighborRow, neighborColumn, numRows, numColumns);
275  }
276  else {
277  return -1;
278  }
279 }
280 
281 int BorderExchange::southeast(int row, int column, int numRows, int numColumns) {
282  if (hasSoutheasternNeighbor(row, column, numRows, numColumns)) {
283  int const neighborRow = row + (row < numRows - 1);
284  int const neighborColumn = column + (column < numColumns - 1);
285  return rankFromRowAndColumn(neighborRow, neighborColumn, numRows, numColumns);
286  }
287  else {
288  return -1;
289  }
290 }
291 
292 bool BorderExchange::hasNorthwesternNeighbor(int row, int column, int numRows, int numColumns) {
293  return (row > 0) || (column > 0);
294 }
295 
296 bool BorderExchange::hasNorthernNeighbor(int row, int column, int numRows, int numColumns) {
297  return row > 0;
298 }
299 
300 bool BorderExchange::hasNortheasternNeighbor(int row, int column, int numRows, int numColumns) {
301  return (row > 0) || (column < numColumns - 1);
302 }
303 
304 bool BorderExchange::hasWesternNeighbor(int row, int column, int numRows, int numColumns) {
305  return column > 0;
306 }
307 
308 bool BorderExchange::hasEasternNeighbor(int row, int column, int numRows, int numColumns) {
309  return column < numColumns - 1;
310 }
311 
312 bool BorderExchange::hasSouthwesternNeighbor(int row, int column, int numRows, int numColumns) {
313  return (row < numRows - 1) || (column > 0);
314 }
315 
316 bool BorderExchange::hasSouthernNeighbor(int row, int column, int numRows, int numColumns) {
317  return row < numRows - 1;
318 }
319 
320 bool BorderExchange::hasSoutheasternNeighbor(int row, int column, int numRows, int numColumns) {
321  return (row < numRows - 1) || (column < numColumns - 1);
322 }
323 
324 int BorderExchange::reverseDirection(int commId, int direction) {
325  int neighbor = neighborIndex(commId, direction);
326  if (neighbor == commId) {
327  return LOCAL;
328  }
329  int revdir = 9 - direction; // Correct unless at an edge of the MPI quilt
330  int batchDimension = mMPIBlock->getBatchDimension();
331  int numRows = mMPIBlock->getNumRows();
332  int numCols = mMPIBlock->getNumColumns();
333  int rankRowColumn = commId % (numRows * numCols);
334  int row = rowFromRank(rankRowColumn, numRows, numCols);
335  int col = columnFromRank(rankRowColumn, numRows, numCols);
336  switch (direction) {
337  case LOCAL:
338  pvAssert(0); // Should have neighbor==commId, so should have already returned
339  break;
340  case NORTHWEST:
341  pvAssert(revdir == SOUTHEAST);
342  if (row == 0) {
343  pvAssert(col > 0);
344  revdir = NORTHEAST;
345  }
346  if (col == 0) {
347  pvAssert(row > 0);
348  revdir = SOUTHWEST;
349  }
350  break;
351  case NORTH:
352  pvAssert(row > 0); // If row==0, there is no north neighbor so
353  // should have already returned.
354  break;
355  case NORTHEAST:
356  pvAssert(revdir == SOUTHWEST);
357  if (row == 0) {
358  pvAssert(col < numCols - 1);
359  revdir = NORTHWEST;
360  }
361  if (col == numCols - 1) {
362  pvAssert(row > 0);
363  revdir = SOUTHEAST;
364  }
365  break;
366  case WEST: pvAssert(col > 0); break;
367  case EAST: pvAssert(col < numCols - 1); break;
368  case SOUTHWEST:
369  pvAssert(revdir == NORTHEAST);
370  if (row == numRows - 1) {
371  pvAssert(col > 0);
372  revdir = SOUTHEAST;
373  }
374  if (col == 0) {
375  pvAssert(row < numRows - 1);
376  revdir = NORTHWEST;
377  }
378  break;
379  case SOUTH: pvAssert(row < numRows - 1); break;
380  case SOUTHEAST:
381  pvAssert(revdir == NORTHWEST);
382  if (row == numRows - 1) {
383  pvAssert(col < numCols - 1);
384  revdir = SOUTHWEST;
385  }
386  if (col == numCols - 1) {
387  pvAssert(row < numRows - 1);
388  revdir = NORTHEAST;
389  }
390  break;
391  default:
392  pvAssert(0); // All possible directions handled in above cases.
393  break;
394  }
395  return revdir;
396 }
397 
398 std::size_t BorderExchange::recvOffset(int direction) {
399  // This check should make sure n is a local rank
400  const int nx = mLayerLoc.nx;
401  const int ny = mLayerLoc.ny;
402  const int leftBorder = mLayerLoc.halo.lt;
403  const int topBorder = mLayerLoc.halo.dn;
404 
405  const int sx = strideXExtended(&mLayerLoc);
406  const int sy = strideYExtended(&mLayerLoc);
407 
408  int offset;
409 
410  switch (direction) {
411  case LOCAL: offset = sx * leftBorder + sy * topBorder; break;
412  case NORTHWEST: offset = 0; break;
413  case NORTH: offset = sx * leftBorder; break;
414  case NORTHEAST: offset = sx * leftBorder + sx * nx; break;
415  case WEST: offset = sy * topBorder; break;
416  case EAST: offset = sx * leftBorder + sx * nx + sy * topBorder; break;
417  case SOUTHWEST: offset = sy * (topBorder + ny); break;
418  case SOUTH: offset = sx * leftBorder + sy * (topBorder + ny); break;
419  case SOUTHEAST: offset = sx * leftBorder + sx * nx + sy * (topBorder + ny); break;
420  default:
421  pvAssert(0); /* All allowable directions handled in above cases */
422  break;
423  }
424  return (std::size_t)offset;
425 }
426 
431 std::size_t BorderExchange::sendOffset(int direction) {
432  const size_t nx = mLayerLoc.nx;
433  const size_t ny = mLayerLoc.ny;
434  const size_t leftBorder = mLayerLoc.halo.lt;
435  const size_t topBorder = mLayerLoc.halo.up;
436 
437  const size_t sx = strideXExtended(&mLayerLoc);
438  const size_t sy = strideYExtended(&mLayerLoc);
439 
440  const int numRows = mMPIBlock->getNumRows();
441  const int numCols = mMPIBlock->getNumColumns();
442  const int row = mMPIBlock->getRowIndex();
443  const int col = mMPIBlock->getColumnIndex();
444 
445  bool hasNorthNeighbor = row > 0;
446  bool hasWestNeighbor = col > 0;
447  bool hasEastNeighbor = col < numCols - 1;
448  bool hasSouthNeighbor = row < numRows - 1;
449 
450  int offset;
451  switch (direction) {
452  case LOCAL: offset = sx * leftBorder + sy * topBorder; break;
453  case NORTHWEST:
454  offset = sx * hasWestNeighbor * leftBorder + sy * hasNorthNeighbor * topBorder;
455  break;
456  case NORTH: offset = sx * leftBorder + sy * topBorder; break;
457  case NORTHEAST:
458  offset = sx * (nx + !hasEastNeighbor * leftBorder) + sy * hasNorthNeighbor * topBorder;
459  break;
460  case WEST: offset = sx * leftBorder + sy * topBorder; break;
461  case EAST: offset = sx * nx + sy * topBorder; break;
462  case SOUTHWEST:
463  offset = sx * hasWestNeighbor * leftBorder + sy * (ny + !hasSouthNeighbor * topBorder);
464  break;
465  case SOUTH: offset = sx * leftBorder + sy * ny; break;
466  case SOUTHEAST:
467  return sx * (nx + !hasEastNeighbor * leftBorder)
468  + sy * (ny + !hasSouthNeighbor * topBorder);
469  break;
470  default:
471  pvAssert(0); /* All allowable directions handled in above cases */
472  break;
473  }
474  return (std::size_t)offset;
475 }
476 
477 // NW and SE corners have tag 33; edges have tag 34; NE and SW have tag 35.
478 // In the top row of processes in the hypercolumn, a process is both the
479 // northeast and east neighbor of the process to its left. If there is only
480 // one row, a process is the northeast, east, and southeast neighbor of the
481 // process to its left. The numbering of tags ensures that the
482 // MPI_Send/MPI_Irecv calls between a pair of processes can be distinguished.
483 std::vector<int> const BorderExchange::mTags = {0, 33, 34, 35, 34, 34, 35, 34, 33};
484 
485 int BorderExchange::exchangeCounter = 1024;
486 
487 } // end namespace PV
int getNumColumns() const
Definition: MPIBlock.hpp:130
int getNumRows() const
Definition: MPIBlock.hpp:125
int reverseDirection(int commId, int direction)
std::size_t recvOffset(int direction)
MPI_Comm getComm() const
Definition: MPIBlock.hpp:90
std::size_t sendOffset(int direction)
int getRank() const
Definition: MPIBlock.hpp:100
int neighborIndex(int commId, int direction)
int getColumnIndex() const
Definition: MPIBlock.hpp:166
int getBatchDimension() const
Definition: MPIBlock.hpp:135
int getRowIndex() const
Definition: MPIBlock.hpp:161