42 #ifndef TPETRA_COMPUTEROWANDCOLUMNONENORMS_DEF_HPP
43 #define TPETRA_COMPUTEROWANDCOLUMNONENORMS_DEF_HPP
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"
65 template<
class SC,
class LO,
class GO,
class NT>
70 const LO lclNumRows =
static_cast<LO
> (rowMap.getNodeNumElements ());
72 std::size_t maxNumEnt {0};
73 for (LO lclRow = 0; lclRow < lclNumRows; ++lclRow) {
75 maxNumEnt = numEnt > maxNumEnt ? numEnt : maxNumEnt;
80 template<
class SC,
class LO,
class GO,
class NT>
82 forEachLocalRowMatrixRow (
85 const std::size_t maxNumEnt,
88 const typename Tpetra::RowMatrix<SC, LO, GO, NT>::nonconst_local_inds_host_view_type& ,
89 const typename Tpetra::RowMatrix<SC, LO, GO, NT>::nonconst_values_host_view_type& ,
90 std::size_t )> doForEachRow)
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);
97 for (LO lclRow = 0; lclRow < lclNumRows; ++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));
102 doForEachRow (lclRow, ind, val, numEnt);
106 template<
class SC,
class LO,
class GO,
class NT>
108 forEachLocalRowMatrixRow (
112 const typename Tpetra::RowMatrix<SC, LO, GO, NT>::nonconst_local_inds_host_view_type& ,
113 const typename Tpetra::RowMatrix<SC, LO, GO, NT>::nonconst_values_host_view_type& ,
114 std::size_t )> doForEachRow)
117 const LO lclNumRows =
static_cast<LO
> (rowMap.getNodeNumElements ());
118 const std::size_t maxNumEnt = lclMaxNumEntriesRowMatrix (A);
120 forEachLocalRowMatrixRow<SC, LO, GO, NT> (A, lclNumRows, maxNumEnt, doForEachRow);
126 template<
class SC,
class LO,
class GO,
class NT>
129 typename NT::device_type>& result,
132 using KAT = Kokkos::ArithTraits<SC>;
133 using mag_type =
typename KAT::mag_type;
135 auto rowNorms_h = Kokkos::create_mirror_view (result.rowNorms);
137 auto rowScaledColNorms_h = Kokkos::create_mirror_view (result.rowScaledColNorms);
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];
149 rowScaledColNorms_h[lclCol] += matrixAbsVal / rowNorm;
157 template<
class SC,
class LO,
class GO,
class NT>
158 EquilibrationInfo<typename Kokkos::ArithTraits<SC>::val_type,
typename NT::device_type>
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;
170 const LO lclNumRows =
static_cast<LO
> (rowMap.getNodeNumElements ());
171 const LO lclNumCols = 0;
172 constexpr
bool assumeSymmetric =
false;
173 equib_info_type result (lclNumRows, lclNumCols, assumeSymmetric);
174 auto result_h = result.createMirrorView ();
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);
185 const GO lclDiagColInd = colMap.getLocalElement (gblRow);
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;
192 if (KAT::isNan (matrixVal)) {
193 result_h.foundNan =
true;
195 const mag_type matrixAbsVal = KAT::abs (matrixVal);
196 rowNorm += matrixAbsVal;
197 const LO lclCol = ind[k];
198 if (lclCol == lclDiagColInd) {
205 if (diagVal == KAT::zero ()) {
206 result_h.foundZeroDiag =
true;
208 if (rowNorm == KAM::zero ()) {
209 result_h.foundZeroRowNorm =
true;
215 result_h.rowDiagonalEntries[lclRow] += diagVal;
216 result_h.rowNorms[lclRow] = rowNorm;
219 result.assign (result_h);
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)
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;
238 const LO lclNumRows =
static_cast<LO
> (rowMap.getNodeNumElements ());
239 const LO lclNumCols =
static_cast<LO
> (colMap.getNodeNumElements ());
242 (lclNumRows, lclNumCols, assumeSymmetric);
243 auto result_h = result.createMirrorView ();
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);
254 const GO lclDiagColInd = colMap.getLocalElement (gblRow);
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;
261 if (KAT::isNan (matrixVal)) {
262 result_h.foundNan =
true;
264 const mag_type matrixAbsVal = KAT::abs (matrixVal);
265 rowNorm += matrixAbsVal;
266 const LO lclCol = ind[k];
267 if (lclCol == lclDiagColInd) {
270 if (! assumeSymmetric) {
271 result_h.colNorms[lclCol] += matrixAbsVal;
277 if (diagVal == KAT::zero ()) {
278 result_h.foundZeroDiag =
true;
280 if (rowNorm == KAM::zero ()) {
281 result_h.foundZeroRowNorm =
true;
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;
299 template<
class SC,
class LO,
class GO,
class NT>
300 class ComputeLocalRowScaledColumnNorms {
303 using val_type =
typename Kokkos::ArithTraits<SC>::val_type;
304 using mag_type =
typename Kokkos::ArithTraits<val_type>::mag_type;
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 ())
315 KOKKOS_INLINE_FUNCTION
void operator () (
const LO lclRow)
const {
316 using KAT = Kokkos::ArithTraits<val_type>;
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);
325 Kokkos::atomic_add (&rowScaledColNorms_[lclCol], matrixAbsVal / rowNorm);
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)
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>;
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);
346 Kokkos::View<mag_type*, device_type> rowScaledColNorms_;
347 Kokkos::View<const mag_type*, device_type> rowNorms_;
350 local_matrix_device_type A_lcl_;
353 template<
class SC,
class LO,
class GO,
class NT>
355 computeLocalRowScaledColumnNorms_CrsMatrix (EquilibrationInfo<
typename Kokkos::ArithTraits<SC>::val_type,
356 typename NT::device_type>& result,
359 using impl_type = ComputeLocalRowScaledColumnNorms<SC, LO, GO, NT>;
360 impl_type::run (result.rowScaledColNorms, result.rowNorms, A);
363 template<
class SC,
class LO,
class GO,
class NT>
365 computeLocalRowScaledColumnNorms (EquilibrationInfo<
typename Kokkos::ArithTraits<SC>::val_type,
366 typename NT::device_type>& result,
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;
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);
386 const crs_matrix_type* A_crs =
dynamic_cast<const crs_matrix_type*
> (&A);
387 if (A_crs ==
nullptr) {
391 computeLocalRowScaledColumnNorms_CrsMatrix (result, *A_crs);
397 template<
class SC,
class LO,
class GO,
class NT>
398 class ComputeLocalRowOneNorms {
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;
406 ComputeLocalRowOneNorms (
const equib_info_type& equib,
407 const local_matrix_device_type& A_lcl,
408 const local_map_type& rowMap,
409 const local_map_type& colMap) :
422 using value_type = int;
424 KOKKOS_INLINE_FUNCTION
void init (value_type& dst)
const
429 KOKKOS_INLINE_FUNCTION
void
430 join (
volatile value_type& dst,
431 const volatile value_type& src)
const
436 KOKKOS_INLINE_FUNCTION
void
437 operator () (
const LO lclRow, value_type& dst)
const
439 using KAT = Kokkos::ArithTraits<val_type>;
440 using mag_type =
typename KAT::mag_type;
441 using KAM = Kokkos::ArithTraits<mag_type>;
443 const GO gblRow = rowMap_.getGlobalElement (lclRow);
445 const GO lclDiagColInd = colMap_.getLocalElement (gblRow);
447 const auto curRow = A_lcl_.rowConst (lclRow);
448 const LO numEnt = curRow.length;
450 mag_type rowNorm {0.0};
451 val_type diagVal {0.0};
453 for (LO k = 0; k < numEnt; ++k) {
454 const val_type matrixVal = curRow.value (k);
455 if (KAT::isInf (matrixVal)) {
458 if (KAT::isNan (matrixVal)) {
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);
471 if (diagVal == KAT::zero ()) {
474 if (rowNorm == KAM::zero ()) {
477 equib_.rowDiagonalEntries[lclRow] = diagVal;
478 equib_.rowNorms[lclRow] = rowNorm;
482 equib_info_type equib_;
483 local_matrix_device_type A_lcl_;
484 local_map_type rowMap_;
485 local_map_type colMap_;
490 template<
class SC,
class LO,
class GO,
class NT>
491 class ComputeLocalRowAndColumnOneNorms {
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;
499 ComputeLocalRowAndColumnOneNorms (
const equib_info_type& equib,
500 const local_matrix_device_type& A_lcl,
501 const local_map_type& rowMap,
502 const local_map_type& colMap) :
515 using value_type = int;
517 KOKKOS_INLINE_FUNCTION
void init (value_type& dst)
const
522 KOKKOS_INLINE_FUNCTION
void
523 join (
volatile value_type& dst,
524 const volatile value_type& src)
const
529 KOKKOS_INLINE_FUNCTION
void
530 operator () (
const LO lclRow, value_type& dst)
const
532 using KAT = Kokkos::ArithTraits<val_type>;
533 using mag_type =
typename KAT::mag_type;
534 using KAM = Kokkos::ArithTraits<mag_type>;
536 const GO gblRow = rowMap_.getGlobalElement (lclRow);
538 const GO lclDiagColInd = colMap_.getLocalElement (gblRow);
540 const auto curRow = A_lcl_.rowConst (lclRow);
541 const LO numEnt = curRow.length;
543 mag_type rowNorm {0.0};
544 val_type diagVal {0.0};
546 for (LO k = 0; k < numEnt; ++k) {
547 const val_type matrixVal = curRow.value (k);
548 if (KAT::isInf (matrixVal)) {
551 if (KAT::isNan (matrixVal)) {
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);
560 if (! equib_.assumeSymmetric) {
561 Kokkos::atomic_add (&(equib_.colNorms[lclCol]), matrixAbsVal);
567 if (diagVal == KAT::zero ()) {
570 if (rowNorm == KAM::zero ()) {
577 equib_.rowDiagonalEntries[lclRow] = diagVal;
578 equib_.rowNorms[lclRow] = rowNorm;
579 if (! equib_.assumeSymmetric &&
580 lclDiagColInd != Tpetra::Details::OrdinalTraits<LO>::invalid ()) {
583 equib_.colDiagonalEntries[lclDiagColInd] += diagVal;
588 equib_info_type equib_;
589 local_matrix_device_type A_lcl_;
590 local_map_type rowMap_;
591 local_map_type colMap_;
596 template<
class SC,
class LO,
class GO,
class NT>
597 EquilibrationInfo<typename Kokkos::ArithTraits<SC>::val_type,
typename NT::device_type>
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;
607 const LO lclNumRows =
static_cast<LO
> (A.
getRowMap ()->getNodeNumElements ());
608 const LO lclNumCols = 0;
609 constexpr
bool assumeSymmetric =
false;
610 equib_info_type equib (lclNumRows, lclNumCols, assumeSymmetric);
616 Kokkos::parallel_reduce (
"computeLocalRowOneNorms",
617 range_type (0, lclNumRows), functor,
619 equib.foundInf = (result & 1) != 0;
620 equib.foundNan = (result & 2) != 0;
621 equib.foundZeroDiag = (result & 4) != 0;
622 equib.foundZeroRowNorm = (result & 8) != 0;
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)
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;
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);
648 Kokkos::parallel_reduce (
"computeLocalRowAndColumnOneNorms",
649 range_type (0, lclNumRows), functor,
651 equib.foundInf = (result & 1) != 0;
652 equib.foundNan = (result & 2) != 0;
653 equib.foundZeroDiag = (result & 4) != 0;
654 equib.foundZeroRowNorm = (result & 8) != 0;
662 template<
class SC,
class LO,
class GO,
class NT>
663 EquilibrationInfo<typename Kokkos::ArithTraits<SC>::val_type,
664 typename NT::device_type>
668 const crs_matrix_type* A_crs =
dynamic_cast<const crs_matrix_type*
> (&A);
670 if (A_crs ==
nullptr) {
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)
705 const crs_matrix_type* A_crs =
dynamic_cast<const crs_matrix_type*
> (&A);
707 if (A_crs ==
nullptr) {
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))
722 if (X.isConstantStride ()) {
723 return Kokkos::subview (X.getLocalViewDevice(Access::ReadOnly),
724 Kokkos::ALL (), whichColumn);
727 auto X_whichColumn = X.getVector (whichColumn);
728 return Kokkos::subview (X_whichColumn->getLocalViewDevice(Access::ReadOnly),
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))
740 if (X.isConstantStride ()) {
741 return Kokkos::subview (X.getLocalViewDevice(Access::ReadWrite),
742 Kokkos::ALL (), whichColumn);
745 auto X_whichColumn = X.getVectorNonConst (whichColumn);
746 return Kokkos::subview(X_whichColumn->getLocalViewDevice(Access::ReadWrite),
751 template<
class SC,
class LO,
class GO,
class NT,
class ViewValueType>
753 copy1DViewIntoMultiVectorColumn (
755 const LO whichColumn,
756 const Kokkos::View<ViewValueType*, typename NT::device_type>& view)
758 auto X_lcl = getLocalView_1d_writeOnly (X, whichColumn);
762 template<
class SC,
class LO,
class GO,
class NT,
class ViewValueType>
764 copyMultiVectorColumnInto1DView (
765 const Kokkos::View<ViewValueType*, typename NT::device_type>& view,
767 const LO whichColumn)
769 auto X_lcl = getLocalView_1d_readOnly (X, whichColumn);
773 template<
class OneDViewType,
class IndexType>
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) {}
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;
792 template<
class OneDViewType>
793 bool findZero (
const OneDViewType& x)
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>;
800 Kokkos::RangePolicy<execution_space, size_type> range (0, x.extent (0));
801 range.set (Kokkos::ChunkSize (500));
804 Kokkos::parallel_reduce (
"findZero", range, functor_type (x), foundZero);
805 return foundZero == 1;
808 template<
class SC,
class LO,
class GO,
class NT>
810 globalizeRowOneNorms (EquilibrationInfo<
typename Kokkos::ArithTraits<SC>::val_type,
811 typename NT::device_type>& equib,
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.");
825 auto exp = G->getExporter ();
826 if (! exp.is_null ()) {
836 mv_type rowMapMV (G->getRowMap (), 2,
false);
838 copy1DViewIntoMultiVectorColumn (rowMapMV, 0, equib.rowNorms);
839 copy1DViewIntoMultiVectorColumn (rowMapMV, 1, equib.rowDiagonalEntries);
841 mv_type rangeMapMV (G->getRangeMap (), 2,
true);
845 copyMultiVectorColumnInto1DView (equib.rowNorms, rowMapMV, 0);
846 copyMultiVectorColumnInto1DView (equib.rowDiagonalEntries, rowMapMV, 1);
851 equib.foundZeroDiag = findZero (equib.rowDiagonalEntries);
852 equib.foundZeroRowNorm = findZero (equib.rowNorms);
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;
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);
870 equib.foundInf = gblNaughtyMatrix[0] == 1;
871 equib.foundNan = gblNaughtyMatrix[1] == 1;
872 equib.foundZeroDiag = gblNaughtyMatrix[2] == 1;
873 equib.foundZeroRowNorm = gblNaughtyMatrix[3] == 1;
876 template<
class SC,
class LO,
class GO,
class NT>
878 globalizeColumnOneNorms (EquilibrationInfo<
typename Kokkos::ArithTraits<SC>::val_type,
879 typename NT::device_type>& equib,
881 const bool assumeSymmetric)
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;
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.");
897 auto imp = G->getImporter ();
898 if (assumeSymmetric) {
899 const LO numCols = 2;
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);
913 mv_type rowNorms_rowMap (G->getRowMap (), numCols,
true);
914 copy1DViewIntoMultiVectorColumn (rowNorms_rowMap, 0, equib.rowNorms);
915 copy1DViewIntoMultiVectorColumn (rowNorms_rowMap, 1, equib.rowDiagonalEntries);
921 std::unique_ptr<mv_type> rowNorms_colMap;
922 if (imp.is_null ()) {
925 std::unique_ptr<mv_type> (
new mv_type (rowNorms_domMap, * (G->getColMap ())));
929 std::unique_ptr<mv_type> (
new mv_type (G->getColMap (), numCols,
true));
934 const LO lclNumCols =
935 static_cast<LO
> (G->getColMap ()->getNodeNumElements ());
936 if (
static_cast<LO
> (equib.colNorms.extent (0)) != lclNumCols) {
938 Kokkos::View<mag_type*, device_type> (
"colNorms", lclNumCols);
940 if (
static_cast<LO
> (equib.colDiagonalEntries.extent (0)) != lclNumCols) {
941 equib.colDiagonalEntries =
942 Kokkos::View<val_type*, device_type> (
"colDiagonalEntries", lclNumCols);
947 copyMultiVectorColumnInto1DView (equib.colNorms, *rowNorms_colMap, 0);
948 copyMultiVectorColumnInto1DView (equib.colDiagonalEntries, *rowNorms_colMap, 1);
951 if (! imp.is_null ()) {
952 const LO numCols = 3;
961 mv_type colMapMV (G->getColMap (), numCols,
false);
963 copy1DViewIntoMultiVectorColumn (colMapMV, 0, equib.colNorms);
964 copy1DViewIntoMultiVectorColumn (colMapMV, 1, equib.colDiagonalEntries);
965 copy1DViewIntoMultiVectorColumn (colMapMV, 2, equib.rowScaledColNorms);
967 mv_type domainMapMV (G->getDomainMap (), numCols,
true);
968 domainMapMV.doExport (colMapMV, *imp,
Tpetra::ADD);
971 copyMultiVectorColumnInto1DView (equib.colNorms, colMapMV, 0);
972 copyMultiVectorColumnInto1DView (equib.colDiagonalEntries, colMapMV, 1);
973 copyMultiVectorColumnInto1DView (equib.rowScaledColNorms, colMapMV, 2);
980 template<
class SC,
class LO,
class GO,
class NT>
981 Details::EquilibrationInfo<typename Kokkos::ArithTraits<SC>::val_type,
982 typename NT::device_type>
985 TEUCHOS_TEST_FOR_EXCEPTION
987 "computeRowOneNorms: Input matrix A must be fillComplete.");
990 Details::globalizeRowOneNorms (result, A);
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)
1000 TEUCHOS_TEST_FOR_EXCEPTION
1002 "computeRowAndColumnOneNorms: Input matrix A must be fillComplete.");
1005 Details::globalizeRowOneNorms (result, A);
1006 if (! assumeSymmetric) {
1010 Details::computeLocalRowScaledColumnNorms (result, A);
1012 Details::globalizeColumnOneNorms (result, A, assumeSymmetric);
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); \
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);
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.
Struct storing results of Tpetra::computeRowAndColumnOneNorms.
void assign(const EquilibrationInfo< ScalarType, SrcDeviceType > &src)
Deep-copy src into *this.