10#ifndef TPETRA_LEFTANDORRIGHTSCALECRSMATRIX_DECL_HPP
11#define TPETRA_LEFTANDORRIGHTSCALECRSMATRIX_DECL_HPP
16#include "TpetraCore_config.h"
17#if KOKKOS_VERSION >= 40799
18#include "KokkosKernels_ArithTraits.hpp"
20#include "Kokkos_ArithTraits.hpp"
22#include "Kokkos_Core.hpp"
63template <
class SC,
class LO,
class GO,
class NT>
66#
if KOKKOS_VERSION >= 40799
67 const typename KokkosKernels::ArithTraits<SC>::mag_type*,
69 const typename Kokkos::ArithTraits<SC>::mag_type*,
71 typename NT::device_type>& rowScalingFactors,
73#
if KOKKOS_VERSION >= 40799
74 const typename KokkosKernels::ArithTraits<SC>::mag_type*,
76 const typename Kokkos::ArithTraits<SC>::mag_type*,
78 typename NT::device_type>& colScalingFactors,
80 const bool rightScale,
81 const bool assumeSymmetric,
111template <
class SC,
class LO,
class GO,
class NT>
113#
if KOKKOS_VERSION >= 40799
114 const Tpetra::Vector<
typename KokkosKernels::ArithTraits<SC>::mag_type,
118 LO, GO, NT>& rowScalingFactors,
119#
if KOKKOS_VERSION >= 40799
120 const Tpetra::Vector<
typename KokkosKernels::ArithTraits<SC>::mag_type,
124 LO, GO, NT>& colScalingFactors,
125 const bool leftScale,
126 const bool rightScale,
127 const bool assumeSymmetric,
Forward declaration of Tpetra::CrsMatrix.
Forward declaration of Tpetra::Vector.
Struct that holds views of the contents of a CrsMatrix.
A distributed dense vector.
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.