template<typename MatrixType_>
class Eigen::IterScaling< MatrixType_ >
iterative scaling algorithm to equilibrate rows and column norms in matrices
This class can be used as a preprocessing tool to accelerate the convergence of iterative methods
This feature is useful to limit the pivoting amount during LU/ILU factorization The scaling strategy as presented here preserves the symmetry of the problem NOTE It is assumed that the matrix does not have empty row or column,
Example with key steps
IterScaling<SparseMatrix<double> >
scal;
b =
scal.LeftScaling().cwiseProduct(
b);
x =
scal.RightScaling().cwiseProduct(
x);
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
Scalar * b
Definition: benchVecAdd.cpp:17
Matrix< SCALARA, Dynamic, Dynamic, opt_A > A
Definition: bench_gemm.cpp:47
EIGEN_BLAS_FUNC() scal(int *n, RealScalar *palpha, RealScalar *px, int *incx)
Definition: level1_impl.h:105
list x
Definition: plotDoE.py:28
- Template Parameters
-
MatrixType_ | the type of the matrix. It should be a real square sparsematrix |
References : D. Ruiz and B. Ucar, A Symmetry Preserving Algorithm for Matrix Scaling, INRIA Research report RR-7552
- See also
- IncompleteLUT
template<typename MatrixType_ >
Compute the left and right diagonal matrices to scale the input matrix mat
FIXME This algorithm will be modified such that the diagonal elements are permuted on the diagonal.
- See also
- LeftScaling() RightScaling()
84 VectorXd Dr, Dc, DrRes, DcRes;
89 double EpsRow = 1.0, EpsCol = 1.0;
96 for (
typename MatrixType::InnerIterator it(
m_matrix,
k); it; ++it) {
97 if (Dr(it.row()) <
abs(it.value())) Dr(it.row()) =
abs(it.value());
99 if (Dc(it.col()) <
abs(it.value())) Dc(it.col()) =
abs(it.value());
102 for (
int i = 0;
i <
m; ++
i) {
105 for (
int i = 0;
i <
n; ++
i) {
109 for (
int i = 0;
i <
m; ++
i) {
112 for (
int i = 0;
i <
n; ++
i) {
119 for (
typename MatrixType::InnerIterator it(
m_matrix,
k); it; ++it) {
120 it.valueRef() = it.value() / (Dr(it.row()) * Dc(it.col()));
122 if (DrRes(it.row()) <
abs(it.value())) DrRes(it.row()) =
abs(it.value());
124 if (DcRes(it.col()) <
abs(it.value())) DcRes(it.col()) =
abs(it.value());
127 DrRes.array() = (1 - DrRes.array()).
abs();
128 EpsRow = DrRes.maxCoeff();
129 DcRes.array() = (1 - DcRes.array()).
abs();
130 EpsCol = DcRes.maxCoeff();
AnnoyingScalar abs(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:135
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Eigen::SparseMatrix< double > mat
Definition: EigenUnitTest.cpp:10
#define eigen_assert(x)
Definition: Macros.h:910
VectorXd m_left
Definition: unsupported/Eigen/src/IterativeSolvers/Scaling.h:166
MatrixType m_matrix
Definition: unsupported/Eigen/src/IterativeSolvers/Scaling.h:163
bool m_isInitialized
Definition: unsupported/Eigen/src/IterativeSolvers/Scaling.h:165
VectorXd m_right
Definition: unsupported/Eigen/src/IterativeSolvers/Scaling.h:167
int m_maxits
Definition: unsupported/Eigen/src/IterativeSolvers/Scaling.h:169
double m_tol
Definition: unsupported/Eigen/src/IterativeSolvers/Scaling.h:168
Index cols() const
Definition: SparseMatrix.h:161
void resize(Index rows, Index cols)
Definition: SparseMatrix.h:734
Index rows() const
Definition: SparseMatrix.h:159
int * m
Definition: level2_cplx_impl.h:294
char char char int int * k
Definition: level2_impl.h:374
References abs(), Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::cols(), eigen_assert, i, k, m, Eigen::IterScaling< MatrixType_ >::m_isInitialized, Eigen::IterScaling< MatrixType_ >::m_left, Eigen::IterScaling< MatrixType_ >::m_matrix, Eigen::IterScaling< MatrixType_ >::m_maxits, Eigen::IterScaling< MatrixType_ >::m_right, Eigen::IterScaling< MatrixType_ >::m_tol, n, Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::resize(), Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::rows(), and sqrt().
Referenced by Eigen::IterScaling< MatrixType_ >::computeRef(), and Eigen::IterScaling< MatrixType_ >::IterScaling().