Tpetra parallel linear algebra
Version of the Day
|
Struct storing results of Tpetra::computeRowAndColumnOneNorms. More...
#include <Tpetra_Details_EquilibrationInfo.hpp>
Public Member Functions | |
template<class SrcDeviceType > | |
void | assign (const EquilibrationInfo< ScalarType, SrcDeviceType > &src) |
Deep-copy src into *this. More... | |
Public Attributes | |
Kokkos::View< mag_type *, device_type > | rowNorms |
One-norms of the matrix's rows, distributed via the row Map. More... | |
Kokkos::View< val_type *, device_type > | rowDiagonalEntries |
Diagonal entries of the matrix, distributed via the row Map. More... | |
Kokkos::View< mag_type *, device_type > | colNorms |
One-norms of the matrix's columns, distributed via the column Map. More... | |
Kokkos::View< val_type *, device_type > | colDiagonalEntries |
Diagonal entries of the matrix, distributed via the column Map. More... | |
Kokkos::View< mag_type *, device_type > | rowScaledColNorms |
One-norms of the matrix's columns, after the matrix's rows have been scaled by rowNorms. More... | |
bool | assumeSymmetric |
Whether to assume that the matrix is (globally) symmetric. More... | |
bool | foundInf |
Found an Inf somewhere in the matrix. More... | |
bool | foundNan |
Found a NaN somewhere in the matrix. More... | |
bool | foundZeroDiag |
Found a zero diagonal entry somewhere in the matrix. More... | |
bool | foundZeroRowNorm |
At least one row of the matrix has a zero norm. More... | |
Struct storing results of Tpetra::computeRowAndColumnOneNorms.
ScalarType | Type of the entries in the Tpetra::RowMatrix / Tpetra::CrsMatrix. This must be the Kokkos-ized version, that is, Kokkos::ArithTraits<ScalarType>::val_type . |
DeviceType | Kokkos::Device specialization. |
Tpetra users may NOT do anything with this struct other than get it as the return value of Tpetra::computeRowAndColumnOneNorms, and pass it into Tpetra::leftAndOrRightScaleCrsMatrix.
Tpetra developers may use the fields in this struct, to implement equilibration, balancing (in the symmetric / Hermitian positive definite case only), or to analyze the matrix (e.g., whether it is diagonally dominant).
If this is the return value of computeRowAndColumnOneNorms, results are always global. computeLocalRowAndColumnOneNorms only computes local results, where "local" means "to an MPI process." globalizeRowOneNorms and globalizeColumnOneNorms "globalize" the results, so they refer to the entire matrix, globally distributed over an MPI communicator.
Definition at line 51 of file Tpetra_Details_EquilibrationInfo.hpp.
|
inline |
Deep-copy src into *this.
Definition at line 112 of file Tpetra_Details_EquilibrationInfo.hpp.
Kokkos::View<mag_type*, device_type> Tpetra::Details::EquilibrationInfo< ScalarType, DeviceType >::rowNorms |
One-norms of the matrix's rows, distributed via the row Map.
Definition at line 162 of file Tpetra_Details_EquilibrationInfo.hpp.
Kokkos::View<val_type*, device_type> Tpetra::Details::EquilibrationInfo< ScalarType, DeviceType >::rowDiagonalEntries |
Diagonal entries of the matrix, distributed via the row Map.
Definition at line 165 of file Tpetra_Details_EquilibrationInfo.hpp.
Kokkos::View<mag_type*, device_type> Tpetra::Details::EquilibrationInfo< ScalarType, DeviceType >::colNorms |
One-norms of the matrix's columns, distributed via the column Map.
If assumeSymmetric is true, this is just a redistributed version of rowNorms.
Definition at line 172 of file Tpetra_Details_EquilibrationInfo.hpp.
Kokkos::View<val_type*, device_type> Tpetra::Details::EquilibrationInfo< ScalarType, DeviceType >::colDiagonalEntries |
Diagonal entries of the matrix, distributed via the column Map.
Only use this if assumeSymmetric is false.
Definition at line 177 of file Tpetra_Details_EquilibrationInfo.hpp.
Kokkos::View<mag_type*, device_type> Tpetra::Details::EquilibrationInfo< ScalarType, DeviceType >::rowScaledColNorms |
One-norms of the matrix's columns, after the matrix's rows have been scaled by rowNorms.
This is only valid if assumeSymmetric is false.
For the nonsymmetric case, we imitate LAPACK's DGEEQU, in doing the row scaling first, then the column scaling. Thus, the column norms are "scaled" by the row norms (above). We still keep the unscaled column norms (colNorms; see above) in that case, because they are diagnostic.
Definition at line 189 of file Tpetra_Details_EquilibrationInfo.hpp.
bool Tpetra::Details::EquilibrationInfo< ScalarType, DeviceType >::assumeSymmetric |
Whether to assume that the matrix is (globally) symmetric.
This affects whether colDiagonalEntries and rowScaledColNorms are valid.
Definition at line 195 of file Tpetra_Details_EquilibrationInfo.hpp.
bool Tpetra::Details::EquilibrationInfo< ScalarType, DeviceType >::foundInf |
Found an Inf somewhere in the matrix.
Definition at line 198 of file Tpetra_Details_EquilibrationInfo.hpp.
bool Tpetra::Details::EquilibrationInfo< ScalarType, DeviceType >::foundNan |
Found a NaN somewhere in the matrix.
Definition at line 201 of file Tpetra_Details_EquilibrationInfo.hpp.
bool Tpetra::Details::EquilibrationInfo< ScalarType, DeviceType >::foundZeroDiag |
Found a zero diagonal entry somewhere in the matrix.
Definition at line 204 of file Tpetra_Details_EquilibrationInfo.hpp.
bool Tpetra::Details::EquilibrationInfo< ScalarType, DeviceType >::foundZeroRowNorm |
At least one row of the matrix has a zero norm.
Definition at line 207 of file Tpetra_Details_EquilibrationInfo.hpp.