10 #ifndef TPETRA_DETAILS_EQUILIBRATIONINFO_HPP
11 #define TPETRA_DETAILS_EQUILIBRATIONINFO_HPP
16 #include "TpetraCore_config.h"
17 #if KOKKOS_VERSION >= 40799
18 #include "KokkosKernels_ArithTraits.hpp"
20 #include "Kokkos_ArithTraits.hpp"
22 #include "Kokkos_Core.hpp"
50 template <
class ScalarType,
class DeviceType>
52 #if KOKKOS_VERSION >= 40799
53 using val_type =
typename KokkosKernels::ArithTraits<ScalarType>::val_type;
55 using val_type =
typename Kokkos::ArithTraits<ScalarType>::val_type;
57 #if KOKKOS_VERSION >= 40799
58 using mag_type =
typename KokkosKernels::ArithTraits<val_type>::mag_type;
60 using mag_type =
typename Kokkos::ArithTraits<val_type>::mag_type;
62 using device_type =
typename DeviceType::device_type;
63 using host_device_type =
typename Kokkos::View<mag_type*, device_type>::host_mirror_type::device_type;
73 const std::size_t lclNumCols,
74 const bool assumeSymmetric_)
75 :
rowNorms(Kokkos::View<mag_type*, device_type>(
"rowNorms", lclNumRows))
76 ,
rowDiagonalEntries(Kokkos::View<val_type*, device_type>(
"rowDiagonalEntries", lclNumRows))
77 ,
colNorms(Kokkos::View<mag_type*, device_type>(
"colNorms", lclNumCols))
79 assumeSymmetric_ ? std::size_t(0) : lclNumCols))
81 assumeSymmetric_ ? std::size_t(0) : lclNumCols))
89 const Kokkos::View<val_type*, device_type>& rowDiagonalEntries_,
90 const Kokkos::View<mag_type*, device_type>& colNorms_,
91 const Kokkos::View<val_type*, device_type>& colDiagonalEntries_,
92 const Kokkos::View<mag_type*, device_type>& rowScaledColNorms_,
93 const bool assumeSymmetric_,
96 const bool foundZeroDiag_,
97 const bool foundZeroRowNorm_)
110 template <
class SrcDeviceType>
113 using execution_space =
typename device_type::execution_space;
122 Kokkos::View<val_type*, device_type>(
"colDiagonalEntries", 0);
129 Kokkos::View<mag_type*, device_type>(
"rowScaledColNorms", 0);
144 auto rowNorms_h = Kokkos::create_mirror_view(
rowNorms);
146 auto colNorms_h = Kokkos::create_mirror_view(
colNorms);
150 return HostMirror{rowNorms_h, rowDiagonalEntries_h, colNorms_h,
213 #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.