10#ifndef MUELU_CLASSICALDROPPING_HPP
11#define MUELU_CLASSICALDROPPING_HPP
14#include "Kokkos_Core.hpp"
15#include "KokkosKernels_ArithTraits.hpp"
16#include "Xpetra_Matrix.hpp"
17#include "MueLu_Utilities.hpp"
45template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node, Misc::StrengthMeasure measure>
48 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
49 using diag_vec_type = Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
54 using diag_view_type =
typename Kokkos::DualView<const scalar_type*, Kokkos::LayoutStride, typename Node::device_type, Kokkos::MemoryUnmanaged>::t_dev;
58 using ATS = KokkosKernels::ArithTraits<scalar_type>;
60 using mATS = KokkosKernels::ArithTraits<magnitudeType>;
72 :
A(A_.getLocalMatrixDevice())
78 auto lclDiag2d =
diagVec->getLocalViewDevice(Tpetra::Access::ReadOnly);
79 diag = Kokkos::subview(lclDiag2d, Kokkos::ALL(), 0);
82 auto lcl2d =
diagVec->getLocalViewDevice(Tpetra::Access::ReadOnly);
83 diag = Kokkos::subview(lcl2d, Kokkos::ALL(), 0);
87 KOKKOS_FORCEINLINE_FUNCTION
89 auto row =
A.rowConst(rlid);
90 size_t offset =
A.graph.row_map(rlid);
92#ifdef MUELU_COALESCE_DROP_DEBUG
94 Kokkos::printf(
"SoC: ");
96 auto clid = row.colidx(k);
98 auto val = row.value(k);
101 auto aiiajj = ATS::magnitude(
diag(rlid)) * ATS::magnitude(
diag(clid));
102 auto aij2 = ATS::magnitude(val) * ATS::magnitude(val);
104 Kokkos::printf(
"%5f ", ATS::sqrt(aij2 / aiiajj));
106 auto neg_aij = -ATS::real(val);
107 auto max_neg_aik =
eps * ATS::real(
diag(rlid));
109 Kokkos::printf(
"%5f ", neg_aij / max_neg_aik);
111 auto aiiajj = ATS::magnitude(
diag(rlid)) * ATS::magnitude(
diag(clid));
112 const bool is_nonpositive = ATS::real(val) <= mATS::zero();
113 magnitudeType aij2 = ATS::magnitude(val) * ATS::magnitude(val);
117 Kokkos::printf(
"%5f ", ATS::sqrt(aij2 / aiiajj));
120 Kokkos::printf(
"\n");
125 auto clid = row.colidx(k);
127 auto val = row.value(k);
130 auto aiiajj = ATS::magnitude(
diag(rlid)) * ATS::magnitude(
diag(clid));
131 auto aij2 = ATS::magnitude(val) * ATS::magnitude(val);
136 auto neg_aij = -ATS::real(val);
137 auto max_neg_aik =
eps * ATS::real(
diag(rlid));
138 results(offset + k) = Kokkos::max((neg_aij < max_neg_aik) ?
DROP :
KEEP,
141 auto aiiajj = ATS::magnitude(
diag(rlid)) * ATS::magnitude(
diag(clid));
142 const bool is_nonpositive = ATS::real(val) <= mATS::zero();
143 magnitudeType aij2 = ATS::magnitude(val) * ATS::magnitude(val);
155template <Misc::StrengthMeasure measure,
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
Classical dropping criterion.
typename matrix_type::local_matrix_device_type local_matrix_type
DropFunctor(matrix_type &A_, magnitudeType threshold, results_view &results_)
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > matrix_type
typename ATS::magnitudeType magnitudeType
typename local_matrix_type::value_type scalar_type
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
typename local_matrix_type::memory_space memory_space
typename local_matrix_type::ordinal_type local_ordinal_type
KokkosKernels::ArithTraits< scalar_type > ATS
KokkosKernels::ArithTraits< magnitudeType > mATS
Kokkos::View< const bool *, memory_space > boundary_nodes_view
KOKKOS_FORCEINLINE_FUNCTION void operator()(const local_ordinal_type rlid) const
Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > diag_vec_type
Teuchos::RCP< diag_vec_type > diagVec
static Teuchos::RCP< Vector > GetMatrixMaxMinusOffDiagonal(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
Return vector containing: max_{i\not=k}(-a_ik), for each for i in the matrix.
static RCP< Vector > GetMatrixOverlappedDiagonal(const Matrix &A)
Extract Overlapped Matrix Diagonal.
auto make_drop_functor(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A_, typename DropFunctor< Scalar, LocalOrdinal, GlobalOrdinal, Node, measure >::magnitudeType threshold, typename DropFunctor< Scalar, LocalOrdinal, GlobalOrdinal, Node, measure >::results_view &results_)
@ SignedSmoothedAggregationMeasure
@ SignedRugeStuebenMeasure
@ SmoothedAggregationMeasure