MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_ScalarDroppingBase.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// MueLu: A package for multigrid based preconditioning
4//
5// Copyright 2012 NTESS and the MueLu contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef MUELU_SCALARDROPPINGBASE_HPP
11#define MUELU_SCALARDROPPINGBASE_HPP
12
13#include "Xpetra_Matrix.hpp"
14
17#include "MueLu_Factory.hpp"
18
19namespace MueLu {
20
21template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
23 public:
24 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
25 using crs_matrix_type = Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
26 using GraphType = Xpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node>;
27 using local_matrix_type = typename crs_matrix_type::local_matrix_device_type;
28 using local_graph_type = typename GraphType::local_graph_device_type;
29 using rowptr_type = typename local_graph_type::row_map_type::non_const_type;
30 using device_type = typename Node::device_type;
31 using memory_space = typename device_type::memory_space;
32 using execution_space = typename local_matrix_type::execution_space;
33 using results_view = Kokkos::View<DecisionType*, memory_space>;
35 using range_type = Kokkos::RangePolicy<LocalOrdinal, execution_space>;
36 using LocalOrdinalVector = Xpetra::Vector<LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node>;
37
38 template <class... Functors>
39 static void runDroppingFunctorsImpl(local_matrix_type& lclA, results_view& results, rowptr_type& filtered_rowptr, LocalOrdinal& nnz_filtered, Functors&... functors) {
40 auto range = range_type(0, lclA.numRows());
41#if !defined(HAVE_MUELU_DEBUG)
42 auto countingFunctor = MatrixConstruction::PointwiseCountingFunctor(lclA, results, filtered_rowptr, functors...);
43
44#else
45 auto debug = Misc::DebugFunctor(lclA, results);
46 auto countingFunctor = MatrixConstruction::PointwiseCountingFunctor(lclA, results, filtered_rowptr, functors..., debug);
47#endif
48 Kokkos::parallel_scan("MueLu::CoalesceDrop::CountEntries", range, countingFunctor, nnz_filtered);
49 }
50
51 template <class... Functors>
52 static void runDroppingFunctors(matrix_type& A, results_view& results, rowptr_type& filtered_rowptr, LocalOrdinal& nnz_filtered, const bool useBlocking, Level& level, const Factory& factory, Functors&... functors) {
53 auto lclA = A.getLocalMatrixDevice();
54 if (useBlocking) {
55 auto BlockNumber = level.template Get<Teuchos::RCP<LocalOrdinalVector>>("BlockNumber", factory.GetFactory("BlockNumber").get());
56 auto block_diagonalize = Misc::BlockDiagonalizeFunctor(A, *BlockNumber, results);
57
58 runDroppingFunctorsImpl(lclA, results, filtered_rowptr, nnz_filtered, block_diagonalize, functors...);
59 } else
60 runDroppingFunctorsImpl(lclA, results, filtered_rowptr, nnz_filtered, functors...);
61 }
62};
63
64} // namespace MueLu
65#endif
MueLu::DefaultLocalOrdinal LocalOrdinal
const RCP< const FactoryBase > GetFactory(const std::string &varName) const
Default implementation of FactoryAcceptor::GetFactory()
Class that holds all level-specific information.
Functor that executes a sequence of sub-functors on each row for a problem with blockSize == 1.
Functor that drops all entries that are not on the block diagonal.
Functor that checks that all entries have been marked.
typename device_type::memory_space memory_space
Kokkos::View< DecisionType *, memory_space > results_view
typename local_graph_type::row_map_type::non_const_type rowptr_type
Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > crs_matrix_type
Xpetra::Vector< LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node > LocalOrdinalVector
Kokkos::RangePolicy< LocalOrdinal, execution_space > range_type
typename crs_matrix_type::local_matrix_device_type local_matrix_type
static void runDroppingFunctors(matrix_type &A, results_view &results, rowptr_type &filtered_rowptr, LocalOrdinal &nnz_filtered, const bool useBlocking, Level &level, const Factory &factory, Functors &... functors)
static void runDroppingFunctorsImpl(local_matrix_type &lclA, results_view &results, rowptr_type &filtered_rowptr, LocalOrdinal &nnz_filtered, Functors &... functors)
Xpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > GraphType
typename GraphType::local_graph_device_type local_graph_type
typename Node::device_type device_type
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > matrix_type
typename local_matrix_type::execution_space execution_space
Namespace for MueLu classes and methods.