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