Tpetra parallel linear algebra Version of the Day
Loading...
Searching...
No Matches
Tpetra_leftAndOrRightScaleCrsMatrix_decl.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Tpetra: Templated Linear Algebra Services Package
4//
5// Copyright 2008 NTESS and the Tpetra contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef TPETRA_LEFTANDORRIGHTSCALECRSMATRIX_DECL_HPP
11#define TPETRA_LEFTANDORRIGHTSCALECRSMATRIX_DECL_HPP
12
15
16#include "TpetraCore_config.h"
17#if KOKKOS_VERSION >= 40799
18#include "KokkosKernels_ArithTraits.hpp"
19#else
20#include "Kokkos_ArithTraits.hpp"
21#endif
22#include "Kokkos_Core.hpp"
24#include "Tpetra_Vector_fwd.hpp"
25
26namespace Tpetra {
27
33 SCALING_MULTIPLY,
34 SCALING_DIVIDE
35};
36
63template <class SC, class LO, class GO, class NT>
65 const Kokkos::View<
66#if KOKKOS_VERSION >= 40799
67 const typename KokkosKernels::ArithTraits<SC>::mag_type*,
68#else
69 const typename Kokkos::ArithTraits<SC>::mag_type*,
70#endif
71 typename NT::device_type>& rowScalingFactors,
72 const Kokkos::View<
73#if KOKKOS_VERSION >= 40799
74 const typename KokkosKernels::ArithTraits<SC>::mag_type*,
75#else
76 const typename Kokkos::ArithTraits<SC>::mag_type*,
77#endif
78 typename NT::device_type>& colScalingFactors,
79 const bool leftScale,
80 const bool rightScale,
81 const bool assumeSymmetric,
82 const EScaling scaling);
83
111template <class SC, class LO, class GO, class NT>
113#if KOKKOS_VERSION >= 40799
114 const Tpetra::Vector<typename KokkosKernels::ArithTraits<SC>::mag_type,
115#else
116 const Tpetra::Vector<typename Kokkos::ArithTraits<SC>::mag_type,
117#endif
118 LO, GO, NT>& rowScalingFactors,
119#if KOKKOS_VERSION >= 40799
120 const Tpetra::Vector<typename KokkosKernels::ArithTraits<SC>::mag_type,
121#else
122 const Tpetra::Vector<typename Kokkos::ArithTraits<SC>::mag_type,
123#endif
124 LO, GO, NT>& colScalingFactors,
125 const bool leftScale,
126 const bool rightScale,
127 const bool assumeSymmetric,
128 const EScaling scaling);
129
130} // namespace Tpetra
131
132#endif // TPETRA_LEFTANDORRIGHTSCALECRSMATRIX_DECL_HPP
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.