Zoltan2
Zoltan2_TpetraCrsColorer.hpp
Go to the documentation of this file.
1 #pragma once
2 
3 #include <stdexcept>
4 
5 #include "Teuchos_RCP.hpp"
6 #include "Teuchos_ParameterList.hpp"
7 #include "Teuchos_TestForException.hpp"
8 
9 #include "Tpetra_Map.hpp"
10 #include "Tpetra_MultiVector.hpp"
11 #include "Tpetra_CrsMatrix.hpp"
12 #include "Tpetra_CrsGraph.hpp"
13 #include "Tpetra_BlockMultiVector.hpp"
14 #include "Tpetra_BlockCrsMatrix.hpp"
15 
19 
20 namespace Zoltan2
21 {
22 
23 // Base class for coloring Tpetra::CrsMatrix for use in column compression
24 template <typename CrsMatrixType>
26 {
27 public:
28  typedef CrsMatrixType matrix_t;
29  typedef typename matrix_t::crs_graph_type graph_t;
30  typedef typename matrix_t::scalar_type scalar_t;
31  typedef typename matrix_t::local_ordinal_type lno_t;
32  typedef typename matrix_t::global_ordinal_type gno_t;
33  typedef typename matrix_t::node_type node_t;
34  typedef typename node_t::device_type device_t;
35  typedef typename device_t::execution_space execution_space;
36  typedef Kokkos::View<int *, device_t> list_of_colors_t;
37  typedef typename list_of_colors_t::HostMirror list_of_colors_host_t;
38 
39  // Constructor
40  TpetraCrsColorer(const Teuchos::RCP<matrix_t> &matrix_);
41 
42  // Destructor
44 
45  // Compute coloring data
46  void
47  computeColoring(Teuchos::ParameterList &coloring_params);
48 
49  // Compute seed matrix
50  template <typename MultiVectorType>
51  void
52  computeSeedMatrix(MultiVectorType &V) const;
53 
54  // Compute seed matrix with distribution fitted to the graph's column map
55  template <typename MultiVectorType>
56  void
57  computeSeedMatrixFitted(MultiVectorType &V) const;
58 
59  // Reconstruct matrix from supplied compressed vector
60  template <typename MultiVectorType>
61  void
62  reconstructMatrix(MultiVectorType &W) const;
63  template <typename MultiVectorType>
64  void
65  reconstructMatrix(MultiVectorType &W, matrix_t &mat) const;
66 
67  // Reconstruct matrix from supplied compressed vector fitted to the graph's
68  // row map
69  template <typename MultiVectorType>
70  void
71  reconstructMatrixFitted(MultiVectorType &W) const;
72  template <typename MultiVectorType>
73  void
74  reconstructMatrixFitted(MultiVectorType &W, matrix_t &mat) const;
75 
76  // Return number of colors
77  // KDD should num_colors be int or size_t?
78  int
79  getNumColors() const
80  {
81  return num_colors;
82  }
83 
84  // Get color given a column index
85  int
86  getColor(const size_t col) const
87  {
88  return list_of_colors_host(col) - 1;
89  }
90 
91  // Check coloring is valid, i.e., no row has two nonzero columns with
92  // same color
93  bool
94  checkColoring() const
95  {
97  }
98 
99 protected:
100  Teuchos::RCP<matrix_t> matrix;
101  Teuchos::RCP<const graph_t> graph;
105 };
106 
108 // Specialization of TpetraCrsColorer for BlockCrsMatrix
109 template <typename SC, typename LO, typename GO, typename NO>
110 class TpetraCrsColorer<Tpetra::BlockCrsMatrix<SC,LO,GO,NO> >
111 {
112 public:
113  typedef Tpetra::BlockCrsMatrix<SC, LO, GO, NO> matrix_t;
114  typedef Tpetra::BlockMultiVector<SC, LO, GO, NO> multivector_t;
115  typedef typename matrix_t::crs_graph_type graph_t;
116  typedef typename matrix_t::scalar_type scalar_t;
117  typedef typename matrix_t::local_ordinal_type lno_t;
118  typedef typename matrix_t::global_ordinal_type gno_t;
119  typedef typename matrix_t::node_type node_t;
120  typedef typename node_t::device_type device_t;
121  typedef typename device_t::execution_space execution_space;
122  typedef Kokkos::View<int *, device_t> list_of_colors_t;
123  typedef typename list_of_colors_t::HostMirror list_of_colors_host_t;
124 
125  // Constructor
126  TpetraCrsColorer(const Teuchos::RCP<matrix_t> &matrix_);
127 
128  // Destructor
130 
131  // Compute coloring data
132  void
133  computeColoring(Teuchos::ParameterList &coloring_params);
134 
135  // Compute seed matrix
136  void
138 
139  // Compute seed matrix with distribution fitted to the graph's column map
140  void
142 
143  // Reconstruct matrix from supplied compressed vector
144  void
146  void
147  reconstructMatrix(multivector_t &W, matrix_t &mat) const;
148 
149  // Reconstruct matrix from supplied compressed vector fitted to the graph's
150  // row map
151  void
153  void
155 
156  // Return number of colors
157  int
158  getNumColors() const
159  {
160  return num_colors * matrix->getBlockSize();
161  }
162 
163  // Get color given a column index
164  int
165  getColor(const size_t col) const
166  {
167  const int block_size = matrix->getBlockSize();
168  const size_t block_col = col / block_size;
169  const int offset = col - block_col * block_size;
170  return (list_of_colors_host(block_col) - 1) * block_size + offset;
171  }
172 
173  // Check coloring is valid, i.e., no row has two nonzero columns with
174  // same color
175  bool
177  {
179  }
180 
181 protected:
182  Teuchos::RCP<matrix_t> matrix;
183  Teuchos::RCP<const graph_t> graph;
187 };
188 
190 
191 template <typename CrsMatrixType>
193  const Teuchos::RCP<matrix_t> &matrix_
194 )
195  : matrix(matrix_),
196  graph(matrix->getCrsGraph()),
197  list_of_colors(),
198  list_of_colors_host(),
199  num_colors(0)
200 {
201 
202 }
203 
205 template <typename CrsMatrixType>
206 void
208  Teuchos::ParameterList &coloring_params
209 )
210 {
211  const std::string library = coloring_params.get("library", "zoltan");
212 
213  if (library == "zoltan") {
214  // Use Zoltan's coloring
215  ZoltanCrsColorer<matrix_t> zz(matrix);
216  zz.computeColoring(coloring_params,
217  num_colors, list_of_colors_host, list_of_colors);
218  }
219  else {
220  // Use Zoltan2's coloring when it is ready
221  TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error,
222  "Zoltan2CrsColorer not yet ready; use parameter library = zoltan");
223 
224  Zoltan2CrsColorer<matrix_t> zz2(matrix);
225  zz2.computeColoring(coloring_params,
226  num_colors, list_of_colors_host, list_of_colors);
227  }
228 }
229 
231 template <typename CrsMatrixType>
232 template <typename MultiVectorType>
233 void
235 {
236  MultiVectorType V_fitted(graph->getColMap(), num_colors);
237 
238  computeSeedMatrixFitted(V_fitted);
239 
240  Tpetra::Import<lno_t, gno_t, node_t>
241  importer(graph->getColMap(), V.getMap());
242 
243  V.doImport(V_fitted, importer, Tpetra::INSERT);
244 }
245 
247 template <typename CrsMatrixType>
248 template <typename MultiVectorType>
249 void
251  MultiVectorType &V) const
252 {
253  // Check V's map is locally fitted to the graph's column map
254  TEUCHOS_TEST_FOR_EXCEPTION(
255  !(graph->getColMap()->isLocallyFitted(*(V.getMap()))), std::domain_error,
256  "Map for supplied vector is not locally fitted to the column map of the"
257  " Jacobian graph. "
258  "You must call the non-fitted version of this function.");
259 
260  auto V_view_dev = V.getLocalViewDevice(Tpetra::Access::OverwriteAll);
261  const size_t num_local_cols = graph->getNodeNumCols();
262  list_of_colors_t my_list_of_colors = list_of_colors;
263 
264  Kokkos::parallel_for(
265  Kokkos::RangePolicy<execution_space>(0, num_local_cols),
266  KOKKOS_LAMBDA(const size_t i) {
267  V_view_dev(i, my_list_of_colors[i] - 1) = scalar_t(1.0); },
268  "TpetraCrsColorer::computeSeedMatrixFitted()");
269 
270 }
271 
273 template <typename CrsMatrixType>
274 template <typename MultiVectorType>
275 void
277 {
278  reconstructMatrix(W, *matrix);
279 }
280 
282 template <typename CrsMatrixType>
283 template <typename MultiVectorType>
284 void
286  MultiVectorType &W,
287  matrix_t &mat) const
288 {
289  MultiVectorType W_fitted(graph->getRowMap(), num_colors);
290 
291  Tpetra::Import<lno_t, gno_t, node_t>
292  importer(W.getMap(), graph->getRowMap());
293 
294  W_fitted.doImport(W, importer, Tpetra::INSERT);
295 
296  reconstructMatrixFitted(W_fitted, mat);
297 }
298 
300 template <typename CrsMatrixType>
301 template <typename MultiVectorType>
302 void
304  MultiVectorType &W) const
305 {
306  reconstructMatrixFitted(W, *matrix);
307 }
308 
310 template <typename CrsMatrixType>
311 template <typename MultiVectorType>
312 void
314  MultiVectorType &W,
315  matrix_t &mat) const
316 {
317  // Check the graph's row map is locally fitted to W's map
318  TEUCHOS_TEST_FOR_EXCEPTION(
319  !(W.getMap()->isLocallyFitted(*(graph->getRowMap()))), std::domain_error,
320  "Row map of the Jacobian graph is not locally fitted to the vector's map."
321  " You must call the non-fitted version of this function.");
322 
323  auto W_view_dev = W.getLocalViewDevice(Tpetra::Access::ReadOnly);
324  auto local_matrix = mat.getLocalMatrixDevice();
325  auto local_graph = local_matrix.graph;
326  const size_t num_local_rows = graph->getNodeNumRows();
327  list_of_colors_t my_list_of_colors = list_of_colors;
328 
329  Kokkos::parallel_for(
330  Kokkos::RangePolicy<execution_space>(0, num_local_rows),
331  KOKKOS_LAMBDA(const size_t row) {
332  const size_t entry_begin = local_graph.row_map(row);
333  const size_t entry_end = local_graph.row_map(row + 1);
334  for (size_t entry = entry_begin; entry < entry_end; entry++)
335  {
336  const size_t col = local_graph.entries(entry);
337  local_matrix.values(entry) = W_view_dev(row,my_list_of_colors[col]-1);
338  }
339  },
340  "TpetraCrsColorer::reconstructMatrixFitted()");
341 }
342 
344 template <typename SC, typename LO, typename GO, typename NO>
346  const Teuchos::RCP<matrix_t> &matrix_
347 )
348  : matrix(matrix_)
349  , graph(Teuchos::rcp(&(matrix->getCrsGraph()), false))
350  , list_of_colors()
351  , list_of_colors_host()
352  , num_colors(0)
353 {
354 }
355 
357 template <typename SC, typename LO, typename GO, typename NO>
358 void
360  multivector_t &block_V) const
361 {
362  multivector_t block_V_fitted(*(graph->getColMap()), matrix->getBlockSize(),
363  num_colors * matrix->getBlockSize());
364 
365  computeSeedMatrixFitted(block_V_fitted);
366 
367  const lno_t block_size = block_V.getBlockSize();
368  auto col_point_map = multivector_t::makePointMap(*graph->getColMap(),
369  block_size);
370  auto blockV_point_map = multivector_t::makePointMap(*block_V.getMap(),
371  block_size);
372  const auto col_point_map_rcp =
373  Teuchos::rcp(new Tpetra::Map<LO, GO, NO>(col_point_map));
374  const auto blockV_point_map_rcp =
375  Teuchos::rcp(new Tpetra::Map<LO, GO, NO>(blockV_point_map));
376 
377  Tpetra::Import<LO, GO, NO> importer_point(col_point_map_rcp,
378  blockV_point_map_rcp);
379 
380  block_V.getMultiVectorView().doImport(block_V_fitted.getMultiVectorView(),
381  importer_point, Tpetra::INSERT);
382 
383  // Tpetra::Import<LO,GO,NO> importer(graph->getColMap(), block_V.getMap());
384  // block_V.doImport(block_V_fitted, importer, Tpetra::INSERT);
385 }
386 
388 template <typename SC, typename LO, typename GO, typename NO>
389 void
391  multivector_t &blockV) const
392 {
393  // Check blockV's map is locally fitted to the graph's column map
394  TEUCHOS_TEST_FOR_EXCEPTION(
395  !(graph->getColMap()->isLocallyFitted(*(blockV.getMap()))),
396  std::domain_error,
397  "Map for supplied vector is not locally fitted to the column map of the"
398  " Jacobian graph. "
399  "You must call the non-fitted version of this function.");
400 
401  const lno_t block_size = blockV.getBlockSize();
402  auto V = blockV.getMultiVectorView();
403 
404  V.sync_device();
405  auto V_view_dev = V.getLocalViewDevice();
406  const size_t num_local_cols = V_view_dev.extent(0) / block_size;
407  list_of_colors_t my_list_of_colors = list_of_colors;
408 
409  Kokkos::parallel_for(
410  Kokkos::RangePolicy<execution_space>(0, num_local_cols),
411  KOKKOS_LAMBDA(const size_t i) {
412  for (lno_t j = 0; j < block_size; ++j)
413  V_view_dev(i*block_size+j, (my_list_of_colors[i]-1)*block_size+j) =
414  scalar_t(1.0);
415  },
416  "TpetraCrsColorer::computeSeedMatrixFitted()");
417 
418 }
419 
421 template <typename SC, typename LO, typename GO, typename NO>
422 void
424  multivector_t &W) const
425 {
426  reconstructMatrix(W, *matrix);
427 }
428 
430 template <typename SC, typename LO, typename GO, typename NO>
431 void
433  multivector_t &block_W,
434  matrix_t &mat) const
435 {
436  multivector_t block_W_fitted(*(graph->getRowMap()), matrix->getBlockSize(),
437  num_colors * matrix->getBlockSize());
438 
439  const lno_t block_size = block_W.getBlockSize();
440  auto row_point_map = multivector_t::makePointMap(*graph->getRowMap(),
441  block_size);
442  auto blockW_point_map = multivector_t::makePointMap(*block_W.getMap(),
443  block_size);
444  const auto row_point_map_rcp =
445  Teuchos::rcp(new Tpetra::Map<LO, GO, NO>(row_point_map));
446  const auto blockW_point_map_rcp =
447  Teuchos::rcp(new Tpetra::Map<LO, GO, NO>(blockW_point_map));
448 
449  Tpetra::Import<lno_t, gno_t, node_t>
450  importer_point(blockW_point_map_rcp, row_point_map_rcp);
451 
452  block_W_fitted.getMultiVectorView().doImport(block_W.getMultiVectorView(),
453  importer_point, Tpetra::INSERT);
454  reconstructMatrixFitted(block_W_fitted, mat);
455 }
456 
458 template <typename SC, typename LO, typename GO, typename NO>
459 void
461  multivector_t &W) const
462 {
463  reconstructMatrixFitted(W, *matrix);
464 }
465 
467 template <typename SC, typename LO, typename GO, typename NO>
468 void
470  multivector_t &block_W,
471  matrix_t &mat) const
472 {
473  // Check the graph's row map is locally fitted to W's map
474  TEUCHOS_TEST_FOR_EXCEPTION(
475  !(block_W.getMap()->isLocallyFitted(*(graph->getRowMap()))),
476  std::domain_error,
477  "Row map of the Jacobian graph is not locally fitted to the vector's map."
478  " You must call the non-fitted version of this function.");
479 
480  // Blocks are block_size x block_size stored with LayoutRight
481  const lno_t block_size = block_W.getBlockSize();
482  const lno_t block_stride = block_size * block_size;
483  const lno_t block_row_stride = block_size;
484 
485  auto W = block_W.getMultiVectorView();
486  auto W_view_dev = W.getLocalViewDevice();
487  auto matrix_vals = matrix->getValuesDevice();
488  auto local_graph = graph->getLocalGraph();
489  const size_t num_local_rows = graph->getNodeNumRows();
490  list_of_colors_t my_list_of_colors = list_of_colors;
491 
492  Kokkos::parallel_for(
493  Kokkos::RangePolicy<execution_space>(0, num_local_rows),
494  KOKKOS_LAMBDA(const size_t block_row) {
495  const size_t entry_begin = local_graph.row_map(block_row);
496  const size_t entry_end = local_graph.row_map(block_row + 1);
497  for (size_t block_entry = entry_begin; block_entry < entry_end;
498  block_entry++)
499  {
500  const size_t block_col = local_graph.entries(block_entry);
501  const size_t block_offset = block_stride * block_entry;
502  const int block_color = my_list_of_colors[block_col] - 1;
503  for (lno_t i = 0; i < block_size; ++i)
504  {
505  const size_t row = block_row * block_size + i;
506  for (lno_t j = 0; j < block_size; ++j)
507  {
508  const size_t entry = block_offset + block_row_stride * i + j;
509  matrix_vals(entry) = W_view_dev(row, block_color*block_size+j);
510  }
511  }
512  }
513  },
514  "TpetraCrsColorer::reconstructMatrix()");
515 }
516 } // namespace Tpetra
void computeColoring(Teuchos::ParameterList &coloring_params)
void computeSeedMatrix(MultiVectorType &V) const
list_of_colors_t::HostMirror list_of_colors_host_t
void computeSeedMatrixFitted(MultiVectorType &V) const
void reconstructMatrixFitted(MultiVectorType &W) const
void reconstructMatrix(MultiVectorType &W) const
TpetraCrsColorer(const Teuchos::RCP< matrix_t > &matrix_)
Kokkos::View< int *, device_t > list_of_colors_t
device_t::execution_space execution_space
Teuchos::RCP< const graph_t > graph
void computeColoring(Teuchos::ParameterList &coloring_params)
list_of_colors_host_t list_of_colors_host
matrix_t::crs_graph_type graph_t
matrix_t::local_ordinal_type lno_t
int getColor(const size_t col) const
matrix_t::global_ordinal_type gno_t
void computeColoring(Teuchos::ParameterList &coloring_params, int &num_colors, list_of_colors_host_t &list_of_colors_host, list_of_colors_t &list_of_colors)
void computeColoring(Teuchos::ParameterList &coloring_params, int &num_colors, list_of_colors_host_t &list_of_colors_host, list_of_colors_t &list_of_colors) const
bool check_coloring(const Tpetra::CrsGraph< LO, GO, NO > &graph, const list_of_colors_t &list_of_colors)
Created by mbenlioglu on Aug 31, 2020.