10#ifndef IFPACK2_DETAILS_CHEBYSHEVKERNEL_DEF_HPP
11#define IFPACK2_DETAILS_CHEBYSHEVKERNEL_DEF_HPP
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"
22#include "KokkosSparse_spmv_impl.hpp"
32template <
class WVector,
42 static_assert(
static_cast<int>(WVector::rank) == 1,
43 "WVector must be a rank 1 View.");
44 static_assert(
static_cast<int>(DVector::rank) == 1,
45 "DVector must be a rank 1 View.");
46 static_assert(
static_cast<int>(BVector::rank) == 1,
47 "BVector must be a rank 1 View.");
48 static_assert(
static_cast<int>(XVector_colMap::rank) == 1,
49 "XVector_colMap must be a rank 1 View.");
50 static_assert(
static_cast<int>(XVector_domMap::rank) == 1,
51 "XVector_domMap must be a rank 1 View.");
53 using execution_space =
typename AMatrix::execution_space;
54 using LO =
typename AMatrix::non_const_ordinal_type;
55 using value_type =
typename AMatrix::non_const_value_type;
56 using team_policy =
typename Kokkos::TeamPolicy<execution_space>;
57 using team_member =
typename team_policy::member_type;
58 using ATV = KokkosKernels::ArithTraits<value_type>;
69 const LO rows_per_team;
89 const size_t numRows = m_A.numRows();
90 const size_t numCols = m_A.numCols();
100 void operator()(
const team_member&
dev)
const {
102 using KAT = KokkosKernels::ArithTraits<residual_value_type>;
104 Kokkos::parallel_for(Kokkos::TeamThreadRange(
dev, 0, rows_per_team),
105 [&](
const LO&
loop) {
107 static_cast<LO
>(
dev.league_rank()) * rows_per_team +
loop;
108 if (
lclRow >= m_A.numRows()) {
111 const KokkosSparse::SparseRowViewConst<AMatrix>
A_row = m_A.rowConst(
lclRow);
115 Kokkos::parallel_reduce(
123 Kokkos::single(Kokkos::PerThread(
dev),
148chebyshev_kernel_vector(
const Scalar& alpha,
157 using execution_space =
typename AMatrix::execution_space;
159 if (A.numRows() == 0) {
164 int vector_length = -1;
167 const int64_t rows_per_team = KokkosSparse::Impl::spmv_launch_parameters<execution_space>(A.numRows(), A.nnz(),
rows_per_thread, team_size, vector_length);
170 using Kokkos::Dynamic;
171 using Kokkos::Schedule;
172 using Kokkos::Static;
173 using Kokkos::TeamPolicy;
183 policyDynamic = policy_type_dynamic(worksets, team_size, vector_length);
184 policyStatic = policy_type_static(worksets, team_size, vector_length);
188 using w_vec_type =
typename WVector::non_const_type;
189 using d_vec_type =
typename DVector::const_type;
190 using b_vec_type =
typename BVector::const_type;
191 using matrix_type = AMatrix;
192 using x_colMap_vec_type =
typename XVector_colMap::const_type;
193 using x_domMap_vec_type =
typename XVector_domMap::non_const_type;
194 using scalar_type =
typename KokkosKernels::ArithTraits<Scalar>::val_type;
196 if (beta == KokkosKernels::ArithTraits<Scalar>::zero()) {
197 constexpr bool use_beta =
false;
200 ChebyshevKernelVectorFunctor<w_vec_type, d_vec_type,
201 b_vec_type, matrix_type,
202 x_colMap_vec_type, x_domMap_vec_type,
206 functor_type func(alpha, w, d, b, A, x_colMap, x_domMap, beta, rows_per_team);
207 if (A.nnz() > 10000000)
208 Kokkos::parallel_for(kernel_label, policyDynamic, func);
210 Kokkos::parallel_for(kernel_label, policyStatic, func);
213 ChebyshevKernelVectorFunctor<w_vec_type, d_vec_type,
214 b_vec_type, matrix_type,
215 x_colMap_vec_type, x_domMap_vec_type,
219 functor_type func(alpha, w, d, b, A, x_colMap, x_domMap, beta, rows_per_team);
220 if (A.nnz() > 10000000)
221 Kokkos::parallel_for(kernel_label, policyDynamic, func);
223 Kokkos::parallel_for(kernel_label, policyStatic, func);
226 constexpr bool use_beta =
true;
229 ChebyshevKernelVectorFunctor<w_vec_type, d_vec_type,
230 b_vec_type, matrix_type,
231 x_colMap_vec_type, x_domMap_vec_type,
235 functor_type func(alpha, w, d, b, A, x_colMap, x_domMap, beta, rows_per_team);
236 if (A.nnz() > 10000000)
237 Kokkos::parallel_for(kernel_label, policyDynamic, func);
239 Kokkos::parallel_for(kernel_label, policyStatic, func);
242 ChebyshevKernelVectorFunctor<w_vec_type, d_vec_type,
243 b_vec_type, matrix_type,
244 x_colMap_vec_type, x_domMap_vec_type,
248 functor_type func(alpha, w, d, b, A, x_colMap, x_domMap, beta, rows_per_team);
249 if (A.nnz() > 10000000)
250 Kokkos::parallel_for(kernel_label, policyDynamic, func);
252 Kokkos::parallel_for(kernel_label, policyStatic, func);
259template <
class TpetraOperatorType>
260ChebyshevKernel<TpetraOperatorType>::
261 ChebyshevKernel(
const Teuchos::RCP<const operator_type>& A,
262 const bool useNativeSpMV)
263 : useNativeSpMV_(useNativeSpMV) {
267template <
class TpetraOperatorType>
268void ChebyshevKernel<TpetraOperatorType>::
269 setMatrix(
const Teuchos::RCP<const operator_type>& A) {
270 if (A_op_.get() != A.get()) {
274 V1_ = std::unique_ptr<multivector_type>(
nullptr);
276 using Teuchos::rcp_dynamic_cast;
277 Teuchos::RCP<const crs_matrix_type> A_crs =
278 rcp_dynamic_cast<const crs_matrix_type>(A);
279 if (A_crs.is_null()) {
280 A_crs_ = Teuchos::null;
281 imp_ = Teuchos::null;
282 exp_ = Teuchos::null;
285 TEUCHOS_ASSERT(A_crs->isFillComplete());
287 auto G = A_crs->getCrsGraph();
288 imp_ = G->getImporter();
289 exp_ = G->getExporter();
290 if (!imp_.is_null()) {
291 if (X_colMap_.get() ==
nullptr ||
292 !X_colMap_->getMap()->isSameAs(*(imp_->getTargetMap()))) {
293 X_colMap_ = std::unique_ptr<multivector_type>(
new multivector_type(imp_->getTargetMap(), 1));
301template <
class TpetraOperatorType>
302void ChebyshevKernel<TpetraOperatorType>::
303 setAuxiliaryVectors(
size_t numVectors) {
304 if ((V1_.get() ==
nullptr) || V1_->getNumVectors() != numVectors) {
305 using MV = multivector_type;
306 V1_ = std::unique_ptr<MV>(
new MV(A_op_->getRangeMap(), numVectors));
310template <
class TpetraOperatorType>
311void ChebyshevKernel<TpetraOperatorType>::
312 compute(multivector_type& W,
322 TEUCHOS_ASSERT(!A_crs_.is_null());
323 fusedCase(W, alpha, D_inv, B, *A_crs_, X, beta);
325 TEUCHOS_ASSERT(!A_op_.is_null());
326 unfusedCase(W, alpha, D_inv, B, *A_op_, X, beta);
330template <
class TpetraOperatorType>
331typename ChebyshevKernel<TpetraOperatorType>::multivector_type&
332ChebyshevKernel<TpetraOperatorType>::
333 importVector(multivector_type& X_domMap) {
334 if (imp_.is_null()) {
337 X_colMap_->doImport(X_domMap, *imp_, Tpetra::REPLACE);
342template <
class TpetraOperatorType>
343bool ChebyshevKernel<TpetraOperatorType>::
344 canFuse(
const multivector_type& B)
const {
350 return B.getNumVectors() == size_t(1) &&
355template <
class TpetraOperatorType>
356void ChebyshevKernel<TpetraOperatorType>::
357 unfusedCase(multivector_type& W,
361 const operator_type& A,
364 using STS = Teuchos::ScalarTraits<SC>;
365 setAuxiliaryVectors(B.getNumVectors());
367 const SC one = Teuchos::ScalarTraits<SC>::one();
370 Tpetra::deep_copy(*V1_, B);
371 A.apply(X, *V1_, Teuchos::NO_TRANS, -one, one);
374 W.elementWiseMultiply(alpha, D_inv, *V1_, beta);
377 X.update(STS::one(), W, STS::one());
380template <
class TpetraOperatorType>
381void ChebyshevKernel<TpetraOperatorType>::
382 fusedCase(multivector_type& W,
384 multivector_type& D_inv,
386 const crs_matrix_type& A,
389 multivector_type& X_colMap = importVector(X);
391 using Impl::chebyshev_kernel_vector;
392 using STS = Teuchos::ScalarTraits<SC>;
394 auto A_lcl = A.getLocalMatrixDevice();
396 auto Dinv_lcl = Kokkos::subview(D_inv.getLocalViewDevice(Tpetra::Access::ReadOnly), Kokkos::ALL(), 0);
397 auto B_lcl = Kokkos::subview(B.getLocalViewDevice(Tpetra::Access::ReadOnly), Kokkos::ALL(), 0);
398 auto X_domMap_lcl = Kokkos::subview(X.getLocalViewDevice(Tpetra::Access::ReadWrite), Kokkos::ALL(), 0);
399 auto X_colMap_lcl = Kokkos::subview(X_colMap.getLocalViewDevice(Tpetra::Access::ReadOnly), Kokkos::ALL(), 0);
401 const bool do_X_update = !imp_.is_null();
402 if (beta == STS::zero()) {
403 auto W_lcl = Kokkos::subview(W.getLocalViewDevice(Tpetra::Access::OverwriteAll), Kokkos::ALL(), 0);
404 chebyshev_kernel_vector(alpha, W_lcl, Dinv_lcl,
406 X_colMap_lcl, X_domMap_lcl,
410 auto W_lcl = Kokkos::subview(W.getLocalViewDevice(Tpetra::Access::ReadWrite), Kokkos::ALL(), 0);
411 chebyshev_kernel_vector(alpha, W_lcl, Dinv_lcl,
413 X_colMap_lcl, X_domMap_lcl,
418 X.update(STS::one(), W, STS::one());
424#define IFPACK2_DETAILS_CHEBYSHEVKERNEL_INSTANT(SC, LO, GO, NT) \
425 template class Ifpack2::Details::ChebyshevKernel<Tpetra::Operator<SC, LO, GO, NT> >;
Ifpack2's implementation of Trilinos::Details::LinearSolver interface.
Definition Ifpack2_Details_LinearSolver_decl.hpp:75
Ifpack2 implementation details.
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:40
Functor for computing W := alpha * D * (B - A*X) + beta * W and X := X+W.
Definition Ifpack2_Details_ChebyshevKernel_def.hpp:41