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 KokkosKernels::ArithTraits<SC>::mag_type*,
35 const typename Kokkos::ArithTraits<SC>::mag_type*,
40 const typename KokkosKernels::ArithTraits<SC>::mag_type*,
42 const typename Kokkos::ArithTraits<SC>::mag_type*,
46 const bool rightScale,
47 const bool assumeSymmetric,
49 if (!leftScale && !rightScale) {
60 auto A_lcl =
A.getLocalMatrixDevice();
62 static_cast<LO
>(
A.getRowMap()->getLocalNumElements());
64 "leftAndOrRightScaleCrsMatrix: Local matrix is not valid. "
65 "This means that A was not created with a local matrix, "
66 "and that fillComplete has never yet been called on A before. "
67 "Please call fillComplete on A at least once first "
68 "before calling this method.");
88 Teuchos::RCP<Teuchos::ParameterList>
params = Teuchos::parameterList();
89 params->set(
"No Nonlocal Changes",
true);
90 A.fillComplete(
A.getDomainMap(),
A.getRangeMap(),
params);
94template <
class SC,
class LO,
class GO,
class NT>
98 typename KokkosKernels::ArithTraits<SC>::mag_type,
100 typename Kokkos::ArithTraits<SC>::mag_type,
105 typename KokkosKernels::ArithTraits<SC>::mag_type,
107 typename Kokkos::ArithTraits<SC>::mag_type,
110 const bool leftScale,
111 const bool rightScale,
112 const bool assumeSymmetric,
114 using device_type =
typename NT::device_type;
116#if KOKKOS_VERSION >= 40799
117 using mag_type =
typename KokkosKernels::ArithTraits<SC>::mag_type;
119 using mag_type =
typename Kokkos::ArithTraits<SC>::mag_type;
121 const char prefix[] =
"leftAndOrRightScaleCrsMatrix: ";
124 Kokkos::View<const mag_type*, device_type>
row_lcl;
125 Kokkos::View<const mag_type*, device_type>
col_lcl;
130 "must be the same as the CrsMatrix's row Map. If you see this "
131 "message, it's likely that you are using a range Map Vector and that "
132 "the CrsMatrix's row Map is overlapping.");
144 "must be the same as the CrsMatrix's column Map. If you see this "
145 "message, it's likely that you are using a domain Map Vector.");
166#if KOKKOS_VERSION >= 40799
167#define TPETRA_LEFTANDORRIGHTSCALECRSMATRIX_INSTANT(SC, LO, GO, NT) \
169 leftAndOrRightScaleCrsMatrix( \
170 Tpetra::CrsMatrix<SC, LO, GO, NT>& A, \
171 const Kokkos::View< \
172 const KokkosKernels::ArithTraits<SC>::mag_type*, \
173 NT::device_type>& rowScalingFactors, \
174 const Kokkos::View< \
175 const KokkosKernels::ArithTraits<SC>::mag_type*, \
176 NT::device_type>& colScalingFactors, \
177 const bool leftScale, \
178 const bool rightScale, \
179 const bool assumeSymmetric, \
180 const EScaling scaling); \
183 leftAndOrRightScaleCrsMatrix( \
184 Tpetra::CrsMatrix<SC, LO, GO, NT>& A, \
185 const Tpetra::Vector<KokkosKernels::ArithTraits<SC>::mag_type, LO, GO, NT>& rowScalingFactors, \
186 const Tpetra::Vector<KokkosKernels::ArithTraits<SC>::mag_type, LO, GO, NT>& colScalingFactors, \
187 const bool leftScale, \
188 const bool rightScale, \
189 const bool assumeSymmetric, \
190 const EScaling scaling);
193#define TPETRA_LEFTANDORRIGHTSCALECRSMATRIX_INSTANT(SC, LO, GO, NT) \
195 leftAndOrRightScaleCrsMatrix( \
196 Tpetra::CrsMatrix<SC, LO, GO, NT>& A, \
197 const Kokkos::View< \
198 const Kokkos::ArithTraits<SC>::mag_type*, \
199 NT::device_type>& rowScalingFactors, \
200 const Kokkos::View< \
201 const Kokkos::ArithTraits<SC>::mag_type*, \
202 NT::device_type>& colScalingFactors, \
203 const bool leftScale, \
204 const bool rightScale, \
205 const bool assumeSymmetric, \
206 const EScaling scaling); \
209 leftAndOrRightScaleCrsMatrix( \
210 Tpetra::CrsMatrix<SC, LO, GO, NT>& A, \
211 const Tpetra::Vector<Kokkos::ArithTraits<SC>::mag_type, LO, GO, NT>& rowScalingFactors, \
212 const Tpetra::Vector<Kokkos::ArithTraits<SC>::mag_type, LO, GO, NT>& colScalingFactors, \
213 const bool leftScale, \
214 const bool rightScale, \
215 const bool assumeSymmetric, \
216 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.