Tpetra parallel linear algebra  Version of the Day
Tpetra_computeRowAndColumnOneNorms_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_COMPUTEROWANDCOLUMNONENORMS_DEF_HPP
43 #define TPETRA_COMPUTEROWANDCOLUMNONENORMS_DEF_HPP
44 
51 
53 #include "Tpetra_CrsMatrix.hpp"
54 #include "Tpetra_Export.hpp"
55 #include "Tpetra_Map.hpp"
56 #include "Tpetra_MultiVector.hpp"
57 #include "Tpetra_RowMatrix.hpp"
58 #include "Kokkos_Core.hpp"
59 #include "Teuchos_CommHelpers.hpp"
60 #include <memory>
61 
62 namespace Tpetra {
63 namespace Details {
64 
65 template<class SC, class LO, class GO, class NT>
66 std::size_t
67 lclMaxNumEntriesRowMatrix (const Tpetra::RowMatrix<SC, LO, GO, NT>& A)
68 {
69  const auto& rowMap = * (A.getRowMap ());
70  const LO lclNumRows = static_cast<LO> (rowMap.getNodeNumElements ());
71 
72  std::size_t maxNumEnt {0};
73  for (LO lclRow = 0; lclRow < lclNumRows; ++lclRow) {
74  const std::size_t numEnt = A.getNumEntriesInLocalRow (lclRow);
75  maxNumEnt = numEnt > maxNumEnt ? numEnt : maxNumEnt;
76  }
77  return maxNumEnt;
78 }
79 
80 template<class SC, class LO, class GO, class NT>
81 void
82 forEachLocalRowMatrixRow (
84  const LO lclNumRows,
85  const std::size_t maxNumEnt,
86  std::function<void (
87  const LO lclRow,
88  const typename Tpetra::RowMatrix<SC, LO, GO, NT>::nonconst_local_inds_host_view_type& /*ind*/,
89  const typename Tpetra::RowMatrix<SC, LO, GO, NT>::nonconst_values_host_view_type& /*val*/,
90  std::size_t /*numEnt*/ )> doForEachRow)
91 {
92  using lids_type = typename Tpetra::RowMatrix<SC, LO, GO, NT>::nonconst_local_inds_host_view_type;
93  using vals_type = typename Tpetra::RowMatrix<SC, LO, GO, NT>::nonconst_values_host_view_type;
94  lids_type indBuf("indices",maxNumEnt);
95  vals_type valBuf("values",maxNumEnt);
96 
97  for (LO lclRow = 0; lclRow < lclNumRows; ++lclRow) {
98  std::size_t numEnt = A.getNumEntriesInLocalRow (lclRow);
99  lids_type ind = Kokkos::subview(indBuf,std::make_pair((size_t)0, numEnt));
100  vals_type val = Kokkos::subview(valBuf,std::make_pair((size_t)0, numEnt));
101  A.getLocalRowCopy (lclRow, ind, val, numEnt);
102  doForEachRow (lclRow, ind, val, numEnt);
103  }
104 }
105 
106 template<class SC, class LO, class GO, class NT>
107 void
108 forEachLocalRowMatrixRow (
110  std::function<void (
111  const LO lclRow,
112  const typename Tpetra::RowMatrix<SC, LO, GO, NT>::nonconst_local_inds_host_view_type& /*ind*/,
113  const typename Tpetra::RowMatrix<SC, LO, GO, NT>::nonconst_values_host_view_type& /*val*/,
114  std::size_t /*numEnt*/ )> doForEachRow)
115 {
116  const auto& rowMap = * (A.getRowMap ());
117  const LO lclNumRows = static_cast<LO> (rowMap.getNodeNumElements ());
118  const std::size_t maxNumEnt = lclMaxNumEntriesRowMatrix (A);
119 
120  forEachLocalRowMatrixRow<SC, LO, GO, NT> (A, lclNumRows, maxNumEnt, doForEachRow);
121 }
122 
126 template<class SC, class LO, class GO, class NT>
127 void
128 computeLocalRowScaledColumnNorms_RowMatrix (EquilibrationInfo<typename Kokkos::ArithTraits<SC>::val_type,
129  typename NT::device_type>& result,
131 {
132  using KAT = Kokkos::ArithTraits<SC>;
133  using mag_type = typename KAT::mag_type;
134 
135  auto rowNorms_h = Kokkos::create_mirror_view (result.rowNorms);
136  Kokkos::deep_copy (rowNorms_h, result.rowNorms);
137  auto rowScaledColNorms_h = Kokkos::create_mirror_view (result.rowScaledColNorms);
138 
139  forEachLocalRowMatrixRow<SC, LO, GO, NT> (A,
140  [&] (const LO lclRow,
141  const typename Tpetra::RowMatrix<SC, LO, GO, NT>::nonconst_local_inds_host_view_type& ind,
142  const typename Tpetra::RowMatrix<SC, LO, GO, NT>::nonconst_values_host_view_type& val,
143  std::size_t numEnt) {
144  const mag_type rowNorm = rowNorms_h[lclRow];
145  for (std::size_t k = 0; k < numEnt; ++k) {
146  const mag_type matrixAbsVal = KAT::abs (val[k]);
147  const LO lclCol = ind[k];
148 
149  rowScaledColNorms_h[lclCol] += matrixAbsVal / rowNorm;
150  }
151  });
152  Kokkos::deep_copy (result.rowScaledColNorms, rowScaledColNorms_h);
153 }
154 
157 template<class SC, class LO, class GO, class NT>
158 EquilibrationInfo<typename Kokkos::ArithTraits<SC>::val_type, typename NT::device_type>
160 {
161  using KAT = Kokkos::ArithTraits<SC>;
162  using val_type = typename KAT::val_type;
163  using mag_type = typename KAT::mag_type;
164  using KAM = Kokkos::ArithTraits<mag_type>;
165  using device_type = typename NT::device_type;
166  using equib_info_type = EquilibrationInfo<val_type, device_type>;
167 
168  const auto& rowMap = * (A.getRowMap ());
169  const auto& colMap = * (A.getColMap ());
170  const LO lclNumRows = static_cast<LO> (rowMap.getNodeNumElements ());
171  const LO lclNumCols = 0; // don't allocate column-related Views
172  constexpr bool assumeSymmetric = false; // doesn't matter here
173  equib_info_type result (lclNumRows, lclNumCols, assumeSymmetric);
174  auto result_h = result.createMirrorView ();
175 
176  forEachLocalRowMatrixRow<SC, LO, GO, NT> (A,
177  [&] (const LO lclRow,
178  const typename Tpetra::RowMatrix<SC, LO, GO, NT>::nonconst_local_inds_host_view_type& ind,
179  const typename Tpetra::RowMatrix<SC, LO, GO, NT>::nonconst_values_host_view_type& val,
180  std::size_t numEnt) {
181  mag_type rowNorm {0.0};
182  val_type diagVal {0.0};
183  const GO gblRow = rowMap.getGlobalElement (lclRow);
184  // OK if invalid(); then we simply won't find the diagonal entry.
185  const GO lclDiagColInd = colMap.getLocalElement (gblRow);
186 
187  for (std::size_t k = 0; k < numEnt; ++k) {
188  const val_type matrixVal = val[k];
189  if (KAT::isInf (matrixVal)) {
190  result_h.foundInf = true;
191  }
192  if (KAT::isNan (matrixVal)) {
193  result_h.foundNan = true;
194  }
195  const mag_type matrixAbsVal = KAT::abs (matrixVal);
196  rowNorm += matrixAbsVal;
197  const LO lclCol = ind[k];
198  if (lclCol == lclDiagColInd) {
199  diagVal += val[k]; // repeats count additively
200  }
201  } // for each entry in row
202 
203  // This is a local result. If the matrix has an overlapping
204  // row Map, then the global result might differ.
205  if (diagVal == KAT::zero ()) {
206  result_h.foundZeroDiag = true;
207  }
208  if (rowNorm == KAM::zero ()) {
209  result_h.foundZeroRowNorm = true;
210  }
211  // NOTE (mfh 24 May 2018) We could actually compute local
212  // rowScaledColNorms in situ at this point, if ! assumeSymmetric
213  // and row Map is the same as range Map (so that the local row
214  // norms are the same as the global row norms).
215  result_h.rowDiagonalEntries[lclRow] += diagVal;
216  result_h.rowNorms[lclRow] = rowNorm;
217  });
218 
219  result.assign (result_h);
220  return result;
221 }
222 
225 template<class SC, class LO, class GO, class NT>
226 EquilibrationInfo<typename Kokkos::ArithTraits<SC>::val_type, typename NT::device_type>
228  const bool assumeSymmetric)
229 {
230  using KAT = Kokkos::ArithTraits<SC>;
231  using val_type = typename KAT::val_type;
232  using mag_type = typename KAT::mag_type;
233  using KAM = Kokkos::ArithTraits<mag_type>;
234  using device_type = typename NT::device_type;
235 
236  const auto& rowMap = * (A.getRowMap ());
237  const auto& colMap = * (A.getColMap ());
238  const LO lclNumRows = static_cast<LO> (rowMap.getNodeNumElements ());
239  const LO lclNumCols = static_cast<LO> (colMap.getNodeNumElements ());
240 
242  (lclNumRows, lclNumCols, assumeSymmetric);
243  auto result_h = result.createMirrorView ();
244 
245  forEachLocalRowMatrixRow<SC, LO, GO, NT> (A,
246  [&] (const LO lclRow,
247  const typename Tpetra::RowMatrix<SC, LO, GO, NT>::nonconst_local_inds_host_view_type& ind,
248  const typename Tpetra::RowMatrix<SC, LO, GO, NT>::nonconst_values_host_view_type& val,
249  std::size_t numEnt) {
250  mag_type rowNorm {0.0};
251  val_type diagVal {0.0};
252  const GO gblRow = rowMap.getGlobalElement (lclRow);
253  // OK if invalid(); then we simply won't find the diagonal entry.
254  const GO lclDiagColInd = colMap.getLocalElement (gblRow);
255 
256  for (std::size_t k = 0; k < numEnt; ++k) {
257  const val_type matrixVal = val[k];
258  if (KAT::isInf (matrixVal)) {
259  result_h.foundInf = true;
260  }
261  if (KAT::isNan (matrixVal)) {
262  result_h.foundNan = true;
263  }
264  const mag_type matrixAbsVal = KAT::abs (matrixVal);
265  rowNorm += matrixAbsVal;
266  const LO lclCol = ind[k];
267  if (lclCol == lclDiagColInd) {
268  diagVal += val[k]; // repeats count additively
269  }
270  if (! assumeSymmetric) {
271  result_h.colNorms[lclCol] += matrixAbsVal;
272  }
273  } // for each entry in row
274 
275  // This is a local result. If the matrix has an overlapping
276  // row Map, then the global result might differ.
277  if (diagVal == KAT::zero ()) {
278  result_h.foundZeroDiag = true;
279  }
280  if (rowNorm == KAM::zero ()) {
281  result_h.foundZeroRowNorm = true;
282  }
283  // NOTE (mfh 24 May 2018) We could actually compute local
284  // rowScaledColNorms in situ at this point, if ! assumeSymmetric
285  // and row Map is the same as range Map (so that the local row
286  // norms are the same as the global row norms).
287  result_h.rowDiagonalEntries[lclRow] += diagVal;
288  result_h.rowNorms[lclRow] = rowNorm;
289  if (! assumeSymmetric &&
290  lclDiagColInd != Tpetra::Details::OrdinalTraits<LO>::invalid ()) {
291  result_h.colDiagonalEntries[lclDiagColInd] += diagVal;
292  }
293  });
294 
295  result.assign (result_h);
296  return result;
297 }
298 
299 template<class SC, class LO, class GO, class NT>
300 class ComputeLocalRowScaledColumnNorms {
301 public:
302  using crs_matrix_type = ::Tpetra::CrsMatrix<SC, LO, GO, NT>;
303  using val_type = typename Kokkos::ArithTraits<SC>::val_type;
304  using mag_type = typename Kokkos::ArithTraits<val_type>::mag_type;
305  using device_type = typename crs_matrix_type::device_type;
306 
307  ComputeLocalRowScaledColumnNorms (const Kokkos::View<mag_type*, device_type>& rowScaledColNorms,
308  const Kokkos::View<const mag_type*, device_type>& rowNorms,
309  const crs_matrix_type& A) :
310  rowScaledColNorms_ (rowScaledColNorms),
311  rowNorms_ (rowNorms),
312  A_lcl_ (A.getLocalMatrixDevice ())
313  {}
314 
315  KOKKOS_INLINE_FUNCTION void operator () (const LO lclRow) const {
316  using KAT = Kokkos::ArithTraits<val_type>;
317 
318  const auto curRow = A_lcl_.rowConst (lclRow);
319  const mag_type rowNorm = rowNorms_[lclRow];
320  const LO numEnt = curRow.length;
321  for (LO k = 0; k < numEnt; ++k) {
322  const mag_type matrixAbsVal = KAT::abs (curRow.value(k));
323  const LO lclCol = curRow.colidx(k);
324 
325  Kokkos::atomic_add (&rowScaledColNorms_[lclCol], matrixAbsVal / rowNorm);
326  }
327  }
328 
329  static void
330  run (const Kokkos::View<mag_type*, device_type>& rowScaledColNorms,
331  const Kokkos::View<const mag_type*, device_type>& rowNorms,
332  const crs_matrix_type& A)
333  {
334  using execution_space = typename device_type::execution_space;
335  using range_type = Kokkos::RangePolicy<execution_space, LO>;
336  using functor_type = ComputeLocalRowScaledColumnNorms<SC, LO, GO, NT>;
337 
338  functor_type functor (rowScaledColNorms, rowNorms, A);
339  const LO lclNumRows =
340  static_cast<LO> (A.getRowMap ()->getNodeNumElements ());
341  Kokkos::parallel_for ("computeLocalRowScaledColumnNorms",
342  range_type (0, lclNumRows), functor);
343  }
344 
345 private:
346  Kokkos::View<mag_type*, device_type> rowScaledColNorms_;
347  Kokkos::View<const mag_type*, device_type> rowNorms_;
348 
349  using local_matrix_device_type = typename crs_matrix_type::local_matrix_device_type;
350  local_matrix_device_type A_lcl_;
351 };
352 
353 template<class SC, class LO, class GO, class NT>
354 void
355 computeLocalRowScaledColumnNorms_CrsMatrix (EquilibrationInfo<typename Kokkos::ArithTraits<SC>::val_type,
356  typename NT::device_type>& result,
358 {
359  using impl_type = ComputeLocalRowScaledColumnNorms<SC, LO, GO, NT>;
360  impl_type::run (result.rowScaledColNorms, result.rowNorms, A);
361 }
362 
363 template<class SC, class LO, class GO, class NT>
364 void
365 computeLocalRowScaledColumnNorms (EquilibrationInfo<typename Kokkos::ArithTraits<SC>::val_type,
366  typename NT::device_type>& result,
368 {
369  using crs_matrix_type = Tpetra::CrsMatrix<SC, LO, GO, NT>;
370  using val_type = typename Kokkos::ArithTraits<SC>::val_type;
371  using mag_type = typename Kokkos::ArithTraits<val_type>::mag_type;
372  using device_type = typename NT::device_type;
373 
374  auto colMapPtr = A.getColMap ();
375  TEUCHOS_TEST_FOR_EXCEPTION
376  (colMapPtr.get () == nullptr, std::invalid_argument,
377  "computeLocalRowScaledColumnNorms: "
378  "Input matrix A must have a nonnull column Map.");
379  const LO lclNumCols = static_cast<LO> (colMapPtr->getNodeNumElements ());
380  if (static_cast<std::size_t> (result.rowScaledColNorms.extent (0)) !=
381  static_cast<std::size_t> (lclNumCols)) {
382  result.rowScaledColNorms =
383  Kokkos::View<mag_type*, device_type> ("rowScaledColNorms", lclNumCols);
384  }
385 
386  const crs_matrix_type* A_crs = dynamic_cast<const crs_matrix_type*> (&A);
387  if (A_crs == nullptr) {
389  }
390  else {
391  computeLocalRowScaledColumnNorms_CrsMatrix (result, *A_crs);
392  }
393 }
394 
395 // Kokkos::parallel_reduce functor that is part of the implementation
396 // of computeLocalRowOneNorms_CrsMatrix.
397 template<class SC, class LO, class GO, class NT>
398 class ComputeLocalRowOneNorms {
399 public:
400  using val_type = typename Kokkos::ArithTraits<SC>::val_type;
401  using equib_info_type = EquilibrationInfo<val_type, typename NT::device_type>;
402  using local_matrix_device_type =
403  typename ::Tpetra::CrsMatrix<SC, LO, GO, NT>::local_matrix_device_type;
404  using local_map_type = typename ::Tpetra::Map<LO, GO, NT>::local_map_type;
405 
406  ComputeLocalRowOneNorms (const equib_info_type& equib, // in/out
407  const local_matrix_device_type& A_lcl, // in
408  const local_map_type& rowMap, // in
409  const local_map_type& colMap) : // in
410  equib_ (equib),
411  A_lcl_ (A_lcl),
412  rowMap_ (rowMap),
413  colMap_ (colMap)
414  {}
415 
416  // (result & 1) != 0 means "found Inf."
417  // (result & 2) != 0 means "found NaN."
418  // (result & 4) != 0 means "found zero diag."
419  // (result & 8) != 0 means "found zero row norm."
420  // Pack into a single int so the reduction is cheaper,
421  // esp. on GPU.
422  using value_type = int;
423 
424  KOKKOS_INLINE_FUNCTION void init (value_type& dst) const
425  {
426  dst = 0;
427  }
428 
429  KOKKOS_INLINE_FUNCTION void
430  join (volatile value_type& dst,
431  const volatile value_type& src) const
432  {
433  dst |= src;
434  }
435 
436  KOKKOS_INLINE_FUNCTION void
437  operator () (const LO lclRow, value_type& dst) const
438  {
439  using KAT = Kokkos::ArithTraits<val_type>;
440  using mag_type = typename KAT::mag_type;
441  using KAM = Kokkos::ArithTraits<mag_type>;
442 
443  const GO gblRow = rowMap_.getGlobalElement (lclRow);
444  // OK if invalid(); then we simply won't find the diagonal entry.
445  const GO lclDiagColInd = colMap_.getLocalElement (gblRow);
446 
447  const auto curRow = A_lcl_.rowConst (lclRow);
448  const LO numEnt = curRow.length;
449 
450  mag_type rowNorm {0.0};
451  val_type diagVal {0.0};
452 
453  for (LO k = 0; k < numEnt; ++k) {
454  const val_type matrixVal = curRow.value (k);
455  if (KAT::isInf (matrixVal)) {
456  dst |= 1;
457  }
458  if (KAT::isNan (matrixVal)) {
459  dst |= 2;
460  }
461  const mag_type matrixAbsVal = KAT::abs (matrixVal);
462  rowNorm += matrixAbsVal;
463  const LO lclCol = curRow.colidx (k);
464  if (lclCol == lclDiagColInd) {
465  diagVal = curRow.value (k); // assume no repeats
466  }
467  } // for each entry in row
468 
469  // This is a local result. If the matrix has an overlapping
470  // row Map, then the global result might differ.
471  if (diagVal == KAT::zero ()) {
472  dst |= 4;
473  }
474  if (rowNorm == KAM::zero ()) {
475  dst |= 8;
476  }
477  equib_.rowDiagonalEntries[lclRow] = diagVal;
478  equib_.rowNorms[lclRow] = rowNorm;
479  }
480 
481 private:
482  equib_info_type equib_;
483  local_matrix_device_type A_lcl_;
484  local_map_type rowMap_;
485  local_map_type colMap_;
486 };
487 
488 // Kokkos::parallel_reduce functor that is part of the implementation
489 // of computeLocalRowAndColumnOneNorms_CrsMatrix.
490 template<class SC, class LO, class GO, class NT>
491 class ComputeLocalRowAndColumnOneNorms {
492 public:
493  using val_type = typename Kokkos::ArithTraits<SC>::val_type;
494  using equib_info_type = EquilibrationInfo<val_type, typename NT::device_type>;
495  using local_matrix_device_type = typename ::Tpetra::CrsMatrix<SC, LO, GO, NT>::local_matrix_device_type;
496  using local_map_type = typename ::Tpetra::Map<LO, GO, NT>::local_map_type;
497 
498 public:
499  ComputeLocalRowAndColumnOneNorms (const equib_info_type& equib, // in/out
500  const local_matrix_device_type& A_lcl, // in
501  const local_map_type& rowMap, // in
502  const local_map_type& colMap) : // in
503  equib_ (equib),
504  A_lcl_ (A_lcl),
505  rowMap_ (rowMap),
506  colMap_ (colMap)
507  {}
508 
509  // (result & 1) != 0 means "found Inf."
510  // (result & 2) != 0 means "found NaN."
511  // (result & 4) != 0 means "found zero diag."
512  // (result & 8) != 0 means "found zero row norm."
513  // Pack into a single int so the reduction is cheaper,
514  // esp. on GPU.
515  using value_type = int;
516 
517  KOKKOS_INLINE_FUNCTION void init (value_type& dst) const
518  {
519  dst = 0;
520  }
521 
522  KOKKOS_INLINE_FUNCTION void
523  join (volatile value_type& dst,
524  const volatile value_type& src) const
525  {
526  dst |= src;
527  }
528 
529  KOKKOS_INLINE_FUNCTION void
530  operator () (const LO lclRow, value_type& dst) const
531  {
532  using KAT = Kokkos::ArithTraits<val_type>;
533  using mag_type = typename KAT::mag_type;
534  using KAM = Kokkos::ArithTraits<mag_type>;
535 
536  const GO gblRow = rowMap_.getGlobalElement (lclRow);
537  // OK if invalid(); then we simply won't find the diagonal entry.
538  const GO lclDiagColInd = colMap_.getLocalElement (gblRow);
539 
540  const auto curRow = A_lcl_.rowConst (lclRow);
541  const LO numEnt = curRow.length;
542 
543  mag_type rowNorm {0.0};
544  val_type diagVal {0.0};
545 
546  for (LO k = 0; k < numEnt; ++k) {
547  const val_type matrixVal = curRow.value (k);
548  if (KAT::isInf (matrixVal)) {
549  dst |= 1;
550  }
551  if (KAT::isNan (matrixVal)) {
552  dst |= 2;
553  }
554  const mag_type matrixAbsVal = KAT::abs (matrixVal);
555  rowNorm += matrixAbsVal;
556  const LO lclCol = curRow.colidx (k);
557  if (lclCol == lclDiagColInd) {
558  diagVal = curRow.value (k); // assume no repeats
559  }
560  if (! equib_.assumeSymmetric) {
561  Kokkos::atomic_add (&(equib_.colNorms[lclCol]), matrixAbsVal);
562  }
563  } // for each entry in row
564 
565  // This is a local result. If the matrix has an overlapping
566  // row Map, then the global result might differ.
567  if (diagVal == KAT::zero ()) {
568  dst |= 4;
569  }
570  if (rowNorm == KAM::zero ()) {
571  dst |= 8;
572  }
573  // NOTE (mfh 24 May 2018) We could actually compute local
574  // rowScaledColNorms in situ at this point, if ! assumeSymmetric
575  // and row Map is the same as range Map (so that the local row
576  // norms are the same as the global row norms).
577  equib_.rowDiagonalEntries[lclRow] = diagVal;
578  equib_.rowNorms[lclRow] = rowNorm;
579  if (! equib_.assumeSymmetric &&
580  lclDiagColInd != Tpetra::Details::OrdinalTraits<LO>::invalid ()) {
581  // Don't need an atomic update here, since this lclDiagColInd is
582  // a one-to-one function of lclRow.
583  equib_.colDiagonalEntries[lclDiagColInd] += diagVal;
584  }
585  }
586 
587 private:
588  equib_info_type equib_;
589  local_matrix_device_type A_lcl_;
590  local_map_type rowMap_;
591  local_map_type colMap_;
592 };
593 
596 template<class SC, class LO, class GO, class NT>
597 EquilibrationInfo<typename Kokkos::ArithTraits<SC>::val_type, typename NT::device_type>
599 {
600  using execution_space = typename NT::device_type::execution_space;
601  using range_type = Kokkos::RangePolicy<execution_space, LO>;
602  using functor_type = ComputeLocalRowOneNorms<SC, LO, GO, NT>;
603  using val_type = typename Kokkos::ArithTraits<SC>::val_type;
604  using device_type = typename NT::device_type;
605  using equib_info_type = EquilibrationInfo<val_type, device_type>;
606 
607  const LO lclNumRows = static_cast<LO> (A.getRowMap ()->getNodeNumElements ());
608  const LO lclNumCols = 0; // don't allocate column-related Views
609  constexpr bool assumeSymmetric = false; // doesn't matter here
610  equib_info_type equib (lclNumRows, lclNumCols, assumeSymmetric);
611 
612  functor_type functor (equib, A.getLocalMatrixDevice (),
613  A.getRowMap ()->getLocalMap (),
614  A.getColMap ()->getLocalMap ());
615  int result = 0;
616  Kokkos::parallel_reduce ("computeLocalRowOneNorms",
617  range_type (0, lclNumRows), functor,
618  result);
619  equib.foundInf = (result & 1) != 0;
620  equib.foundNan = (result & 2) != 0;
621  equib.foundZeroDiag = (result & 4) != 0;
622  equib.foundZeroRowNorm = (result & 8) != 0;
623  return equib;
624 }
625 
628 template<class SC, class LO, class GO, class NT>
629 EquilibrationInfo<typename Kokkos::ArithTraits<SC>::val_type, typename NT::device_type>
631  const bool assumeSymmetric)
632 {
633  using execution_space = typename NT::device_type::execution_space;
634  using range_type = Kokkos::RangePolicy<execution_space, LO>;
635  using functor_type = ComputeLocalRowAndColumnOneNorms<SC, LO, GO, NT>;
636  using val_type = typename Kokkos::ArithTraits<SC>::val_type;
637  using device_type = typename NT::device_type;
638  using equib_info_type = EquilibrationInfo<val_type, device_type>;
639 
640  const LO lclNumRows = static_cast<LO> (A.getRowMap ()->getNodeNumElements ());
641  const LO lclNumCols = static_cast<LO> (A.getColMap ()->getNodeNumElements ());
642  equib_info_type equib (lclNumRows, lclNumCols, assumeSymmetric);
643 
644  functor_type functor (equib, A.getLocalMatrixDevice (),
645  A.getRowMap ()->getLocalMap (),
646  A.getColMap ()->getLocalMap ());
647  int result = 0;
648  Kokkos::parallel_reduce ("computeLocalRowAndColumnOneNorms",
649  range_type (0, lclNumRows), functor,
650  result);
651  equib.foundInf = (result & 1) != 0;
652  equib.foundNan = (result & 2) != 0;
653  equib.foundZeroDiag = (result & 4) != 0;
654  equib.foundZeroRowNorm = (result & 8) != 0;
655  return equib;
656 }
657 
662 template<class SC, class LO, class GO, class NT>
663 EquilibrationInfo<typename Kokkos::ArithTraits<SC>::val_type,
664  typename NT::device_type>
666 {
667  using crs_matrix_type = Tpetra::CrsMatrix<SC, LO, GO, NT>;
668  const crs_matrix_type* A_crs = dynamic_cast<const crs_matrix_type*> (&A);
669 
670  if (A_crs == nullptr) {
672  }
673  else {
674  return computeLocalRowOneNorms_CrsMatrix (*A_crs);
675  }
676 }
677 
699 template<class SC, class LO, class GO, class NT>
700 EquilibrationInfo<typename Kokkos::ArithTraits<SC>::val_type, typename NT::device_type>
702  const bool assumeSymmetric)
703 {
704  using crs_matrix_type = Tpetra::CrsMatrix<SC, LO, GO, NT>;
705  const crs_matrix_type* A_crs = dynamic_cast<const crs_matrix_type*> (&A);
706 
707  if (A_crs == nullptr) {
708  return computeLocalRowAndColumnOneNorms_RowMatrix (A, assumeSymmetric);
709  }
710  else {
711  return computeLocalRowAndColumnOneNorms_CrsMatrix (*A_crs, assumeSymmetric);
712  }
713 }
714 
715 template<class SC, class LO, class GO, class NT>
716 auto getLocalView_1d_readOnly (
718  const LO whichColumn)
719 -> decltype (Kokkos::subview (X.getLocalViewDevice(Access::ReadOnly),
720  Kokkos::ALL (), whichColumn))
721 {
722  if (X.isConstantStride ()) {
723  return Kokkos::subview (X.getLocalViewDevice(Access::ReadOnly),
724  Kokkos::ALL (), whichColumn);
725  }
726  else {
727  auto X_whichColumn = X.getVector (whichColumn);
728  return Kokkos::subview (X_whichColumn->getLocalViewDevice(Access::ReadOnly),
729  Kokkos::ALL (), 0);
730  }
731 }
732 
733 template<class SC, class LO, class GO, class NT>
734 auto getLocalView_1d_writeOnly (
736  const LO whichColumn)
737 -> decltype (Kokkos::subview (X.getLocalViewDevice(Access::ReadWrite),
738  Kokkos::ALL (), whichColumn))
739 {
740  if (X.isConstantStride ()) {
741  return Kokkos::subview (X.getLocalViewDevice(Access::ReadWrite),
742  Kokkos::ALL (), whichColumn);
743  }
744  else {
745  auto X_whichColumn = X.getVectorNonConst (whichColumn);
746  return Kokkos::subview(X_whichColumn->getLocalViewDevice(Access::ReadWrite),
747  Kokkos::ALL (), 0);
748  }
749 }
750 
751 template<class SC, class LO, class GO, class NT, class ViewValueType>
752 void
753 copy1DViewIntoMultiVectorColumn (
755  const LO whichColumn,
756  const Kokkos::View<ViewValueType*, typename NT::device_type>& view)
757 {
758  auto X_lcl = getLocalView_1d_writeOnly (X, whichColumn);
759  Tpetra::Details::copyConvert (X_lcl, view);
760 }
761 
762 template<class SC, class LO, class GO, class NT, class ViewValueType>
763 void
764 copyMultiVectorColumnInto1DView (
765  const Kokkos::View<ViewValueType*, typename NT::device_type>& view,
767  const LO whichColumn)
768 {
769  auto X_lcl = getLocalView_1d_readOnly (X, whichColumn);
770  Tpetra::Details::copyConvert (view, X_lcl);
771 }
772 
773 template<class OneDViewType, class IndexType>
774 class FindZero {
775 public:
776  static_assert (OneDViewType::Rank == 1,
777  "OneDViewType must be a rank-1 Kokkos::View.");
778  static_assert (std::is_integral<IndexType>::value,
779  "IndexType must be a built-in integer type.");
780  FindZero (const OneDViewType& x) : x_ (x) {}
781  // Kokkos historically didn't like bool reduction results on CUDA,
782  // so we use int as the reduction result type.
783  KOKKOS_INLINE_FUNCTION void
784  operator () (const IndexType i, int& result) const {
785  using val_type = typename OneDViewType::non_const_value_type;
786  result = (x_(i) == Kokkos::ArithTraits<val_type>::zero ()) ? 1 : result;
787  }
788 private:
789  OneDViewType x_;
790 };
791 
792 template<class OneDViewType>
793 bool findZero (const OneDViewType& x)
794 {
795  using view_type = typename OneDViewType::const_type;
796  using execution_space = typename view_type::execution_space;
797  using size_type = typename view_type::size_type;
798  using functor_type = FindZero<view_type, size_type>;
799 
800  Kokkos::RangePolicy<execution_space, size_type> range (0, x.extent (0));
801  range.set (Kokkos::ChunkSize (500)); // adjust as needed
802 
803  int foundZero = 0;
804  Kokkos::parallel_reduce ("findZero", range, functor_type (x), foundZero);
805  return foundZero == 1;
806 }
807 
808 template<class SC, class LO, class GO, class NT>
809 void
810 globalizeRowOneNorms (EquilibrationInfo<typename Kokkos::ArithTraits<SC>::val_type,
811  typename NT::device_type>& equib,
813 {
814  using mv_type = Tpetra::MultiVector<SC, LO, GO, NT>;
815 
816  auto G = A.getGraph ();
817  TEUCHOS_TEST_FOR_EXCEPTION
818  (G.get () == nullptr, std::invalid_argument,
819  "globalizeRowOneNorms: Input RowMatrix A must have a nonnull graph "
820  "(that is, getGraph() must return nonnull).");
821  TEUCHOS_TEST_FOR_EXCEPTION
822  (! G->isFillComplete (), std::invalid_argument,
823  "globalizeRowOneNorms: Input CrsGraph G must be fillComplete.");
824 
825  auto exp = G->getExporter ();
826  if (! exp.is_null ()) {
827  // If the matrix has an overlapping row Map, first Export the
828  // local row norms with ADD CombineMode to a range Map Vector to
829  // get the global row norms, then reverse them back with REPLACE
830  // CombineMode to the row Map Vector. Ditto for the local row
831  // diagonal entries. Use SC instead of mag_type, so we can
832  // communicate both row norms and row diagonal entries at once.
833 
834  // FIXME (mfh 16 May 2018) Clever DualView tricks could possibly
835  // avoid the local copy here.
836  mv_type rowMapMV (G->getRowMap (), 2, false);
837 
838  copy1DViewIntoMultiVectorColumn (rowMapMV, 0, equib.rowNorms);
839  copy1DViewIntoMultiVectorColumn (rowMapMV, 1, equib.rowDiagonalEntries);
840  {
841  mv_type rangeMapMV (G->getRangeMap (), 2, true);
842  rangeMapMV.doExport (rowMapMV, *exp, Tpetra::ADD); // forward mode
843  rowMapMV.doImport (rangeMapMV, *exp, Tpetra::REPLACE); // reverse mode
844  }
845  copyMultiVectorColumnInto1DView (equib.rowNorms, rowMapMV, 0);
846  copyMultiVectorColumnInto1DView (equib.rowDiagonalEntries, rowMapMV, 1);
847 
848  // It's not common for users to solve linear systems with a
849  // nontrival Export, so it's OK for this to cost an additional
850  // pass over rowDiagonalEntries.
851  equib.foundZeroDiag = findZero (equib.rowDiagonalEntries);
852  equib.foundZeroRowNorm = findZero (equib.rowNorms);
853  }
854 
855  constexpr int allReduceCount = 4;
856  int lclNaughtyMatrix[allReduceCount];
857  lclNaughtyMatrix[0] = equib.foundInf ? 1 : 0;
858  lclNaughtyMatrix[1] = equib.foundNan ? 1 : 0;
859  lclNaughtyMatrix[2] = equib.foundZeroDiag ? 1 : 0;
860  lclNaughtyMatrix[3] = equib.foundZeroRowNorm ? 1 : 0;
861 
862  using Teuchos::outArg;
863  using Teuchos::REDUCE_MAX;
864  using Teuchos::reduceAll;
865  auto comm = G->getComm ();
866  int gblNaughtyMatrix[allReduceCount];
867  reduceAll<int, int> (*comm, REDUCE_MAX, allReduceCount,
868  lclNaughtyMatrix, gblNaughtyMatrix);
869 
870  equib.foundInf = gblNaughtyMatrix[0] == 1;
871  equib.foundNan = gblNaughtyMatrix[1] == 1;
872  equib.foundZeroDiag = gblNaughtyMatrix[2] == 1;
873  equib.foundZeroRowNorm = gblNaughtyMatrix[3] == 1;
874 }
875 
876 template<class SC, class LO, class GO, class NT>
877 void
878 globalizeColumnOneNorms (EquilibrationInfo<typename Kokkos::ArithTraits<SC>::val_type,
879  typename NT::device_type>& equib,
881  const bool assumeSymmetric) // if so, use row norms
882 {
883  using val_type = typename Kokkos::ArithTraits<SC>::val_type;
884  using mag_type = typename Kokkos::ArithTraits<val_type>::mag_type;
886  using device_type = typename NT::device_type;
887 
888  auto G = A.getGraph ();
889  TEUCHOS_TEST_FOR_EXCEPTION
890  (G.get () == nullptr, std::invalid_argument,
891  "globalizeColumnOneNorms: Input RowMatrix A must have a nonnull graph "
892  "(that is, getGraph() must return nonnull).");
893  TEUCHOS_TEST_FOR_EXCEPTION
894  (! G->isFillComplete (), std::invalid_argument,
895  "globalizeColumnOneNorms: Input CrsGraph G must be fillComplete.");
896 
897  auto imp = G->getImporter ();
898  if (assumeSymmetric) {
899  const LO numCols = 2;
900  // Redistribute local row info to global column info.
901 
902  // Get the data into a MultiVector on the domain Map.
903  mv_type rowNorms_domMap (G->getDomainMap (), numCols, false);
904  const bool rowMapSameAsDomainMap = G->getRowMap ()->isSameAs (* (G->getDomainMap ()));
905  if (rowMapSameAsDomainMap) {
906  copy1DViewIntoMultiVectorColumn (rowNorms_domMap, 0, equib.rowNorms);
907  copy1DViewIntoMultiVectorColumn (rowNorms_domMap, 1, equib.rowDiagonalEntries);
908  }
909  else {
910  // This is not a common case; it would normally arise when the
911  // matrix has an overlapping row Map.
912  Tpetra::Export<LO, GO, NT> rowToDom (G->getRowMap (), G->getDomainMap ());
913  mv_type rowNorms_rowMap (G->getRowMap (), numCols, true);
914  copy1DViewIntoMultiVectorColumn (rowNorms_rowMap, 0, equib.rowNorms);
915  copy1DViewIntoMultiVectorColumn (rowNorms_rowMap, 1, equib.rowDiagonalEntries);
916  rowNorms_domMap.doExport (rowNorms_rowMap, rowToDom, Tpetra::REPLACE);
917  }
918 
919  // Use the existing Import to redistribute the row norms from the
920  // domain Map to the column Map.
921  std::unique_ptr<mv_type> rowNorms_colMap;
922  if (imp.is_null ()) {
923  // Shallow copy of rowNorms_domMap.
924  rowNorms_colMap =
925  std::unique_ptr<mv_type> (new mv_type (rowNorms_domMap, * (G->getColMap ())));
926  }
927  else {
928  rowNorms_colMap =
929  std::unique_ptr<mv_type> (new mv_type (G->getColMap (), numCols, true));
930  rowNorms_colMap->doImport (rowNorms_domMap, *imp, Tpetra::REPLACE);
931  }
932 
933  // Make sure the result has allocations of the right size.
934  const LO lclNumCols =
935  static_cast<LO> (G->getColMap ()->getNodeNumElements ());
936  if (static_cast<LO> (equib.colNorms.extent (0)) != lclNumCols) {
937  equib.colNorms =
938  Kokkos::View<mag_type*, device_type> ("colNorms", lclNumCols);
939  }
940  if (static_cast<LO> (equib.colDiagonalEntries.extent (0)) != lclNumCols) {
941  equib.colDiagonalEntries =
942  Kokkos::View<val_type*, device_type> ("colDiagonalEntries", lclNumCols);
943  }
944 
945  // Copy row norms and diagonal entries, appropriately
946  // redistributed, into column norms resp. diagonal entries.
947  copyMultiVectorColumnInto1DView (equib.colNorms, *rowNorms_colMap, 0);
948  copyMultiVectorColumnInto1DView (equib.colDiagonalEntries, *rowNorms_colMap, 1);
949  }
950  else {
951  if (! imp.is_null ()) {
952  const LO numCols = 3;
953  // If the matrix has an overlapping column Map (this is usually
954  // the case), first Export (reverse-mode Import) the local info
955  // to a domain Map Vector to get the global info, then Import
956  // them back with REPLACE CombineMode to the column Map Vector.
957  // Ditto for the row-scaled column norms.
958 
959  // FIXME (mfh 16 May 2018) Clever DualView tricks could possibly
960  // avoid the local copy here.
961  mv_type colMapMV (G->getColMap (), numCols, false);
962 
963  copy1DViewIntoMultiVectorColumn (colMapMV, 0, equib.colNorms);
964  copy1DViewIntoMultiVectorColumn (colMapMV, 1, equib.colDiagonalEntries);
965  copy1DViewIntoMultiVectorColumn (colMapMV, 2, equib.rowScaledColNorms);
966  {
967  mv_type domainMapMV (G->getDomainMap (), numCols, true);
968  domainMapMV.doExport (colMapMV, *imp, Tpetra::ADD); // reverse mode
969  colMapMV.doImport (domainMapMV, *imp, Tpetra::REPLACE); // forward mode
970  }
971  copyMultiVectorColumnInto1DView (equib.colNorms, colMapMV, 0);
972  copyMultiVectorColumnInto1DView (equib.colDiagonalEntries, colMapMV, 1);
973  copyMultiVectorColumnInto1DView (equib.rowScaledColNorms, colMapMV, 2);
974  }
975  }
976 }
977 
978 } // namespace Details
979 
980 template<class SC, class LO, class GO, class NT>
981 Details::EquilibrationInfo<typename Kokkos::ArithTraits<SC>::val_type,
982  typename NT::device_type>
984 {
985  TEUCHOS_TEST_FOR_EXCEPTION
986  (! A.isFillComplete (), std::invalid_argument,
987  "computeRowOneNorms: Input matrix A must be fillComplete.");
988  auto result = Details::computeLocalRowOneNorms (A);
989 
990  Details::globalizeRowOneNorms (result, A);
991  return result;
992 }
993 
994 template<class SC, class LO, class GO, class NT>
995 Details::EquilibrationInfo<typename Kokkos::ArithTraits<SC>::val_type,
996  typename NT::device_type>
998  const bool assumeSymmetric)
999 {
1000  TEUCHOS_TEST_FOR_EXCEPTION
1001  (! A.isFillComplete (), std::invalid_argument,
1002  "computeRowAndColumnOneNorms: Input matrix A must be fillComplete.");
1003  auto result = Details::computeLocalRowAndColumnOneNorms (A, assumeSymmetric);
1004 
1005  Details::globalizeRowOneNorms (result, A);
1006  if (! assumeSymmetric) {
1007  // Row-norm-scaled column norms are trivial if the matrix is
1008  // symmetric, since the row norms and column norms are the same in
1009  // that case.
1010  Details::computeLocalRowScaledColumnNorms (result, A);
1011  }
1012  Details::globalizeColumnOneNorms (result, A, assumeSymmetric);
1013  return result;
1014 }
1015 
1016 } // namespace Tpetra
1017 
1018 //
1019 // Explicit instantiation macro
1020 //
1021 // Must be expanded from within the Tpetra namespace!
1022 //
1023 
1024 #define TPETRA_COMPUTEROWANDCOLUMNONENORMS_INSTANT(SC,LO,GO,NT) \
1025  template Details::EquilibrationInfo<Kokkos::ArithTraits<SC>::val_type, NT::device_type> \
1026  computeRowOneNorms (const Tpetra::RowMatrix<SC, LO, GO, NT>& A); \
1027  \
1028  template Details::EquilibrationInfo<Kokkos::ArithTraits<SC>::val_type, NT::device_type> \
1029  computeRowAndColumnOneNorms (const Tpetra::RowMatrix<SC, LO, GO, NT>& A, \
1030  const bool assumeSymmetric);
1031 
1032 #endif // TPETRA_COMPUTEROWANDCOLUMNONENORMS_DEF_HPP
Declare and define Tpetra::Details::copyConvert, an implementation detail of Tpetra (in particular,...
typename Node::device_type device_type
The Kokkos device type.
KokkosSparse::CrsMatrix< impl_scalar_type, local_ordinal_type, device_type, void, typename local_graph_device_type::size_type > local_matrix_device_type
The specialization of Kokkos::CrsMatrix that represents the part of the sparse matrix on each MPI pro...
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.
Teuchos::RCP< const map_type > getRowMap() const override
The Map that describes the row distribution in this matrix.
Communication plan for data redistribution from a (possibly) multiply-owned to a uniquely-owned distr...
One or more distributed dense vectors.
virtual void getLocalRowCopy(LocalOrdinal LocalRow, nonconst_local_inds_host_view_type &Indices, nonconst_values_host_view_type &Values, size_t &NumEntries) const =0
Get a copy of the given local row's entries.
virtual Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const =0
The Map that describes the distribution of rows over processes.
virtual Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const =0
The Map that describes the distribution of columns over processes.
virtual Teuchos::RCP< const RowGraph< LocalOrdinal, GlobalOrdinal, Node > > getGraph() const =0
The RowGraph associated with this matrix.
virtual size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const =0
The current number of entries on the calling process in the specified local row.
virtual bool isFillComplete() const =0
Whether fillComplete() has been called.
Implementation details of Tpetra.
EquilibrationInfo< typename Kokkos::ArithTraits< SC >::val_type, typename NT::device_type > computeLocalRowOneNorms(const Tpetra::RowMatrix< SC, LO, GO, NT > &A)
Compute LOCAL row one-norms ("row sums" etc.) of the input sparse matrix A.
EquilibrationInfo< typename Kokkos::ArithTraits< SC >::val_type, typename NT::device_type > computeLocalRowAndColumnOneNorms_CrsMatrix(const Tpetra::CrsMatrix< SC, LO, GO, NT > &A, const bool assumeSymmetric)
Implementation of computeLocalRowAndColumnOneNorms for a Tpetra::CrsMatrix.
void copyConvert(const OutputViewType &dst, const InputViewType &src)
Copy values from the 1-D Kokkos::View src, to the 1-D Kokkos::View dst, of the same length....
void computeLocalRowScaledColumnNorms_RowMatrix(EquilibrationInfo< typename Kokkos::ArithTraits< SC >::val_type, typename NT::device_type > &result, const Tpetra::RowMatrix< SC, LO, GO, NT > &A)
For a given Tpetra::RowMatrix that is not a Tpetra::CrsMatrix, assume that result....
EquilibrationInfo< typename Kokkos::ArithTraits< SC >::val_type, typename NT::device_type > computeLocalRowAndColumnOneNorms(const Tpetra::RowMatrix< SC, LO, GO, NT > &A, const bool assumeSymmetric)
Compute LOCAL row and column one-norms ("row sums" etc.) of the input sparse matrix A....
EquilibrationInfo< typename Kokkos::ArithTraits< SC >::val_type, typename NT::device_type > computeLocalRowOneNorms_CrsMatrix(const Tpetra::CrsMatrix< SC, LO, GO, NT > &A)
Implementation of computeLocalRowOneNorms for a Tpetra::CrsMatrix.
EquilibrationInfo< typename Kokkos::ArithTraits< SC >::val_type, typename NT::device_type > computeLocalRowOneNorms_RowMatrix(const Tpetra::RowMatrix< SC, LO, GO, NT > &A)
Implementation of computeLocalRowOneNorms for a Tpetra::RowMatrix that is NOT a Tpetra::CrsMatrix.
EquilibrationInfo< typename Kokkos::ArithTraits< SC >::val_type, typename NT::device_type > computeLocalRowAndColumnOneNorms_RowMatrix(const Tpetra::RowMatrix< SC, LO, GO, NT > &A, const bool assumeSymmetric)
Implementation of computeLocalRowAndColumnOneNorms for a Tpetra::RowMatrix that is NOT a Tpetra::CrsM...
Namespace Tpetra contains the class and methods constituting the Tpetra library.
void deep_copy(MultiVector< DS, DL, DG, DN > &dst, const MultiVector< SS, SL, SG, SN > &src)
Copy the contents of the MultiVector src into dst.
Details::EquilibrationInfo< typename Kokkos::ArithTraits< SC >::val_type, typename NT::device_type > computeRowOneNorms(const Tpetra::RowMatrix< SC, LO, GO, NT > &A)
Compute global row one-norms ("row sums") of the input sparse matrix A, in a way suitable for one-sid...
Details::EquilibrationInfo< typename Kokkos::ArithTraits< SC >::val_type, typename NT::device_type > computeRowAndColumnOneNorms(const Tpetra::RowMatrix< SC, LO, GO, NT > &A, const bool assumeSymmetric)
Compute global row and column one-norms ("row sums" and "column sums") of the input sparse matrix A,...
@ REPLACE
Replace existing values with new values.
@ ADD
Sum new values.
Struct storing results of Tpetra::computeRowAndColumnOneNorms.
void assign(const EquilibrationInfo< ScalarType, SrcDeviceType > &src)
Deep-copy src into *this.