10#ifndef TPETRA_LEFTANDORRIGHTSCALECRSMATRIX_DEF_HPP
11#define TPETRA_LEFTANDORRIGHTSCALECRSMATRIX_DEF_HPP
20#include "Tpetra_CrsMatrix.hpp"
21#include "Tpetra_Vector.hpp"
25#include "Teuchos_TestForException.hpp"
29template<
class SC,
class LO,
class GO,
class NT>
33 const typename Kokkos::ArithTraits<SC>::mag_type*,
36 const typename Kokkos::ArithTraits<SC>::mag_type*,
39 const bool rightScale,
40 const bool assumeSymmetric,
43 if (! leftScale && ! rightScale) {
54 auto A_lcl =
A.getLocalMatrixDevice ();
56 static_cast<LO
> (
A.getRowMap ()->getLocalNumElements ());
58 (A_lcl.numRows () !=
lclNumRows, std::invalid_argument,
59 "leftAndOrRightScaleCrsMatrix: Local matrix is not valid. "
60 "This means that A was not created with a local matrix, "
61 "and that fillComplete has never yet been called on A before. "
62 "Please call fillComplete on A at least once first "
63 "before calling this method.");
84 Teuchos::RCP<Teuchos::ParameterList>
params = Teuchos::parameterList ();
85 params->set (
"No Nonlocal Changes",
true);
86 A.fillComplete (
A.getDomainMap (),
A.getRangeMap (),
params);
90template<
class SC,
class LO,
class GO,
class NT>
94 typename Kokkos::ArithTraits<SC>::mag_type,
97 typename Kokkos::ArithTraits<SC>::mag_type,
100 const bool rightScale,
101 const bool assumeSymmetric,
104 using device_type =
typename NT::device_type;
106 using mag_type =
typename Kokkos::ArithTraits<SC>::mag_type;
107 const char prefix[] =
"leftAndOrRightScaleCrsMatrix: ";
110 Kokkos::View<const mag_type*, device_type>
row_lcl;
111 Kokkos::View<const mag_type*, device_type>
col_lcl;
116 (!
same, std::invalid_argument,
prefix <<
"rowScalingFactors's Map "
117 "must be the same as the CrsMatrix's row Map. If you see this "
118 "message, it's likely that you are using a range Map Vector and that "
119 "the CrsMatrix's row Map is overlapping.");
131 (!
same, std::invalid_argument,
prefix <<
"colScalingFactors's Map "
132 "must be the same as the CrsMatrix's column Map. If you see this "
133 "message, it's likely that you are using a domain Map Vector.");
154#define TPETRA_LEFTANDORRIGHTSCALECRSMATRIX_INSTANT(SC,LO,GO,NT) \
156 leftAndOrRightScaleCrsMatrix ( \
157 Tpetra::CrsMatrix<SC, LO, GO, NT>& A, \
158 const Kokkos::View< \
159 const Kokkos::ArithTraits<SC>::mag_type*, \
160 NT::device_type>& rowScalingFactors, \
161 const Kokkos::View< \
162 const Kokkos::ArithTraits<SC>::mag_type*, \
163 NT::device_type>& colScalingFactors, \
164 const bool leftScale, \
165 const bool rightScale, \
166 const bool assumeSymmetric, \
167 const EScaling scaling); \
170 leftAndOrRightScaleCrsMatrix ( \
171 Tpetra::CrsMatrix<SC, LO, GO, NT>& A, \
172 const Tpetra::Vector<Kokkos::ArithTraits<SC>::mag_type, LO, GO, NT>& rowScalingFactors, \
173 const Tpetra::Vector<Kokkos::ArithTraits<SC>::mag_type, LO, GO, NT>& colScalingFactors, \
174 const bool leftScale, \
175 const bool rightScale, \
176 const bool assumeSymmetric, \
177 const EScaling scaling);
Declaration of Tpetra::Details::Behavior, a class that describes Tpetra's behavior.
Declaration and definition of Tpetra::Details::leftScaleLocalCrsMatrix.
Declaration and definition of Tpetra::Details::rightScaleLocalCrsMatrix.
Struct that holds views of the contents of a CrsMatrix.
static bool debug()
Whether Tpetra is in debug mode.
A distributed dense vector.
void leftScaleLocalCrsMatrix(const LocalSparseMatrixType &A_lcl, const ScalingFactorsViewType &scalingFactors, const bool assumeSymmetric, const bool divide=true)
Left-scale a KokkosSparse::CrsMatrix.
void rightScaleLocalCrsMatrix(const LocalSparseMatrixType &A_lcl, const ScalingFactorsViewType &scalingFactors, const bool assumeSymmetric, const bool divide=true)
Right-scale a KokkosSparse::CrsMatrix.
Namespace Tpetra contains the class and methods constituting the Tpetra library.
EScaling
Whether "scaling" a matrix means multiplying or dividing its entries.
void leftAndOrRightScaleCrsMatrix(Tpetra::CrsMatrix< SC, LO, GO, NT > &A, const Kokkos::View< const typename Kokkos::ArithTraits< SC >::mag_type *, typename NT::device_type > &rowScalingFactors, const Kokkos::View< const typename Kokkos::ArithTraits< SC >::mag_type *, typename NT::device_type > &colScalingFactors, const bool leftScale, const bool rightScale, const bool assumeSymmetric, const EScaling scaling)
Left-scale and/or right-scale (in that order) the entries of the input Tpetra::CrsMatrix A.