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