Tpetra parallel linear algebra Version of the Day
Loading...
Searching...
No Matches
Tpetra_leftAndOrRightScaleCrsMatrix_def.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_DEF_HPP
11#define TPETRA_LEFTANDORRIGHTSCALECRSMATRIX_DEF_HPP
12
19
20#include "Tpetra_CrsMatrix.hpp"
21#include "Tpetra_Vector.hpp"
25#include "Teuchos_TestForException.hpp"
26
27namespace Tpetra {
28
29template <class SC, class LO, class GO, class NT>
31 const Kokkos::View<
32#if KOKKOS_VERSION >= 40799
33 const typename KokkosKernels::ArithTraits<SC>::mag_type*,
34#else
35 const typename Kokkos::ArithTraits<SC>::mag_type*,
36#endif
37 typename NT::device_type>& rowScalingFactors,
38 const Kokkos::View<
39#if KOKKOS_VERSION >= 40799
40 const typename KokkosKernels::ArithTraits<SC>::mag_type*,
41#else
42 const typename Kokkos::ArithTraits<SC>::mag_type*,
43#endif
44 typename NT::device_type>& colScalingFactors,
45 const bool leftScale,
46 const bool rightScale,
47 const bool assumeSymmetric,
48 const EScaling scaling) {
49 if (!leftScale && !rightScale) {
50 return;
51 }
52
53 const bool A_fillComplete_on_input = A.isFillComplete();
55 // Make sure that A has a valid local matrix. It might not if it
56 // was not created with a local matrix, and if fillComplete has
57 // never been called on it before. A never-initialized (and thus
58 // invalid) local matrix has zero rows, because it was default
59 // constructed.
60 auto A_lcl = A.getLocalMatrixDevice();
61 const LO lclNumRows =
62 static_cast<LO>(A.getRowMap()->getLocalNumElements());
63 TEUCHOS_TEST_FOR_EXCEPTION(A_lcl.numRows() != lclNumRows, std::invalid_argument,
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.");
69 } else {
70 A.resumeFill();
71 }
72
73 const bool divide = scaling == SCALING_DIVIDE;
74 if (leftScale) {
75 Details::leftScaleLocalCrsMatrix(A.getLocalMatrixDevice(),
77 assumeSymmetric,
78 divide);
79 }
80 if (rightScale) {
81 Details::rightScaleLocalCrsMatrix(A.getLocalMatrixDevice(),
83 assumeSymmetric,
84 divide);
85 }
86
87 if (A_fillComplete_on_input) { // put A back how we found it
88 Teuchos::RCP<Teuchos::ParameterList> params = Teuchos::parameterList();
89 params->set("No Nonlocal Changes", true);
90 A.fillComplete(A.getDomainMap(), A.getRangeMap(), params);
91 }
92}
93
94template <class SC, class LO, class GO, class NT>
96 const Tpetra::Vector<
97#if KOKKOS_VERSION >= 40799
98 typename KokkosKernels::ArithTraits<SC>::mag_type,
99#else
100 typename Kokkos::ArithTraits<SC>::mag_type,
101#endif
102 LO, GO, NT>& rowScalingFactors,
103 const Tpetra::Vector<
104#if KOKKOS_VERSION >= 40799
105 typename KokkosKernels::ArithTraits<SC>::mag_type,
106#else
107 typename Kokkos::ArithTraits<SC>::mag_type,
108#endif
109 LO, GO, NT>& colScalingFactors,
110 const bool leftScale,
111 const bool rightScale,
112 const bool assumeSymmetric,
113 const EScaling scaling) {
114 using device_type = typename NT::device_type;
115 using dev_memory_space = typename device_type::memory_space;
116#if KOKKOS_VERSION >= 40799
117 using mag_type = typename KokkosKernels::ArithTraits<SC>::mag_type;
118#else
119 using mag_type = typename Kokkos::ArithTraits<SC>::mag_type;
120#endif
121 const char prefix[] = "leftAndOrRightScaleCrsMatrix: ";
122 const bool debug = ::Tpetra::Details::Behavior::debug();
123
124 Kokkos::View<const mag_type*, device_type> row_lcl;
125 Kokkos::View<const mag_type*, device_type> col_lcl;
126 if (leftScale) {
127 if (debug) {
128 const bool same = rowScalingFactors.getMap()->isSameAs(*(A.getRowMap()));
129 TEUCHOS_TEST_FOR_EXCEPTION(!same, std::invalid_argument, prefix << "rowScalingFactors's Map "
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.");
133 }
134 // if (rowScalingFactors.template need_sync<dev_memory_space> ()) {
135 // const_cast<vec_type&> (rowScalingFactors).template sync<dev_memory_space> ();
136 // }
137 auto row_lcl_2d = rowScalingFactors.template getLocalView<dev_memory_space>(Access::ReadOnly);
138 row_lcl = Kokkos::subview(row_lcl_2d, Kokkos::ALL(), 0);
139 }
140 if (rightScale) {
141 if (debug) {
142 const bool same = colScalingFactors.getMap()->isSameAs(*(A.getColMap()));
143 TEUCHOS_TEST_FOR_EXCEPTION(!same, std::invalid_argument, prefix << "colScalingFactors's Map "
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.");
146 }
147 // if (colScalingFactors.template need_sync<dev_memory_space> ()) {
148 // const_cast<vec_type&> (colScalingFactors).template sync<dev_memory_space> ();
149 // }
150 auto col_lcl_2d = colScalingFactors.template getLocalView<dev_memory_space>(Access::ReadOnly);
151 col_lcl = Kokkos::subview(col_lcl_2d, Kokkos::ALL(), 0);
152 }
153
154 leftAndOrRightScaleCrsMatrix(A, row_lcl, col_lcl, leftScale, rightScale,
155 assumeSymmetric, scaling);
156}
157
158} // namespace Tpetra
159
160//
161// Explicit instantiation macro
162//
163// Must be expanded from within the Tpetra namespace!
164//
165
166#if KOKKOS_VERSION >= 40799
167#define TPETRA_LEFTANDORRIGHTSCALECRSMATRIX_INSTANT(SC, LO, GO, NT) \
168 template void \
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); \
181 \
182 template void \
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);
191
192#else
193#define TPETRA_LEFTANDORRIGHTSCALECRSMATRIX_INSTANT(SC, LO, GO, NT) \
194 template void \
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); \
207 \
208 template void \
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);
217
218#endif
219
220#endif // TPETRA_LEFTANDORRIGHTSCALECRSMATRIX_DEF_HPP
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.