MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_MatrixFreeTentativeP_def.hpp
Go to the documentation of this file.
1 // @HEADER
2 //
3 // ***********************************************************************
4 //
5 // MueLu: A package for multigrid based preconditioning
6 // Copyright 2012 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact
39 // Jonathan Hu (jhu@sandia.gov)
40 // Andrey Prokopenko (aprokop@sandia.gov)
41 // Ray Tuminaro (rstumin@sandia.gov)
42 //
43 // ***********************************************************************
44 //
45 // @HEADER
46 #ifndef MUELU_MATRIXFREETENTATIVEP_DEF_HPP
47 #define MUELU_MATRIXFREETENTATIVEP_DEF_HPP
48 
50 
51 #include "MueLu_Aggregates.hpp"
52 
53 #include <Tpetra_KokkosCompat_ClassicNodeAPI_Wrapper.hpp>
54 
55 #include "Teuchos_ScalarTraits.hpp"
56 
57 #include "Xpetra_Map.hpp"
58 #include "Xpetra_Vector.hpp"
59 #include "Xpetra_MultiVector.hpp"
60 #include "Xpetra_Operator.hpp"
61 
62 namespace MueLu {
63 
64 // compute Y = alpha*R*X + beta*Y
65 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
67  MultiVector &Y,
68  Teuchos::ETransp mode,
69  Scalar alpha,
70  Scalar beta) const {
71  using impl_scalar_type = typename Kokkos::ArithTraits<Scalar>::val_type;
72  impl_scalar_type implAlpha = alpha;
73 
74  // Step 1: Y*=beta*Y, setup necessary structures
75  Y.scale(beta);
76 
77  // TODO: probably smarter to sqrt the whole aggSizes once, but may be slower if it's done in a separate kernel launch?
78  typename Aggregates::aggregates_sizes_type::const_type aggSizes = aggregates_->ComputeAggregateSizes();
79 
80  auto kokkos_view_X = X.getDeviceLocalView(Xpetra::Access::ReadOnly);
81  auto kokkos_view_Y = Y.getDeviceLocalView(Xpetra::Access::ReadWrite);
82  LO numCols = kokkos_view_X.extent(1);
83 
84  if (mode == Teuchos::TRANS) { // if we're in restrictor mode
85  auto vertex2AggId = aggregates_->GetVertex2AggId();
86  auto vertex2AggIdView = vertex2AggId->getDeviceLocalView(Xpetra::Access::ReadOnly);
87  LO numNodes = kokkos_view_X.extent(0);
88 
89  // Step 2: Compute Y=Y+alpha*R*X
90  // recall R*X is an average of X over aggregates
91  Kokkos::parallel_for(
92  "MueLu:MatrixFreeTentativeR_kokkos:apply", md_range_type({0, 0}, {numCols, numNodes}),
93  KOKKOS_LAMBDA(const int colIdx, const int NodeIdx) {
94  LO aggIdx = vertex2AggIdView(NodeIdx, 0);
95  if (aggIdx != -1) { // depending on maps, vertex2agg
96  Kokkos::atomic_add(&kokkos_view_Y(aggIdx, colIdx), implAlpha * kokkos_view_X(NodeIdx, colIdx) / Kokkos::sqrt(aggSizes(aggIdx)));
97  }
98  });
99  } else { // if we're in prolongator mode
100  const auto vertex2Agg = aggregates_->GetVertex2AggId();
101  auto vertex2AggView = vertex2Agg->getDeviceLocalView(Xpetra::Access::ReadOnly);
102  LO numNodes = kokkos_view_Y.extent(0);
103 
104  // Step 2: Compute Y=Y+alpha*P*X
105  // recall P*X is essentially injection of X, but sum if a node belongs to multiple aggregates
106  Kokkos::parallel_for(
107  "MueLu:MatrixFreeTentativeP:apply", md_range_type({0, 0}, {numCols, numNodes}),
108  KOKKOS_LAMBDA(const int colIdx, const int fineIdx) {
109  LO aggIdx = vertex2AggView(fineIdx, 0);
110  kokkos_view_Y(fineIdx, colIdx) += implAlpha * kokkos_view_X(aggIdx, colIdx) / Kokkos::sqrt(aggSizes(aggIdx));
111  });
112  }
113 }
114 
115 // I don't care
116 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
118  TEUCHOS_TEST_FOR_EXCEPTION(true, Exceptions::RuntimeError, "MatrixFreeTentativeP residual would make no sense as the operator is not square!");
119 }
120 
121 } // namespace MueLu
122 
123 #define MUELU_MATRIXFREETENTATIVEP_SHORT
124 #endif // MUELU_MATRIXFREETENTATIVEP_DEF_HPP
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
LocalOrdinal LO
void residual(const MultiVector &X, const MultiVector &B, MultiVector &R) const override
MueLu::DefaultScalar Scalar
virtual void scale(const Scalar &alpha)=0
Kokkos::MDRangePolicy< local_ordinal_type, execution_space, Kokkos::Rank< 2 > > md_range_type
Exception throws to report errors in the internal logical of the program.
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