Tpetra parallel linear algebra  Version of the Day
Tpetra_leftAndOrRightScaleCrsMatrix_def.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Tpetra: Templated Linear Algebra Services Package
5 // Copyright (2008) Sandia Corporation
6 //
7 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
8 // the U.S. Government retains certain rights in this software.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ************************************************************************
40 // @HEADER
41 
42 #ifndef TPETRA_LEFTANDORRIGHTSCALECRSMATRIX_DEF_HPP
43 #define TPETRA_LEFTANDORRIGHTSCALECRSMATRIX_DEF_HPP
44 
51 
52 #include "Tpetra_CrsMatrix.hpp"
53 #include "Tpetra_Vector.hpp"
57 #include "Teuchos_TestForException.hpp"
58 
59 namespace Tpetra {
60 
61 template<class SC, class LO, class GO, class NT>
62 void
64  const Kokkos::View<
65  const typename Kokkos::ArithTraits<SC>::mag_type*,
66  typename NT::device_type>& rowScalingFactors,
67  const Kokkos::View<
68  const typename Kokkos::ArithTraits<SC>::mag_type*,
69  typename NT::device_type>& colScalingFactors,
70  const bool leftScale,
71  const bool rightScale,
72  const bool assumeSymmetric,
73  const EScaling scaling)
74 {
75  if (! leftScale && ! rightScale) {
76  return;
77  }
78 
79  const bool A_fillComplete_on_input = A.isFillComplete ();
80  if (! A_fillComplete_on_input) {
81  // Make sure that A has a valid local matrix. It might not if it
82  // was not created with a local matrix, and if fillComplete has
83  // never been called on it before. A never-initialized (and thus
84  // invalid) local matrix has zero rows, because it was default
85  // constructed.
86  auto A_lcl = A.getLocalMatrixDevice ();
87  const LO lclNumRows =
88  static_cast<LO> (A.getRowMap ()->getNodeNumElements ());
89  TEUCHOS_TEST_FOR_EXCEPTION
90  (A_lcl.numRows () != lclNumRows, std::invalid_argument,
91  "leftAndOrRightScaleCrsMatrix: Local matrix is not valid. "
92  "This means that A was not created with a local matrix, "
93  "and that fillComplete has never yet been called on A before. "
94  "Please call fillComplete on A at least once first "
95  "before calling this method.");
96  }
97  else {
98  A.resumeFill ();
99  }
100 
101  const bool divide = scaling == SCALING_DIVIDE;
102  if (leftScale) {
104  rowScalingFactors,
105  assumeSymmetric,
106  divide);
107  }
108  if (rightScale) {
110  colScalingFactors,
111  assumeSymmetric,
112  divide);
113  }
114 
115  if (A_fillComplete_on_input) { // put A back how we found it
116  Teuchos::RCP<Teuchos::ParameterList> params = Teuchos::parameterList ();
117  params->set ("No Nonlocal Changes", true);
118  A.fillComplete (A.getDomainMap (), A.getRangeMap (), params);
119  }
120 }
121 
122 template<class SC, class LO, class GO, class NT>
123 void
125  const Tpetra::Vector<
126  typename Kokkos::ArithTraits<SC>::mag_type,
127  LO, GO, NT>& rowScalingFactors,
128  const Tpetra::Vector<
129  typename Kokkos::ArithTraits<SC>::mag_type,
130  LO, GO, NT>& colScalingFactors,
131  const bool leftScale,
132  const bool rightScale,
133  const bool assumeSymmetric,
134  const EScaling scaling)
135 {
136  using device_type = typename NT::device_type;
137  using dev_memory_space = typename device_type::memory_space;
138  using mag_type = typename Kokkos::ArithTraits<SC>::mag_type;
139  const char prefix[] = "leftAndOrRightScaleCrsMatrix: ";
140  const bool debug = ::Tpetra::Details::Behavior::debug ();
141 
142  Kokkos::View<const mag_type*, device_type> row_lcl;
143  Kokkos::View<const mag_type*, device_type> col_lcl;
144  if (leftScale) {
145  if (debug) {
146  const bool same = rowScalingFactors.getMap ()->isSameAs (* (A.getRowMap ()));
147  TEUCHOS_TEST_FOR_EXCEPTION
148  (! same, std::invalid_argument, prefix << "rowScalingFactors's Map "
149  "must be the same as the CrsMatrix's row Map. If you see this "
150  "message, it's likely that you are using a range Map Vector and that "
151  "the CrsMatrix's row Map is overlapping.");
152  }
153  // if (rowScalingFactors.template need_sync<dev_memory_space> ()) {
154  // const_cast<vec_type&> (rowScalingFactors).template sync<dev_memory_space> ();
155  // }
156  auto row_lcl_2d = rowScalingFactors.template getLocalView<dev_memory_space> (Access::ReadOnly);
157  row_lcl = Kokkos::subview (row_lcl_2d, Kokkos::ALL (), 0);
158  }
159  if (rightScale) {
160  if (debug) {
161  const bool same = colScalingFactors.getMap ()->isSameAs (* (A.getColMap ()));
162  TEUCHOS_TEST_FOR_EXCEPTION
163  (! same, std::invalid_argument, prefix << "colScalingFactors's Map "
164  "must be the same as the CrsMatrix's column Map. If you see this "
165  "message, it's likely that you are using a domain Map Vector.");
166  }
167  // if (colScalingFactors.template need_sync<dev_memory_space> ()) {
168  // const_cast<vec_type&> (colScalingFactors).template sync<dev_memory_space> ();
169  // }
170  auto col_lcl_2d = colScalingFactors.template getLocalView<dev_memory_space> (Access::ReadOnly);
171  col_lcl = Kokkos::subview (col_lcl_2d, Kokkos::ALL (), 0);
172  }
173 
174  leftAndOrRightScaleCrsMatrix (A, row_lcl, col_lcl, leftScale, rightScale,
175  assumeSymmetric, scaling);
176 }
177 
178 } // namespace Tpetra
179 
180 //
181 // Explicit instantiation macro
182 //
183 // Must be expanded from within the Tpetra namespace!
184 //
185 
186 #define TPETRA_LEFTANDORRIGHTSCALECRSMATRIX_INSTANT(SC,LO,GO,NT) \
187  template void \
188  leftAndOrRightScaleCrsMatrix ( \
189  Tpetra::CrsMatrix<SC, LO, GO, NT>& A, \
190  const Kokkos::View< \
191  const Kokkos::ArithTraits<SC>::mag_type*, \
192  NT::device_type>& rowScalingFactors, \
193  const Kokkos::View< \
194  const Kokkos::ArithTraits<SC>::mag_type*, \
195  NT::device_type>& colScalingFactors, \
196  const bool leftScale, \
197  const bool rightScale, \
198  const bool assumeSymmetric, \
199  const EScaling scaling); \
200  \
201  template void \
202  leftAndOrRightScaleCrsMatrix ( \
203  Tpetra::CrsMatrix<SC, LO, GO, NT>& A, \
204  const Tpetra::Vector<Kokkos::ArithTraits<SC>::mag_type, LO, GO, NT>& rowScalingFactors, \
205  const Tpetra::Vector<Kokkos::ArithTraits<SC>::mag_type, LO, GO, NT>& colScalingFactors, \
206  const bool leftScale, \
207  const bool rightScale, \
208  const bool assumeSymmetric, \
209  const EScaling scaling);
210 
211 #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.
Teuchos::RCP< const map_type > getDomainMap() const override
The domain Map of this matrix.
Teuchos::RCP< const map_type > getRangeMap() const override
The range Map of this matrix.
local_matrix_device_type getLocalMatrixDevice() const
The local sparse matrix.
Teuchos::RCP< const map_type > getColMap() const override
The Map that describes the column distribution in this matrix.
void fillComplete(const Teuchos::RCP< const map_type > &domainMap, const Teuchos::RCP< const map_type > &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Tell the matrix that you are done changing its structure or values, and that you are ready to do comp...
bool isFillComplete() const override
Whether the matrix is fill complete.
Teuchos::RCP< const map_type > getRowMap() const override
The Map that describes the row distribution in this matrix.
void resumeFill(const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Resume operations that may change the values or structure of the matrix.
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.