MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_ClassicalDropping.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_CLASSICALDROPPING_HPP
11#define MUELU_CLASSICALDROPPING_HPP
12
14#include "Kokkos_Core.hpp"
15#include "KokkosKernels_ArithTraits.hpp"
16#include "Xpetra_Matrix.hpp"
17#include "MueLu_Utilities.hpp"
18
20
45template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node, Misc::StrengthMeasure measure>
47 public:
48 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
49 using diag_vec_type = Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
50 using local_matrix_type = typename matrix_type::local_matrix_device_type;
51 using scalar_type = typename local_matrix_type::value_type;
52 using local_ordinal_type = typename local_matrix_type::ordinal_type;
53 using memory_space = typename local_matrix_type::memory_space;
54 using diag_view_type = typename Kokkos::DualView<const scalar_type*, Kokkos::LayoutStride, typename Node::device_type, Kokkos::MemoryUnmanaged>::t_dev;
55
56 using results_view = Kokkos::View<DecisionType*, memory_space>;
57
58 using ATS = KokkosKernels::ArithTraits<scalar_type>;
59 using magnitudeType = typename ATS::magnitudeType;
60 using mATS = KokkosKernels::ArithTraits<magnitudeType>;
61 using boundary_nodes_view = Kokkos::View<const bool*, memory_space>;
62
63 private:
65 Teuchos::RCP<diag_vec_type> diagVec;
66 diag_view_type diag; // corresponds to overlapped diagonal
69
70 public:
72 : A(A_.getLocalMatrixDevice())
73 , eps(threshold)
74 , results(results_) {
75 // Construct ghosted matrix diagonal
76 if constexpr ((measure == Misc::SmoothedAggregationMeasure) || (measure == Misc::SignedSmoothedAggregationMeasure)) {
78 auto lclDiag2d = diagVec->getLocalViewDevice(Tpetra::Access::ReadOnly);
79 diag = Kokkos::subview(lclDiag2d, Kokkos::ALL(), 0);
80 } else if constexpr (measure == Misc::SignedRugeStuebenMeasure) {
82 auto lcl2d = diagVec->getLocalViewDevice(Tpetra::Access::ReadOnly);
83 diag = Kokkos::subview(lcl2d, Kokkos::ALL(), 0);
84 }
85 }
86
87 KOKKOS_FORCEINLINE_FUNCTION
88 void operator()(const local_ordinal_type rlid) const {
89 auto row = A.rowConst(rlid);
90 size_t offset = A.graph.row_map(rlid);
91
92#ifdef MUELU_COALESCE_DROP_DEBUG
93 {
94 Kokkos::printf("SoC: ");
95 for (local_ordinal_type k = 0; k < row.length; ++k) {
96 auto clid = row.colidx(k);
97
98 auto val = row.value(k);
99
100 if constexpr (measure == Misc::SmoothedAggregationMeasure) {
101 auto aiiajj = ATS::magnitude(diag(rlid)) * ATS::magnitude(diag(clid)); // |a_ii|*|a_jj|
102 auto aij2 = ATS::magnitude(val) * ATS::magnitude(val); // |a_ij|^2
103
104 Kokkos::printf("%5f ", ATS::sqrt(aij2 / aiiajj));
105 } else if constexpr (measure == Misc::SignedRugeStuebenMeasure) {
106 auto neg_aij = -ATS::real(val);
107 auto max_neg_aik = eps * ATS::real(diag(rlid));
108
109 Kokkos::printf("%5f ", neg_aij / max_neg_aik);
110 } else if constexpr (measure == Misc::SignedSmoothedAggregationMeasure) {
111 auto aiiajj = ATS::magnitude(diag(rlid)) * ATS::magnitude(diag(clid)); // |a_ii|*|a_jj|
112 const bool is_nonpositive = ATS::real(val) <= mATS::zero();
113 magnitudeType aij2 = ATS::magnitude(val) * ATS::magnitude(val); // |a_ij|^2
114 // + |a_ij|^2, if a_ij < 0, - |a_ij|^2 if a_ij >=0
115 if (!is_nonpositive)
116 aij2 = -aij2;
117 Kokkos::printf("%5f ", ATS::sqrt(aij2 / aiiajj));
118 }
119 }
120 Kokkos::printf("\n");
121 }
122#endif
123
124 for (local_ordinal_type k = 0; k < row.length; ++k) {
125 auto clid = row.colidx(k);
126
127 auto val = row.value(k);
128
129 if constexpr (measure == Misc::SmoothedAggregationMeasure) {
130 auto aiiajj = ATS::magnitude(diag(rlid)) * ATS::magnitude(diag(clid)); // |a_ii|*|a_jj|
131 auto aij2 = ATS::magnitude(val) * ATS::magnitude(val); // |a_ij|^2
132
133 results(offset + k) = Kokkos::max((aij2 <= eps * eps * aiiajj) ? DROP : KEEP,
134 results(offset + k));
135 } else if constexpr (measure == Misc::SignedRugeStuebenMeasure) {
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,
139 results(offset + k));
140 } else if constexpr (measure == Misc::SignedSmoothedAggregationMeasure) {
141 auto aiiajj = ATS::magnitude(diag(rlid)) * ATS::magnitude(diag(clid)); // |a_ii|*|a_jj|
142 const bool is_nonpositive = ATS::real(val) <= mATS::zero();
143 magnitudeType aij2 = ATS::magnitude(val) * ATS::magnitude(val); // |a_ij|^2
144 // + |a_ij|^2, if a_ij < 0, - |a_ij|^2 if a_ij >=0
145 if (!is_nonpositive)
146 aij2 = -aij2;
147 results(offset + k) = Kokkos::max((aij2 <= eps * eps * aiiajj) ? DROP : KEEP,
148 results(offset + k));
149 }
150 }
151 }
152};
153
154// helper function to allow partial template deduction
155template <Misc::StrengthMeasure measure, class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
156auto make_drop_functor(Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& A_,
159 auto functor = DropFunctor<Scalar, LocalOrdinal, GlobalOrdinal, Node, measure>(A_, threshold, results_);
160 return functor;
161}
162
163} // namespace MueLu::ClassicalDropping
164
165#endif
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 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
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_)