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 "Kokkos_ArithTraits.hpp"
16#include "Xpetra_Matrix.hpp"
17#include "MueLu_Utilities.hpp"
18
20
30template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
31class SAFunctor {
32 private:
33 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
34 using diag_vec_type = Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
35 using local_matrix_type = typename matrix_type::local_matrix_type;
36 using scalar_type = typename local_matrix_type::value_type;
37 using local_ordinal_type = typename local_matrix_type::ordinal_type;
38 using memory_space = typename local_matrix_type::memory_space;
39 using diag_view_type = typename Kokkos::DualView<const scalar_type*, Kokkos::LayoutStride, typename Node::device_type, Kokkos::MemoryUnmanaged>::t_dev;
40
41 using results_view = Kokkos::View<DecisionType*, memory_space>;
42
43 using ATS = Kokkos::ArithTraits<scalar_type>;
44 using magnitudeType = typename ATS::magnitudeType;
45 using boundary_nodes_view = Kokkos::View<const bool*, memory_space>;
46
48 Teuchos::RCP<diag_vec_type> diagVec;
49 diag_view_type diag; // corresponds to overlapped diagonal
52
53 public:
55 : A(A_.getLocalMatrixDevice())
56 , eps(threshold)
57 , results(results_) {
59 auto lclDiag2d = diagVec->getDeviceLocalView(Xpetra::Access::ReadOnly);
60 diag = Kokkos::subview(lclDiag2d, Kokkos::ALL(), 0);
61 }
62
63 KOKKOS_FORCEINLINE_FUNCTION
64 void operator()(const local_ordinal_type rlid) const {
65 auto row = A.rowConst(rlid);
66 size_t offset = A.graph.row_map(rlid);
67 for (local_ordinal_type k = 0; k < row.length; ++k) {
68 auto clid = row.colidx(k);
69
70 auto val = row.value(k);
71 auto aiiajj = ATS::magnitude(diag(rlid)) * ATS::magnitude(diag(clid)); // |a_ii|*|a_jj|
72 auto aij2 = ATS::magnitude(val) * ATS::magnitude(val); // |a_ij|^2
73
74 results(offset + k) = Kokkos::max((aij2 <= eps * eps * aiiajj) ? DROP : KEEP,
75 results(offset + k));
76 }
77 }
78};
79
89template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
91 private:
92 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
93 using local_matrix_type = typename matrix_type::local_matrix_type;
94 using scalar_type = typename local_matrix_type::value_type;
95 using local_ordinal_type = typename local_matrix_type::ordinal_type;
96 using memory_space = typename local_matrix_type::memory_space;
97
98 using results_view = Kokkos::View<DecisionType*, memory_space>;
99
100 using ATS = Kokkos::ArithTraits<scalar_type>;
101 using magnitudeType = typename ATS::magnitudeType;
102 using boundary_nodes_view = Kokkos::View<const bool*, memory_space>;
103
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;
106
108 Teuchos::RCP<diag_vec_type> diagVec;
109 diag_view_type diag; // corresponds to overlapped diagonal
112
113 public:
115 : A(A_.getLocalMatrixDevice())
116 , eps(threshold)
117 , results(results_) {
119 auto lclDiag2d = diagVec->getDeviceLocalView(Xpetra::Access::ReadOnly);
120 diag = Kokkos::subview(lclDiag2d, Kokkos::ALL(), 0);
121 }
122
123 KOKKOS_FORCEINLINE_FUNCTION
124 void operator()(const local_ordinal_type rlid) const {
125 auto row = A.rowConst(rlid);
126 size_t offset = A.graph.row_map(rlid);
127 for (local_ordinal_type k = 0; k < row.length; ++k) {
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,
132 results(offset + k));
133 }
134 }
135};
136
146template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
148 private:
149 using matrix_type = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
150 using diag_vec_type = Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
151 using local_matrix_type = typename matrix_type::local_matrix_type;
152 using scalar_type = typename local_matrix_type::value_type;
153 using local_ordinal_type = typename local_matrix_type::ordinal_type;
154 using memory_space = typename local_matrix_type::memory_space;
155 using diag_view_type = typename Kokkos::DualView<const scalar_type*, Kokkos::LayoutStride, typename Node::device_type, Kokkos::MemoryUnmanaged>::t_dev;
156
157 using results_view = Kokkos::View<DecisionType*, memory_space>;
158
159 using ATS = Kokkos::ArithTraits<scalar_type>;
160 using magnitudeType = typename ATS::magnitudeType;
161 using mATS = Kokkos::ArithTraits<magnitudeType>;
162 using boundary_nodes_view = Kokkos::View<const bool*, memory_space>;
163
165 Teuchos::RCP<diag_vec_type> diagVec;
166 diag_view_type diag; // corresponds to overlapped diagonal
169
170 public:
172 : A(A_.getLocalMatrixDevice())
173 , eps(threshold)
174 , results(results_) {
175 // Construct ghosted matrix diagonal
177 auto lclDiag2d = diagVec->getDeviceLocalView(Xpetra::Access::ReadOnly);
178 diag = Kokkos::subview(lclDiag2d, Kokkos::ALL(), 0);
179 }
180
181 KOKKOS_FORCEINLINE_FUNCTION
182 void operator()(const local_ordinal_type rlid) const {
183 auto row = A.rowConst(rlid);
184 size_t offset = A.graph.row_map(rlid);
185 for (local_ordinal_type k = 0; k < row.length; ++k) {
186 auto clid = row.colidx(k);
187
188 auto val = row.value(k);
189 auto aiiajj = ATS::magnitude(diag(rlid)) * ATS::magnitude(diag(clid)); // |a_ii|*|a_jj|
190 const bool is_nonpositive = ATS::real(val) <= mATS::zero();
191 magnitudeType aij2 = ATS::magnitude(val) * ATS::magnitude(val); // |a_ij|^2
192 // + |a_ij|^2, if a_ij < 0, - |a_ij|^2 if a_ij >=0
193 if (is_nonpositive)
194 aij2 = -aij2;
195 results(offset + k) = Kokkos::max((aij2 <= eps * eps * aiiajj) ? DROP : KEEP,
196 results(offset + k));
197 }
198 }
199};
200
201} // namespace MueLu::ClassicalDropping
202
203#endif
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_)
typename local_matrix_type::memory_space memory_space
Kokkos::ArithTraits< scalar_type > ATS
Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > diag_vec_type
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
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
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< magnitudeType > mATS
KOKKOS_FORCEINLINE_FUNCTION void operator()(const local_ordinal_type rlid) const
typename local_matrix_type::ordinal_type local_ordinal_type
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)