10#ifndef MUELU_CLASSICALDROPPING_HPP
11#define MUELU_CLASSICALDROPPING_HPP
14#include "Kokkos_Core.hpp"
15#include "Kokkos_ArithTraits.hpp"
16#include "Xpetra_Matrix.hpp"
17#include "MueLu_Utilities.hpp"
30template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
33 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
34 using diag_vec_type = Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
39 using diag_view_type =
typename Kokkos::DualView<const scalar_type*, Kokkos::LayoutStride, typename Node::device_type, Kokkos::MemoryUnmanaged>::t_dev;
43 using ATS = Kokkos::ArithTraits<scalar_type>;
55 :
A(A_.getLocalMatrixDevice())
59 auto lclDiag2d =
diagVec->getDeviceLocalView(Xpetra::Access::ReadOnly);
60 diag = Kokkos::subview(lclDiag2d, Kokkos::ALL(), 0);
63 KOKKOS_FORCEINLINE_FUNCTION
65 auto row =
A.rowConst(rlid);
66 size_t offset =
A.graph.row_map(rlid);
68 auto clid = row.colidx(k);
70 auto val = row.value(k);
71 auto aiiajj = ATS::magnitude(
diag(rlid)) * ATS::magnitude(
diag(clid));
72 auto aij2 = ATS::magnitude(val) * ATS::magnitude(val);
89template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
92 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
100 using ATS = Kokkos::ArithTraits<scalar_type>;
104 using diag_vec_type = Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
105 using diag_view_type =
typename Kokkos::DualView<const scalar_type*, Kokkos::LayoutStride, typename Node::device_type, Kokkos::MemoryUnmanaged>::t_dev;
115 :
A(A_.getLocalMatrixDevice())
119 auto lclDiag2d =
diagVec->getDeviceLocalView(Xpetra::Access::ReadOnly);
120 diag = Kokkos::subview(lclDiag2d, Kokkos::ALL(), 0);
123 KOKKOS_FORCEINLINE_FUNCTION
125 auto row =
A.rowConst(rlid);
126 size_t offset =
A.graph.row_map(rlid);
128 auto val = row.value(k);
129 auto neg_aij = -ATS::real(val);
130 auto max_neg_aik =
eps * ATS::real(
diag(rlid));
131 results(offset + k) = Kokkos::max((neg_aij <= max_neg_aik) ?
DROP :
KEEP,
146template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
149 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
150 using diag_vec_type = Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
155 using diag_view_type =
typename Kokkos::DualView<const scalar_type*, Kokkos::LayoutStride, typename Node::device_type, Kokkos::MemoryUnmanaged>::t_dev;
159 using ATS = Kokkos::ArithTraits<scalar_type>;
161 using mATS = Kokkos::ArithTraits<magnitudeType>;
172 :
A(A_.getLocalMatrixDevice())
177 auto lclDiag2d =
diagVec->getDeviceLocalView(Xpetra::Access::ReadOnly);
178 diag = Kokkos::subview(lclDiag2d, Kokkos::ALL(), 0);
181 KOKKOS_FORCEINLINE_FUNCTION
183 auto row =
A.rowConst(rlid);
184 size_t offset =
A.graph.row_map(rlid);
186 auto clid = row.colidx(k);
188 auto val = row.value(k);
189 auto aiiajj = ATS::magnitude(
diag(rlid)) * ATS::magnitude(
diag(clid));
190 const bool is_nonpositive = ATS::real(val) <= mATS::zero();
191 magnitudeType aij2 = ATS::magnitude(val) * ATS::magnitude(val);
typename matrix_type::local_matrix_type local_matrix_type
typename Kokkos::DualView< const scalar_type *, Kokkos::LayoutStride, typename Node::device_type, Kokkos::MemoryUnmanaged >::t_dev diag_view_type
typename local_matrix_type::ordinal_type local_ordinal_type
KOKKOS_FORCEINLINE_FUNCTION void operator()(const local_ordinal_type rlid) const
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > matrix_type
SAFunctor(matrix_type &A_, magnitudeType threshold, results_view &results_)
Teuchos::RCP< diag_vec_type > diagVec
typename local_matrix_type::memory_space memory_space
Kokkos::ArithTraits< scalar_type > ATS
Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > diag_vec_type
typename ATS::magnitudeType magnitudeType
typename local_matrix_type::value_type scalar_type
Kokkos::View< const bool *, memory_space > boundary_nodes_view
Kokkos::View< DecisionType *, memory_space > results_view
SignedRSFunctor(matrix_type &A_, magnitudeType threshold, results_view &results_)
typename Kokkos::DualView< const scalar_type *, Kokkos::LayoutStride, typename Node::device_type, Kokkos::MemoryUnmanaged >::t_dev diag_view_type
typename local_matrix_type::memory_space memory_space
Kokkos::View< const bool *, memory_space > boundary_nodes_view
typename local_matrix_type::ordinal_type local_ordinal_type
typename ATS::magnitudeType magnitudeType
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > matrix_type
Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > diag_vec_type
typename matrix_type::local_matrix_type local_matrix_type
Kokkos::View< DecisionType *, memory_space > results_view
Teuchos::RCP< diag_vec_type > diagVec
Kokkos::ArithTraits< scalar_type > ATS
typename local_matrix_type::value_type scalar_type
KOKKOS_FORCEINLINE_FUNCTION void operator()(const local_ordinal_type rlid) const
Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > diag_vec_type
SignedSAFunctor(matrix_type &A_, magnitudeType threshold, results_view &results_)
typename local_matrix_type::value_type scalar_type
Kokkos::View< const bool *, memory_space > boundary_nodes_view
Kokkos::ArithTraits< scalar_type > ATS
Kokkos::ArithTraits< magnitudeType > mATS
Teuchos::RCP< diag_vec_type > diagVec
KOKKOS_FORCEINLINE_FUNCTION void operator()(const local_ordinal_type rlid) const
typename local_matrix_type::ordinal_type local_ordinal_type
typename ATS::magnitudeType magnitudeType
typename matrix_type::local_matrix_type local_matrix_type
typename local_matrix_type::memory_space memory_space
Kokkos::View< DecisionType *, memory_space > results_view
typename Kokkos::DualView< const scalar_type *, Kokkos::LayoutStride, typename Node::device_type, Kokkos::MemoryUnmanaged >::t_dev diag_view_type
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > matrix_type
static Teuchos::RCP< Vector > GetMatrixMaxMinusOffDiagonal(const Xpetra::Matrix< Scalar, DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > &A)
static RCP< Vector > GetMatrixOverlappedDiagonal(const Matrix &A)