Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Tpetra_Details_leftScaleLocalCrsMatrix.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Tpetra: Templated Linear Algebra Services Package
5 // Copyright (2008) Sandia Corporation
6 //
7 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
8 // the U.S. Government retains certain rights in this software.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ************************************************************************
40 // @HEADER
41 
42 #ifndef TPETRA_DETAILS_LEFTSCALELOCALCRSMATRIX_HPP
43 #define TPETRA_DETAILS_LEFTSCALELOCALCRSMATRIX_HPP
44 
51 
52 #include "TpetraCore_config.h"
53 #include "Kokkos_Core.hpp"
54 #include "Kokkos_ArithTraits.hpp"
55 #include <type_traits>
56 
57 namespace Tpetra {
58 namespace Details {
59 
68 template<class LocalSparseMatrixType,
69  class ScalingFactorsViewType,
70  const bool divide>
72 public:
73  using val_type =
74  typename std::remove_const<typename LocalSparseMatrixType::value_type>::type;
75  using mag_type = typename ScalingFactorsViewType::non_const_value_type;
76  static_assert (ScalingFactorsViewType::rank == 1,
77  "scalingFactors must be a rank-1 Kokkos::View.");
78  using device_type = typename LocalSparseMatrixType::device_type;
79  using LO = typename LocalSparseMatrixType::ordinal_type;
80  using policy_type = Kokkos::TeamPolicy<typename device_type::execution_space, LO>;
81 
90  LeftScaleLocalCrsMatrix (const LocalSparseMatrixType& A_lcl,
91  const ScalingFactorsViewType& scalingFactors,
92  const bool assumeSymmetric) :
93  A_lcl_ (A_lcl),
94  scalingFactors_ (scalingFactors),
95  assumeSymmetric_ (assumeSymmetric)
96  {}
97 
98  KOKKOS_INLINE_FUNCTION void
99  operator () (const typename policy_type::member_type & team) const
100  {
101  using KAM = Kokkos::ArithTraits<mag_type>;
102 
103  const LO lclRow = team.league_rank();
104  const mag_type curRowNorm = scalingFactors_(lclRow);
105  // Users are responsible for any divisions or multiplications by
106  // zero.
107  const mag_type scalingFactor = assumeSymmetric_ ?
108  KAM::sqrt (curRowNorm) : curRowNorm;
109  auto curRow = A_lcl_.row (lclRow);
110  const LO numEnt = curRow.length;
111  Kokkos::parallel_for(Kokkos::TeamThreadRange(team, numEnt), [&](const LO k) {
112  if (divide) { // constexpr, so should get compiled out
113  curRow.value (k) = curRow.value(k) / scalingFactor;
114  }
115  else {
116  curRow.value (k) = curRow.value(k) * scalingFactor;
117  }
118  });
119  }
120 
121 private:
122  LocalSparseMatrixType A_lcl_;
123  typename ScalingFactorsViewType::const_type scalingFactors_;
124  bool assumeSymmetric_;
125 };
126 
140 template<class LocalSparseMatrixType, class ScalingFactorsViewType>
141 void
142 leftScaleLocalCrsMatrix (const LocalSparseMatrixType& A_lcl,
143  const ScalingFactorsViewType& scalingFactors,
144  const bool assumeSymmetric,
145  const bool divide = true)
146 {
147  using device_type = typename LocalSparseMatrixType::device_type;
148  using execution_space = typename device_type::execution_space;
149  using LO = typename LocalSparseMatrixType::ordinal_type;
150  using policy_type = Kokkos::TeamPolicy<execution_space, LO>;
151 
152  const LO lclNumRows = A_lcl.numRows ();
153  if (divide) {
154  using functor_type =
155  LeftScaleLocalCrsMatrix<LocalSparseMatrixType,
156  typename ScalingFactorsViewType::const_type, true>;
157  functor_type functor (A_lcl, scalingFactors, assumeSymmetric);
158  Kokkos::parallel_for ("leftScaleLocalCrsMatrix",
159  policy_type (lclNumRows, Kokkos::AUTO), functor);
160  }
161  else {
162  using functor_type =
163  LeftScaleLocalCrsMatrix<LocalSparseMatrixType,
164  typename ScalingFactorsViewType::const_type, false>;
165  functor_type functor (A_lcl, scalingFactors, assumeSymmetric);
166  Kokkos::parallel_for ("leftScaleLocalCrsMatrix",
167  policy_type (lclNumRows, Kokkos::AUTO), functor);
168  }
169 }
170 
171 } // namespace Details
172 } // namespace Tpetra
173 
174 #endif // TPETRA_DETAILS_LEFTSCALELOCALCRSMATRIX_HPP
void leftScaleLocalCrsMatrix(const LocalSparseMatrixType &A_lcl, const ScalingFactorsViewType &scalingFactors, const bool assumeSymmetric, const bool divide=true)
Left-scale a KokkosSparse::CrsMatrix.
Kokkos::parallel_for functor that left-scales a KokkosSparse::CrsMatrix.
LeftScaleLocalCrsMatrix(const LocalSparseMatrixType &A_lcl, const ScalingFactorsViewType &scalingFactors, const bool assumeSymmetric)