MueLu  Version of the Day
MueLu_MatrixFreeTentativeP_kokkos_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_KOKKOS_DEF_HPP
47 #define MUELU_MATRIXFREETENTATIVEP_KOKKOS_DEF_HPP
48 
50 
51 #include "MueLu_Aggregates_kokkos.hpp"
52 #include "MueLu_Level.hpp"
53 #include "MueLu_MasterList.hpp"
54 #include "MueLu_PerfUtils.hpp"
55 #include "MueLu_PFactory.hpp"
56 #include "MueLu_Monitor.hpp"
57 #include "MueLu_Utilities_kokkos.hpp"
58 
59 #include <KokkosCompat_ClassicNodeAPI_Wrapper.hpp>
60 
61 #include "Teuchos_ScalarTraits.hpp"
62 
63 namespace MueLu {
64 
65  // compute Y = alpha*R*X + beta*Y
66  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class DeviceType>
68  MultiVector &Y,
69  Teuchos::ETransp mode,
70  Scalar alpha,
71  Scalar beta) const {
72 
73  using impl_scalar_type = typename Kokkos::ArithTraits<Scalar>::val_type;
74  impl_scalar_type implAlpha = alpha;
75 
76  // Step 1: Y*=beta*Y, setup necessary structures
77  Y.scale(beta);
78 
79  // TODO: probably smarter to sqrt the whole aggSizes once, but may be slower if it's done in a separate kernel launch?
80  typename Aggregates_kokkos::aggregates_sizes_type::const_type aggSizes = aggregates_->ComputeAggregateSizes();
81 
82  auto kokkos_view_X = X.getDeviceLocalView(Xpetra::Access::ReadOnly);
83  auto kokkos_view_Y = Y.getDeviceLocalView(Xpetra::Access::ReadWrite);
84  LO numCols = kokkos_view_X.extent(1);
85 
86  if(mode == Teuchos::TRANS) { // if we're in restrictor mode
87  auto vertex2AggId = aggregates_->GetVertex2AggId();
88  auto vertex2AggIdView = vertex2AggId->getDeviceLocalView(Xpetra::Access::ReadOnly);
89  LO numNodes = kokkos_view_X.extent(0);
90 
91  // Step 2: Compute Y=Y+alpha*R*X
92  // recall R*X is an average of X over aggregates
93  Kokkos::parallel_for("MueLu:MatrixFreeTentativeR_kokkos:apply", md_range_type({0,0},{numCols,numNodes}),
94  KOKKOS_LAMBDA(const int colIdx, const int NodeIdx) {
95  LO aggIdx = vertex2AggIdView(NodeIdx,0);
96  if(aggIdx != -1) { // depending on maps, vertex2agg
97  Kokkos::atomic_add(&kokkos_view_Y(aggIdx,colIdx),implAlpha*kokkos_view_X(NodeIdx,colIdx)/Kokkos::sqrt(aggSizes(aggIdx)));
98  }
99  });
100  } else { // if we're in prolongator mode
101  const auto vertex2Agg = aggregates_->GetVertex2AggId();
102  auto vertex2AggView = vertex2Agg->getDeviceLocalView(Xpetra::Access::ReadOnly);
103  LO numNodes = kokkos_view_Y.extent(0);
104 
105  // Step 2: Compute Y=Y+alpha*P*X
106  // recall P*X is essentially injection of X, but sum if a node belongs to multiple aggregates
107  Kokkos::parallel_for("MueLu:MatrixFreeTentativeP_kokkos: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 DeviceType>
117  void MatrixFreeTentativeP_kokkos<Scalar,LocalOrdinal,GlobalOrdinal,Kokkos::Compat::KokkosDeviceWrapperNode<DeviceType>>::residual(const MultiVector &X, const MultiVector &B, MultiVector &R) const {
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_KOKKOS_SHORT
124 #endif // MUELU_MATRIXFREETENTATIVEP_KOKKOS_DEF_HPP
Namespace for MueLu classes and methods.
MueLu::DefaultScalar Scalar
Exception throws to report errors in the internal logical of the program.