MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_MatrixFreeTentativeP_def.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_MATRIXFREETENTATIVEP_DEF_HPP
11#define MUELU_MATRIXFREETENTATIVEP_DEF_HPP
12
14
15#include "MueLu_Aggregates.hpp"
16
17#include <Tpetra_KokkosCompat_ClassicNodeAPI_Wrapper.hpp>
18
19#include "Teuchos_ScalarTraits.hpp"
20
21#include "Xpetra_Map.hpp"
22#include "Xpetra_Vector.hpp"
23#include "Xpetra_MultiVector.hpp"
24#include "Xpetra_Operator.hpp"
25
26namespace MueLu {
27
28// compute Y = alpha*R*X + beta*Y
29template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
31 MultiVector &Y,
32 Teuchos::ETransp mode,
33 Scalar alpha,
34 Scalar beta) const {
35 using impl_scalar_type = typename KokkosKernels::ArithTraits<Scalar>::val_type;
36 impl_scalar_type implAlpha = alpha;
37
38 // Step 1: Y*=beta*Y, setup necessary structures
39 Y.scale(beta);
40
41 // TODO: probably smarter to sqrt the whole aggSizes once, but may be slower if it's done in a separate kernel launch?
42 typename Aggregates::aggregates_sizes_type::const_type aggSizes = aggregates_->ComputeAggregateSizes();
43
44 auto kokkos_view_X = X.getLocalViewDevice(Tpetra::Access::ReadOnly);
45 auto kokkos_view_Y = Y.getLocalViewDevice(Tpetra::Access::ReadWrite);
46 LO numCols = kokkos_view_X.extent(1);
47
48 if (mode == Teuchos::TRANS) { // if we're in restrictor mode
49 auto vertex2AggId = aggregates_->GetVertex2AggId();
50 auto vertex2AggIdView = vertex2AggId->getLocalViewDevice(Tpetra::Access::ReadOnly);
51 LO numNodes = kokkos_view_X.extent(0);
52
53 // Step 2: Compute Y=Y+alpha*R*X
54 // recall R*X is an average of X over aggregates
55 Kokkos::parallel_for(
56 "MueLu:MatrixFreeTentativeR_kokkos:apply", md_range_type({0, 0}, {numCols, numNodes}),
57 KOKKOS_LAMBDA(const int colIdx, const int NodeIdx) {
58 LO aggIdx = vertex2AggIdView(NodeIdx, 0);
59 if (aggIdx != -1) { // depending on maps, vertex2agg
60 Kokkos::atomic_add(&kokkos_view_Y(aggIdx, colIdx), implAlpha * kokkos_view_X(NodeIdx, colIdx) / Kokkos::sqrt(aggSizes(aggIdx)));
61 }
62 });
63 } else { // if we're in prolongator mode
64 const auto vertex2Agg = aggregates_->GetVertex2AggId();
65 auto vertex2AggView = vertex2Agg->getLocalViewDevice(Tpetra::Access::ReadOnly);
66 LO numNodes = kokkos_view_Y.extent(0);
67
68 // Step 2: Compute Y=Y+alpha*P*X
69 // recall P*X is essentially injection of X, but sum if a node belongs to multiple aggregates
70 Kokkos::parallel_for(
71 "MueLu:MatrixFreeTentativeP:apply", md_range_type({0, 0}, {numCols, numNodes}),
72 KOKKOS_LAMBDA(const int colIdx, const int fineIdx) {
73 LO aggIdx = vertex2AggView(fineIdx, 0);
74 kokkos_view_Y(fineIdx, colIdx) += implAlpha * kokkos_view_X(aggIdx, colIdx) / Kokkos::sqrt(aggSizes(aggIdx));
75 });
76 }
77}
78
79// I don't care
80template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
81void MatrixFreeTentativeP<Scalar, LocalOrdinal, GlobalOrdinal, Node>::residual(const MultiVector &X, const MultiVector &B, MultiVector &R) const {
82 TEUCHOS_TEST_FOR_EXCEPTION(true, Exceptions::RuntimeError, "MatrixFreeTentativeP residual would make no sense as the operator is not square!");
83}
84
85} // namespace MueLu
86
87#define MUELU_MATRIXFREETENTATIVEP_SHORT
88#endif // MUELU_MATRIXFREETENTATIVEP_DEF_HPP
MueLu::DefaultScalar Scalar
Exception throws to report errors in the internal logical of the program.
void residual(const MultiVector &X, const MultiVector &B, MultiVector &R) const override
Kokkos::MDRangePolicy< local_ordinal_type, execution_space, Kokkos::Rank< 2 > > md_range_type
void apply(const MultiVector &X, MultiVector &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=Teuchos::ScalarTraits< Scalar >::one(), Scalar beta=Teuchos::ScalarTraits< Scalar >::zero()) const override
Namespace for MueLu classes and methods.