Tpetra parallel linear algebra
Version of the Day
Toggle main menu visibility
Loading...
Searching...
No Matches
core
src
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>
71
class
LeftScaleLocalCrsMatrix
{
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
88
LeftScaleLocalCrsMatrix
(
const
LocalSparseMatrixType& A_lcl,
89
const
ScalingFactorsViewType& scalingFactors,
90
const
bool
assumeSymmetric) :
91
A_lcl_ (A_lcl),
92
scalingFactors_ (scalingFactors),
93
assumeSymmetric_ (assumeSymmetric)
94
{}
95
96
KOKKOS_INLINE_FUNCTION
void
97
operator () (
const
typename
LocalSparseMatrixType::ordinal_type lclRow)
const
98
{
99
using
LO =
typename
LocalSparseMatrixType::ordinal_type;
100
using
KAM = Kokkos::ArithTraits<mag_type>;
101
102
const
mag_type
curRowNorm = scalingFactors_(lclRow);
103
// Users are responsible for any divisions or multiplications by
104
// zero.
105
const
mag_type
scalingFactor = assumeSymmetric_ ?
106
KAM::sqrt (curRowNorm) : curRowNorm;
107
auto
curRow = A_lcl_.row (lclRow);
108
const
LO numEnt = curRow.length;
109
for
(LO k = 0; k < numEnt; ++k) {
110
if
(divide) {
// constexpr, so should get compiled out
111
curRow.value (k) = curRow.value(k) / scalingFactor;
112
}
113
else
{
114
curRow.value (k) = curRow.value(k) * scalingFactor;
115
}
116
}
117
}
118
119
private
:
120
LocalSparseMatrixType A_lcl_;
121
typename
ScalingFactorsViewType::const_type scalingFactors_;
122
bool
assumeSymmetric_;
123
};
124
138
template
<
class
LocalSparseMatrixType,
class
ScalingFactorsViewType>
139
void
140
leftScaleLocalCrsMatrix
(
const
LocalSparseMatrixType& A_lcl,
141
const
ScalingFactorsViewType& scalingFactors,
142
const
bool
assumeSymmetric,
143
const
bool
divide =
true
)
144
{
145
using
device_type
=
typename
LocalSparseMatrixType::device_type;
146
using
execution_space
=
typename
device_type::execution_space;
147
using
LO =
typename
LocalSparseMatrixType::ordinal_type;
148
using
range_type = Kokkos::RangePolicy<execution_space, LO>;
149
150
const
LO lclNumRows = A_lcl.numRows ();
151
if
(divide) {
152
using
functor_type =
153
LeftScaleLocalCrsMatrix
<LocalSparseMatrixType,
154
typename
ScalingFactorsViewType::const_type,
true
>;
155
functor_type functor (A_lcl, scalingFactors, assumeSymmetric);
156
Kokkos::parallel_for (
"leftScaleLocalCrsMatrix"
,
157
range_type (0, lclNumRows), functor);
158
}
159
else
{
160
using
functor_type =
161
LeftScaleLocalCrsMatrix
<LocalSparseMatrixType,
162
typename
ScalingFactorsViewType::const_type,
false
>;
163
functor_type functor (A_lcl, scalingFactors, assumeSymmetric);
164
Kokkos::parallel_for (
"leftScaleLocalCrsMatrix"
,
165
range_type (0, lclNumRows), functor);
166
}
167
}
168
169
}
// namespace Details
170
}
// namespace Tpetra
171
172
#endif
// TPETRA_DETAILS_LEFTSCALELOCALCRSMATRIX_HPP
Tpetra::device_type
Tpetra::mag_type
Tpetra::Details::LeftScaleLocalCrsMatrix
Kokkos::parallel_for functor that left-scales a KokkosSparse::CrsMatrix.
Definition
Tpetra_Details_leftScaleLocalCrsMatrix.hpp:71
Tpetra::Details::LeftScaleLocalCrsMatrix::LeftScaleLocalCrsMatrix
LeftScaleLocalCrsMatrix(const LocalSparseMatrixType &A_lcl, const ScalingFactorsViewType &scalingFactors, const bool assumeSymmetric)
Definition
Tpetra_Details_leftScaleLocalCrsMatrix.hpp:88
Tpetra::Details
Nonmember function that computes a residual Computes R = B - A * X.
Definition
Tpetra_KokkosRefactor_Details_MultiVectorLocalDeepCopy.hpp:51
Tpetra::Details::leftScaleLocalCrsMatrix
void leftScaleLocalCrsMatrix(const LocalSparseMatrixType &A_lcl, const ScalingFactorsViewType &scalingFactors, const bool assumeSymmetric, const bool divide=true)
Left-scale a KokkosSparse::CrsMatrix.
Definition
Tpetra_Details_leftScaleLocalCrsMatrix.hpp:140
Tpetra
Namespace Tpetra contains the class and methods constituting the Tpetra library.
Generated by
1.17.0