Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
Ifpack2_Details_InverseDiagonalKernel_def.hpp
1// @HEADER
2// *****************************************************************************
3// Ifpack2: Templated Object-Oriented Algebraic Preconditioner Package
4//
5// Copyright 2009 NTESS and the Ifpack2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef IFPACK2_DETAILS_INVERSEDIAGONALKERNEL_DEF_HPP
11#define IFPACK2_DETAILS_INVERSEDIAGONALKERNEL_DEF_HPP
12
13#include "Tpetra_CrsMatrix.hpp"
14#include "Tpetra_MultiVector.hpp"
15#include "Tpetra_Operator.hpp"
16#include "Tpetra_Vector.hpp"
17#include "Tpetra_Export_decl.hpp"
18#include "Tpetra_Import_decl.hpp"
19#include "KokkosKernels_ArithTraits.hpp"
20#include "Teuchos_Assert.hpp"
21#include <type_traits>
22#include "KokkosSparse_spmv_impl.hpp"
23
24namespace Ifpack2 {
25namespace Details {
26namespace Impl {
27
31template <class DVector,
32 class AMatrix,
33 class DiagOffsetType,
34 bool do_L1,
35 bool fix_tiny>
37 using execution_space = typename AMatrix::execution_space;
38 using LO = typename AMatrix::non_const_ordinal_type;
39 using value_type = typename AMatrix::non_const_value_type;
40 using team_policy = typename Kokkos::TeamPolicy<execution_space>;
41 using team_member = typename team_policy::member_type;
42 using ATV = KokkosKernels::ArithTraits<value_type>;
43 // using IST = typename vector_type::impl_scalar_type;
44 using magnitude_type = typename ATV::mag_type;
45 using MATV = KokkosKernels::ArithTraits<magnitude_type>;
46
47 DVector m_d;
48 AMatrix m_A;
49 DiagOffsetType m_offsets;
50 magnitude_type m_L1Eta;
51 magnitude_type m_MinDiagonalValue;
52
54 const AMatrix& m_A_,
56 const magnitude_type m_L1Eta_,
57 const magnitude_type m_MinDiagonalValue_)
58 : m_d(m_d_)
59 , m_A(m_A_)
60 , m_offsets(m_offsets_)
61 , m_L1Eta(m_L1Eta_)
62 , m_MinDiagonalValue(m_MinDiagonalValue_) {
63 const size_t numRows = m_A.numRows();
64
65 TEUCHOS_ASSERT(numRows == size_t(m_d.extent(0)));
66 TEUCHOS_ASSERT(numRows == size_t(m_offsets.extent(0)));
67 }
68
70 void operator()(const LO lclRow) const {
71 const size_t INV = Tpetra::Details::OrdinalTraits<size_t>::invalid();
72 const value_type one = ATV::one();
73
74 // In case the row has no entries
75 m_d(lclRow, 0) = ATV::zero();
76
77 if (m_offsets(lclRow) != INV) {
78 auto curRow = m_A.rowConst(lclRow);
79 value_type d = curRow.value(m_offsets(lclRow));
80
81 if (do_L1) {
82 // Setup for L1 Methods.
83 // Here we add half the value of the off-processor entries in the row,
84 // but only if diagonal isn't sufficiently large.
85 //
86 // This follows from Equation (6.5) in: Baker, Falgout, Kolev and
87 // Yang. "Multigrid Smoothers for Ultraparallel Computing." SIAM
88 // J. Sci. Comput., Vol. 33, No. 5. (2011), pp. 2864-2887.
89
90 const magnitude_type half = MATV::one() / (MATV::one() + MATV::one());
91 const LO numRows = static_cast<LO>(m_A.numRows());
92 const LO row_length = static_cast<LO>(curRow.length);
93 magnitude_type diagonal_boost = MATV::zero();
94 for (LO iEntry = 0; iEntry < row_length; iEntry++) {
95 if (curRow.colidx(iEntry) >= numRows)
96 diagonal_boost += ATV::magnitude(curRow.value(iEntry));
97 }
99 if (ATV::magnitude(d) < m_L1Eta * diagonal_boost)
100 d += diagonal_boost;
101 }
102
103 if (fix_tiny) {
104 // Replace diagonal entries that are too small.
105
106 if (ATV::magnitude(d) <= m_MinDiagonalValue)
107 d = m_MinDiagonalValue;
108 }
109
110 // invert diagonal entries
111 m_d(lclRow, 0) = one / d;
112 }
113 }
114};
115
116} // namespace Impl
117
118template <class TpetraOperatorType>
120 InverseDiagonalKernel(const Teuchos::RCP<const operator_type>& A) {
121 setMatrix(A);
122}
123
124template <class TpetraOperatorType>
125void InverseDiagonalKernel<TpetraOperatorType>::
126 setMatrix(const Teuchos::RCP<const operator_type>& A) {
127 if (A_op_.get() != A.get()) {
128 A_op_ = A;
129
130 using Teuchos::rcp_dynamic_cast;
131 A_crs_ = rcp_dynamic_cast<const crs_matrix_type>(A);
132
133 TEUCHOS_TEST_FOR_EXCEPTION(A_crs_.is_null(), std::logic_error,
134 "Ifpack2::Details::InverseDiagonalKernel: operator A must be a Tpetra::CrsMatrix.");
135
136 const size_t lclNumRows = A_crs_->getRowMap()->getLocalNumElements();
137
138 if (offsets_.extent(0) < lclNumRows) {
139 using Kokkos::view_alloc;
140 using Kokkos::WithoutInitializing;
141 using offsets_view_type = Kokkos::View<size_t*, device_type>;
142
143 offsets_ = offsets_view_type(); // clear 1st to save mem
144 auto howAlloc = view_alloc("offsets", WithoutInitializing);
145 offsets_ = offsets_view_type(howAlloc, lclNumRows);
146 }
147
148 A_crs_->getCrsGraph()->getLocalDiagOffsets(offsets_);
149 }
150}
151
152template <class TpetraOperatorType>
153void InverseDiagonalKernel<TpetraOperatorType>::
154 compute(vector_type& D_inv,
155 bool do_l1, magnitude_type L1Eta,
156 bool fixTinyDiagEntries, magnitude_type MinDiagonalValue) {
157 // Canonicalize template arguments to avoid redundant instantiations.
158 using d_type = typename vector_type::dual_view_type::t_dev;
159 // using h_matrix_type = typename crs_matrix_type::local_matrix_host_type;
160 using d_matrix_type = typename crs_matrix_type::local_matrix_device_type;
161
162 const char kernel_label[] = "inverse_diagonal_kernel";
163 using execution_space = typename NT::execution_space;
164 using range_type = Kokkos::RangePolicy<execution_space, LO>;
165 const size_t lclNumRows = A_crs_->getRowMap()->getLocalNumElements();
166 auto policy = range_type(0, lclNumRows);
167
168 d_type d = D_inv.getLocalViewDevice(Tpetra::Access::OverwriteAll);
169 d_matrix_type a = A_crs_->getLocalMatrixDevice();
170
171 if (do_l1) {
172 constexpr bool do_l1_template = true;
173 if (fixTinyDiagEntries) {
174 constexpr bool fix_tiny_template = true;
175 using functor_type =
176 Impl::InverseDiagonalWithExtraction<d_type,
177 d_matrix_type,
178 offset_type,
179 do_l1_template,
180 fix_tiny_template>;
181 functor_type func(d, a, offsets_, L1Eta, MinDiagonalValue);
182 Kokkos::parallel_for(kernel_label, policy, func);
183 } else {
184 constexpr bool fix_tiny_template = false;
185 using functor_type =
186 Impl::InverseDiagonalWithExtraction<d_type,
187 d_matrix_type,
188 offset_type,
189 do_l1_template,
190 fix_tiny_template>;
191 functor_type func(d, a, offsets_, L1Eta, MinDiagonalValue);
192 Kokkos::parallel_for(kernel_label, policy, func);
193 }
194 } else {
195 constexpr bool do_l1_template = false;
196 if (fixTinyDiagEntries) {
197 constexpr bool fix_tiny_template = true;
198 using functor_type =
199 Impl::InverseDiagonalWithExtraction<d_type,
200 d_matrix_type,
201 offset_type,
202 do_l1_template,
203 fix_tiny_template>;
204 functor_type func(d, a, offsets_, L1Eta, MinDiagonalValue);
205 Kokkos::parallel_for(kernel_label, policy, func);
206 } else {
207 constexpr bool fix_tiny_template = false;
208 using functor_type =
209 Impl::InverseDiagonalWithExtraction<d_type,
210 d_matrix_type,
211 offset_type,
212 do_l1_template,
213 fix_tiny_template>;
214 functor_type func(d, a, offsets_, L1Eta, MinDiagonalValue);
215 Kokkos::parallel_for(kernel_label, policy, func);
216 }
217 }
218}
219
220} // namespace Details
221} // namespace Ifpack2
222
223#define IFPACK2_DETAILS_INVERSEDIAGONALKERNEL_INSTANT(SC, LO, GO, NT) \
224 template class Ifpack2::Details::InverseDiagonalKernel<Tpetra::Operator<SC, LO, GO, NT> >;
225
226#endif // IFPACK2_DETAILS_INVERSEDIAGONALKERNEL_DEF_HPP
Compute scaled damped residual for Chebyshev.
Definition Ifpack2_Details_InverseDiagonalKernel_decl.hpp:45
Ifpack2's implementation of Trilinos::Details::LinearSolver interface.
Definition Ifpack2_Details_LinearSolver_decl.hpp:75
void setMatrix(const Teuchos::RCP< const OP > &A)
Set the solver's matrix.
Definition Ifpack2_Details_LinearSolver_def.hpp:53
Ifpack2 implementation details.
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:40
Functor for extracting the inverse diagonal of a matrix.
Definition Ifpack2_Details_InverseDiagonalKernel_def.hpp:36