5 #include "Teuchos_RCP.hpp"
6 #include "Teuchos_ParameterList.hpp"
7 #include "Teuchos_TestForException.hpp"
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"
24 template <
typename CrsMatrixType>
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;
50 template <
typename MultiVectorType>
55 template <
typename MultiVectorType>
60 template <
typename MultiVectorType>
63 template <
typename MultiVectorType>
69 template <
typename MultiVectorType>
72 template <
typename MultiVectorType>
109 template <
typename SC,
typename LO,
typename GO,
typename NO>
113 typedef Tpetra::BlockCrsMatrix<SC, LO, GO, NO>
matrix_t;
115 typedef typename matrix_t::crs_graph_type
graph_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;
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;
191 template <
typename CrsMatrixType>
193 const Teuchos::RCP<matrix_t> &matrix_
196 graph(matrix->getCrsGraph()),
198 list_of_colors_host(),
205 template <
typename CrsMatrixType>
208 Teuchos::ParameterList &coloring_params
211 const std::string library = coloring_params.get(
"library",
"zoltan");
213 if (library ==
"zoltan") {
217 num_colors, list_of_colors_host, list_of_colors);
221 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::logic_error,
222 "Zoltan2CrsColorer not yet ready; use parameter library = zoltan");
226 num_colors, list_of_colors_host, list_of_colors);
231 template <
typename CrsMatrixType>
232 template <
typename MultiVectorType>
236 MultiVectorType V_fitted(graph->getColMap(), num_colors);
238 computeSeedMatrixFitted(V_fitted);
240 Tpetra::Import<lno_t, gno_t, node_t>
241 importer(graph->getColMap(), V.getMap());
243 V.doImport(V_fitted, importer, Tpetra::INSERT);
247 template <
typename CrsMatrixType>
248 template <
typename MultiVectorType>
251 MultiVectorType &V)
const
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"
258 "You must call the non-fitted version of this function.");
260 auto V_view_dev = V.getLocalViewDevice(Tpetra::Access::OverwriteAll);
261 const size_t num_local_cols = graph->getNodeNumCols();
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()");
273 template <
typename CrsMatrixType>
274 template <
typename MultiVectorType>
278 reconstructMatrix(W, *matrix);
282 template <
typename CrsMatrixType>
283 template <
typename MultiVectorType>
289 MultiVectorType W_fitted(graph->getRowMap(), num_colors);
291 Tpetra::Import<lno_t, gno_t, node_t>
292 importer(W.getMap(), graph->getRowMap());
294 W_fitted.doImport(W, importer, Tpetra::INSERT);
296 reconstructMatrixFitted(W_fitted, mat);
300 template <
typename CrsMatrixType>
301 template <
typename MultiVectorType>
304 MultiVectorType &W)
const
306 reconstructMatrixFitted(W, *matrix);
310 template <
typename CrsMatrixType>
311 template <
typename MultiVectorType>
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.");
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();
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++)
336 const size_t col = local_graph.entries(entry);
337 local_matrix.values(entry) = W_view_dev(row,my_list_of_colors[col]-1);
340 "TpetraCrsColorer::reconstructMatrixFitted()");
344 template <
typename SC,
typename LO,
typename GO,
typename NO>
346 const Teuchos::RCP<matrix_t> &matrix_
349 , graph(
Teuchos::rcp(&(matrix->getCrsGraph()), false))
351 , list_of_colors_host()
357 template <
typename SC,
typename LO,
typename GO,
typename NO>
362 multivector_t block_V_fitted(*(graph->getColMap()), matrix->getBlockSize(),
363 num_colors * matrix->getBlockSize());
365 computeSeedMatrixFitted(block_V_fitted);
367 const lno_t block_size = block_V.getBlockSize();
368 auto col_point_map = multivector_t::makePointMap(*graph->getColMap(),
370 auto blockV_point_map = multivector_t::makePointMap(*block_V.getMap(),
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));
377 Tpetra::Import<LO, GO, NO> importer_point(col_point_map_rcp,
378 blockV_point_map_rcp);
380 block_V.getMultiVectorView().doImport(block_V_fitted.getMultiVectorView(),
381 importer_point, Tpetra::INSERT);
388 template <
typename SC,
typename LO,
typename GO,
typename NO>
394 TEUCHOS_TEST_FOR_EXCEPTION(
395 !(graph->getColMap()->isLocallyFitted(*(blockV.getMap()))),
397 "Map for supplied vector is not locally fitted to the column map of the"
399 "You must call the non-fitted version of this function.");
401 const lno_t block_size = blockV.getBlockSize();
402 auto V = blockV.getMultiVectorView();
405 auto V_view_dev = V.getLocalViewDevice();
406 const size_t num_local_cols = V_view_dev.extent(0) / block_size;
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) =
416 "TpetraCrsColorer::computeSeedMatrixFitted()");
421 template <
typename SC,
typename LO,
typename GO,
typename NO>
426 reconstructMatrix(W, *matrix);
430 template <
typename SC,
typename LO,
typename GO,
typename NO>
436 multivector_t block_W_fitted(*(graph->getRowMap()), matrix->getBlockSize(),
437 num_colors * matrix->getBlockSize());
439 const lno_t block_size = block_W.getBlockSize();
440 auto row_point_map = multivector_t::makePointMap(*graph->getRowMap(),
442 auto blockW_point_map = multivector_t::makePointMap(*block_W.getMap(),
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));
449 Tpetra::Import<lno_t, gno_t, node_t>
450 importer_point(blockW_point_map_rcp, row_point_map_rcp);
452 block_W_fitted.getMultiVectorView().doImport(block_W.getMultiVectorView(),
453 importer_point, Tpetra::INSERT);
454 reconstructMatrixFitted(block_W_fitted, mat);
458 template <
typename SC,
typename LO,
typename GO,
typename NO>
463 reconstructMatrixFitted(W, *matrix);
467 template <
typename SC,
typename LO,
typename GO,
typename NO>
474 TEUCHOS_TEST_FOR_EXCEPTION(
475 !(block_W.getMap()->isLocallyFitted(*(graph->getRowMap()))),
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.");
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;
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();
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;
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)
505 const size_t row = block_row * block_size + i;
506 for (
lno_t j = 0; j < block_size; ++j)
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);
514 "TpetraCrsColorer::reconstructMatrix()");
matrix_t::local_ordinal_type lno_t
list_of_colors_t::HostMirror list_of_colors_host_t
node_t::device_type device_t
matrix_t::global_ordinal_type gno_t
list_of_colors_host_t list_of_colors_host
list_of_colors_t list_of_colors
int getColor(const size_t col) const
void computeColoring(Teuchos::ParameterList &coloring_params)
Teuchos::RCP< matrix_t > matrix
matrix_t::node_type node_t
Tpetra::BlockMultiVector< SC, LO, GO, NO > multivector_t
matrix_t::crs_graph_type graph_t
Tpetra::BlockCrsMatrix< SC, LO, GO, NO > matrix_t
Kokkos::View< int *, device_t > list_of_colors_t
bool checkColoring() const
Teuchos::RCP< const graph_t > graph
matrix_t::scalar_type scalar_t
device_t::execution_space execution_space
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
node_t::device_type device_t
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
matrix_t::node_type node_t
bool checkColoring() const
list_of_colors_t list_of_colors
matrix_t::scalar_type scalar_t
int getColor(const size_t col) const
Teuchos::RCP< matrix_t > matrix
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.