10 #ifndef TPETRA_DETAILS_EQUILIBRATIONINFO_HPP
11 #define TPETRA_DETAILS_EQUILIBRATIONINFO_HPP
16 #include "TpetraCore_config.h"
17 #include "Kokkos_ArithTraits.hpp"
18 #include "Kokkos_Core.hpp"
46 template<
class ScalarType,
class DeviceType>
48 using val_type =
typename Kokkos::ArithTraits<ScalarType>::val_type;
49 using mag_type =
typename Kokkos::ArithTraits<val_type>::mag_type;
50 using device_type =
typename DeviceType::device_type;
51 using host_device_type =
typename Kokkos::View<mag_type*, device_type>::HostMirror::device_type;
62 const std::size_t lclNumCols,
63 const bool assumeSymmetric_) :
64 rowNorms (Kokkos::View<mag_type*, device_type> (
"rowNorms", lclNumRows)),
65 rowDiagonalEntries (Kokkos::View<val_type*, device_type> (
"rowDiagonalEntries", lclNumRows)),
66 colNorms (Kokkos::View<mag_type*, device_type> (
"colNorms", lclNumCols)),
83 const Kokkos::View<val_type*, device_type>& rowDiagonalEntries_,
84 const Kokkos::View<mag_type*, device_type>& colNorms_,
85 const Kokkos::View<val_type*, device_type>& colDiagonalEntries_,
86 const Kokkos::View<mag_type*, device_type>& rowScaledColNorms_,
87 const bool assumeSymmetric_,
90 const bool foundZeroDiag_,
91 const bool foundZeroRowNorm_) :
105 template<
class SrcDeviceType>
109 using execution_space =
typename device_type::execution_space;
118 Kokkos::View<val_type*, device_type> (
"colDiagonalEntries", 0);
126 Kokkos::View<mag_type*, device_type> (
"rowScaledColNorms", 0);
143 auto rowNorms_h = Kokkos::create_mirror_view (
rowNorms);
145 auto colNorms_h = Kokkos::create_mirror_view (
colNorms);
149 return HostMirror {rowNorms_h, rowDiagonalEntries_h, colNorms_h,
212 #endif // TPETRA_DETAILS_EQUILIBRATIONINFO_HPP
Kokkos::View< mag_type *, device_type > colNorms
One-norms of the matrix's columns, distributed via the column Map.
Kokkos::View< mag_type *, device_type > rowNorms
One-norms of the matrix's rows, distributed via the row Map.
bool foundZeroDiag
Found a zero diagonal entry somewhere in the matrix.
bool assumeSymmetric
Whether to assume that the matrix is (globally) symmetric.
Kokkos::View< val_type *, device_type > colDiagonalEntries
Diagonal entries of the matrix, distributed via the column Map.
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.
Kokkos::View< val_type *, device_type > rowDiagonalEntries
Diagonal entries of the matrix, distributed via the row Map.
Kokkos::View< mag_type *, device_type > rowScaledColNorms
One-norms of the matrix's columns, after the matrix's rows have been scaled by rowNorms.
void assign(const EquilibrationInfo< ScalarType, SrcDeviceType > &src)
Deep-copy src into *this.
bool foundNan
Found a NaN somewhere in the matrix.
bool foundInf
Found an Inf somewhere in the matrix.
bool foundZeroRowNorm
At least one row of the matrix has a zero norm.