Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Tpetra_Details_rightScaleLocalCrsMatrix.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_RIGHTSCALELOCALCRSMATRIX_HPP
43 #define TPETRA_DETAILS_RIGHTSCALELOCALCRSMATRIX_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 
90  RightScaleLocalCrsMatrix (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 LocalSparseMatrixType::ordinal_type lclRow) const
100  {
101  using LO = typename LocalSparseMatrixType::ordinal_type;
102  using KAM = Kokkos::ArithTraits<mag_type>;
103 
104  auto curRow = A_lcl_.row (lclRow);
105  const LO numEnt = curRow.length;
106  for (LO k = 0; k < numEnt; ++k) {
107  const LO lclColInd = curRow.colidx(k);
108  const mag_type curColNorm = scalingFactors_(lclColInd);
109  // Users are responsible for any divisions or multiplications by
110  // zero.
111  const mag_type scalingFactor = assumeSymmetric_ ?
112  KAM::sqrt (curColNorm) : curColNorm;
113  if (divide) { // constexpr, so should get compiled out
114  curRow.value(k) = curRow.value(k) / scalingFactor;
115  }
116  else {
117  curRow.value(k) = curRow.value(k) * scalingFactor;
118  }
119  }
120  }
121 
122 private:
123  LocalSparseMatrixType A_lcl_;
124  typename ScalingFactorsViewType::const_type scalingFactors_;
125  bool assumeSymmetric_;
126 };
127 
141 template<class LocalSparseMatrixType, class ScalingFactorsViewType>
142 void
143 rightScaleLocalCrsMatrix (const LocalSparseMatrixType& A_lcl,
144  const ScalingFactorsViewType& scalingFactors,
145  const bool assumeSymmetric,
146  const bool divide = true)
147 {
148  using device_type = typename LocalSparseMatrixType::device_type;
149  using execution_space = typename device_type::execution_space;
150  using LO = typename LocalSparseMatrixType::ordinal_type;
151  using range_type = Kokkos::RangePolicy<execution_space, LO>;
152 
153  const LO lclNumRows = A_lcl.numRows ();
154  if (divide) {
155  using functor_type =
156  RightScaleLocalCrsMatrix<LocalSparseMatrixType,
157  typename ScalingFactorsViewType::const_type, true>;
158  functor_type functor (A_lcl, scalingFactors, assumeSymmetric);
159  Kokkos::parallel_for ("rightScaleLocalCrsMatrix",
160  range_type (0, lclNumRows), functor);
161  }
162  else {
163  using functor_type =
164  RightScaleLocalCrsMatrix<LocalSparseMatrixType,
165  typename ScalingFactorsViewType::const_type, false>;
166  functor_type functor (A_lcl, scalingFactors, assumeSymmetric);
167  Kokkos::parallel_for ("rightScaleLocalCrsMatrix",
168  range_type (0, lclNumRows), functor);
169  }
170 }
171 
172 } // namespace Details
173 } // namespace Tpetra
174 
175 #endif // TPETRA_DETAILS_RIGHTSCALELOCALCRSMATRIX_HPP
RightScaleLocalCrsMatrix(const LocalSparseMatrixType &A_lcl, const ScalingFactorsViewType &scalingFactors, const bool assumeSymmetric)
Kokkos::parallel_for functor that right-scales a KokkosSparse::CrsMatrix.
void rightScaleLocalCrsMatrix(const LocalSparseMatrixType &A_lcl, const ScalingFactorsViewType &scalingFactors, const bool assumeSymmetric, const bool divide=true)
Right-scale a KokkosSparse::CrsMatrix.