mirror of
https://gitlab.com/libeigen/eigen.git
synced 2025-04-19 08:09:36 +08:00
Fix incomplete cholesky.
This commit is contained in:
parent
f1adb0ccc2
commit
352ede96e4
@ -32,17 +32,19 @@ namespace Eigen {
|
|||||||
*
|
*
|
||||||
* \implsparsesolverconcept
|
* \implsparsesolverconcept
|
||||||
*
|
*
|
||||||
* It performs the following incomplete factorization: \f$ S P A P' S \approx L L' \f$
|
* It performs the following incomplete factorization: \f$ S P A P' S + \sigma I \approx L L' \f$
|
||||||
* where L is a lower triangular factor, S is a diagonal scaling matrix, and P is a
|
* where L is a lower triangular factor, S is a diagonal scaling matrix, P is a
|
||||||
* fill-in reducing permutation as computed by the ordering method.
|
* fill-in reducing permutation as computed by the ordering method, and \f$ \sigma \f$ is a shift
|
||||||
|
* for ensuring the decomposed matrix is positive definite.
|
||||||
*
|
*
|
||||||
* \b Shifting \b strategy: Let \f$ B = S P A P' S \f$ be the scaled matrix on which the factorization is carried out,
|
* \b Shifting \b strategy: Let \f$ B = S P A P' S \f$ be the scaled matrix on which the factorization is carried out,
|
||||||
* and \f$ \beta \f$ be the minimum value of the diagonal. If \f$ \beta > 0 \f$ then, the factorization is directly
|
* and \f$ \beta \f$ be the minimum value of the diagonal. If \f$ \beta > 0 \f$ then, the factorization is directly
|
||||||
* performed on the matrix B. Otherwise, the factorization is performed on the shifted matrix \f$ B + (\sigma+|\beta| I
|
* performed on the matrix B, and \sigma = 0. Otherwise, the factorization is performed on the shifted matrix \f$ B +
|
||||||
* \f$ where \f$ \sigma \f$ is the initial shift value as returned and set by setInitialShift() method. The default
|
* \sigma I \f$ for a shifting factor \f$ \sigma \f$. We start with \f$ \sigma = \sigma_0 - \beta \f$, where \f$
|
||||||
* value is \f$ \sigma = 10^{-3} \f$. If the factorization fails, then the shift in doubled until it succeed or a
|
* \sigma_0 \f$ is the initial shift value as returned and set by setInitialShift() method. The default value is \f$
|
||||||
* maximum of ten attempts. If it still fails, as returned by the info() method, then you can either increase the
|
* \sigma_0 = 10^{-3} \f$. If the factorization fails, then the shift in doubled until it succeed or a maximum of ten
|
||||||
* initial shift, or better use another preconditioning technique.
|
* attempts. If it still fails, as returned by the info() method, then you can either increase the initial shift, or
|
||||||
|
* better use another preconditioning technique.
|
||||||
*
|
*
|
||||||
*/
|
*/
|
||||||
template <typename Scalar, int UpLo_ = Lower, typename OrderingType_ = AMDOrdering<int> >
|
template <typename Scalar, int UpLo_ = Lower, typename OrderingType_ = AMDOrdering<int> >
|
||||||
@ -176,6 +178,9 @@ class IncompleteCholesky : public SparseSolverBase<IncompleteCholesky<Scalar, Up
|
|||||||
return m_perm;
|
return m_perm;
|
||||||
}
|
}
|
||||||
|
|
||||||
|
/** \returns the final shift parameter from the computation */
|
||||||
|
RealScalar shift() const { return m_shift; }
|
||||||
|
|
||||||
protected:
|
protected:
|
||||||
FactorType m_L; // The lower part stored in CSC
|
FactorType m_L; // The lower part stored in CSC
|
||||||
VectorRx m_scale; // The vector for scaling the matrix
|
VectorRx m_scale; // The vector for scaling the matrix
|
||||||
@ -184,6 +189,7 @@ class IncompleteCholesky : public SparseSolverBase<IncompleteCholesky<Scalar, Up
|
|||||||
bool m_factorizationIsOk;
|
bool m_factorizationIsOk;
|
||||||
ComputationInfo m_info;
|
ComputationInfo m_info;
|
||||||
PermutationType m_perm;
|
PermutationType m_perm;
|
||||||
|
RealScalar m_shift; // The final shift parameter.
|
||||||
|
|
||||||
private:
|
private:
|
||||||
inline void updateList(Ref<const VectorIx> colPtr, Ref<VectorIx> rowIdx, Ref<VectorSx> vals, const Index& col,
|
inline void updateList(Ref<const VectorIx> colPtr, Ref<VectorIx> rowIdx, Ref<VectorSx> vals, const Index& col,
|
||||||
@ -220,8 +226,9 @@ void IncompleteCholesky<Scalar, UpLo_, OrderingType>::factorize(const MatrixType
|
|||||||
// matrix has a zero on the diagonal.
|
// matrix has a zero on the diagonal.
|
||||||
bool modified = false;
|
bool modified = false;
|
||||||
for (Index i = 0; i < mat.cols(); ++i) {
|
for (Index i = 0; i < mat.cols(); ++i) {
|
||||||
if (numext::is_exactly_zero(m_L.coeff(i, i))) {
|
bool inserted = false;
|
||||||
m_L.insert(i, i) = Scalar(0);
|
m_L.findOrInsertCoeff(i, i, &inserted);
|
||||||
|
if (inserted) {
|
||||||
modified = true;
|
modified = true;
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
@ -270,8 +277,8 @@ void IncompleteCholesky<Scalar, UpLo_, OrderingType>::factorize(const MatrixType
|
|||||||
|
|
||||||
FactorType L_save = m_L;
|
FactorType L_save = m_L;
|
||||||
|
|
||||||
RealScalar shift = 0;
|
m_shift = RealScalar(0);
|
||||||
if (mindiag <= RealScalar(0.)) shift = m_initialShift - mindiag;
|
if (mindiag <= RealScalar(0.)) m_shift = m_initialShift - mindiag;
|
||||||
|
|
||||||
m_info = NumericalIssue;
|
m_info = NumericalIssue;
|
||||||
|
|
||||||
@ -279,7 +286,7 @@ void IncompleteCholesky<Scalar, UpLo_, OrderingType>::factorize(const MatrixType
|
|||||||
int iter = 0;
|
int iter = 0;
|
||||||
do {
|
do {
|
||||||
// Apply the shift to the diagonal elements of the matrix
|
// Apply the shift to the diagonal elements of the matrix
|
||||||
for (Index j = 0; j < n; j++) vals[colPtr[j]] += shift;
|
for (Index j = 0; j < n; j++) vals[colPtr[j]] += m_shift;
|
||||||
|
|
||||||
// jki version of the Cholesky factorization
|
// jki version of the Cholesky factorization
|
||||||
Index j = 0;
|
Index j = 0;
|
||||||
@ -323,7 +330,7 @@ void IncompleteCholesky<Scalar, UpLo_, OrderingType>::factorize(const MatrixType
|
|||||||
if (++iter >= 10) return;
|
if (++iter >= 10) return;
|
||||||
|
|
||||||
// increase shift
|
// increase shift
|
||||||
shift = numext::maxi(m_initialShift, RealScalar(2) * shift);
|
m_shift = numext::maxi(m_initialShift, RealScalar(2) * m_shift);
|
||||||
// restore m_L, col_pattern, and listCol
|
// restore m_L, col_pattern, and listCol
|
||||||
vals = Map<const VectorSx>(L_save.valuePtr(), nnz);
|
vals = Map<const VectorSx>(L_save.valuePtr(), nnz);
|
||||||
rowIdx = Map<const VectorIx>(L_save.innerIndexPtr(), nnz);
|
rowIdx = Map<const VectorIx>(L_save.innerIndexPtr(), nnz);
|
||||||
|
@ -217,15 +217,18 @@ class SparseMatrix : public SparseCompressedBase<SparseMatrix<Scalar_, Options_,
|
|||||||
return m_data.atInRange(m_outerIndex[outer], end, inner);
|
return m_data.atInRange(m_outerIndex[outer], end, inner);
|
||||||
}
|
}
|
||||||
|
|
||||||
/** \returns a non-const reference to the value of the matrix at position \a i, \a j
|
/** \returns a non-const reference to the value of the matrix at position \a i, \a j.
|
||||||
*
|
*
|
||||||
* If the element does not exist then it is inserted via the insert(Index,Index) function
|
* If the element does not exist then it is inserted via the insert(Index,Index) function
|
||||||
* which itself turns the matrix into a non compressed form if that was not the case.
|
* which itself turns the matrix into a non compressed form if that was not the case.
|
||||||
|
* The output parameter `inserted` is set to true.
|
||||||
|
*
|
||||||
|
* Otherwise, if the element does exist, `inserted` will be set to false.
|
||||||
*
|
*
|
||||||
* This is a O(log(nnz_j)) operation (binary search) plus the cost of insert(Index,Index)
|
* This is a O(log(nnz_j)) operation (binary search) plus the cost of insert(Index,Index)
|
||||||
* function if the element does not already exist.
|
* function if the element does not already exist.
|
||||||
*/
|
*/
|
||||||
inline Scalar& coeffRef(Index row, Index col) {
|
inline Scalar& findOrInsertCoeff(Index row, Index col, bool* inserted) {
|
||||||
eigen_assert(row >= 0 && row < rows() && col >= 0 && col < cols());
|
eigen_assert(row >= 0 && row < rows() && col >= 0 && col < cols());
|
||||||
const Index outer = IsRowMajor ? row : col;
|
const Index outer = IsRowMajor ? row : col;
|
||||||
const Index inner = IsRowMajor ? col : row;
|
const Index inner = IsRowMajor ? col : row;
|
||||||
@ -240,17 +243,37 @@ class SparseMatrix : public SparseCompressedBase<SparseMatrix<Scalar_, Options_,
|
|||||||
m_innerNonZeros[outer]++;
|
m_innerNonZeros[outer]++;
|
||||||
m_data.index(end) = StorageIndex(inner);
|
m_data.index(end) = StorageIndex(inner);
|
||||||
m_data.value(end) = Scalar(0);
|
m_data.value(end) = Scalar(0);
|
||||||
|
if (inserted != nullptr) {
|
||||||
|
*inserted = true;
|
||||||
|
}
|
||||||
return m_data.value(end);
|
return m_data.value(end);
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
if ((dst < end) && (m_data.index(dst) == inner))
|
if ((dst < end) && (m_data.index(dst) == inner)) {
|
||||||
// this coefficient exists, return a refernece to it
|
// this coefficient exists, return a refernece to it
|
||||||
|
if (inserted != nullptr) {
|
||||||
|
*inserted = false;
|
||||||
|
}
|
||||||
return m_data.value(dst);
|
return m_data.value(dst);
|
||||||
else
|
} else {
|
||||||
|
if (inserted != nullptr) {
|
||||||
|
*inserted = true;
|
||||||
|
}
|
||||||
// insertion will require reconfiguring the buffer
|
// insertion will require reconfiguring the buffer
|
||||||
return insertAtByOuterInner(outer, inner, dst);
|
return insertAtByOuterInner(outer, inner, dst);
|
||||||
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
|
/** \returns a non-const reference to the value of the matrix at position \a i, \a j
|
||||||
|
*
|
||||||
|
* If the element does not exist then it is inserted via the insert(Index,Index) function
|
||||||
|
* which itself turns the matrix into a non compressed form if that was not the case.
|
||||||
|
*
|
||||||
|
* This is a O(log(nnz_j)) operation (binary search) plus the cost of insert(Index,Index)
|
||||||
|
* function if the element does not already exist.
|
||||||
|
*/
|
||||||
|
inline Scalar& coeffRef(Index row, Index col) { return findOrInsertCoeff(row, col, nullptr); }
|
||||||
|
|
||||||
/** \returns a reference to a novel non zero coefficient with coordinates \a row x \a col.
|
/** \returns a reference to a novel non zero coefficient with coordinates \a row x \a col.
|
||||||
* The non zero coefficient must \b not already exist.
|
* The non zero coefficient must \b not already exist.
|
||||||
*
|
*
|
||||||
|
@ -54,10 +54,28 @@ void bug1150() {
|
|||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
|
void test_non_spd() {
|
||||||
|
Eigen::SparseMatrix<double> A(2, 2);
|
||||||
|
A.insert(0, 0) = 0;
|
||||||
|
A.insert(1, 1) = 3;
|
||||||
|
|
||||||
|
Eigen::IncompleteCholesky<double> solver(A);
|
||||||
|
|
||||||
|
// Recover original matrix.
|
||||||
|
Eigen::MatrixXd M = solver.permutationP().transpose() *
|
||||||
|
(solver.scalingS().asDiagonal().inverse() *
|
||||||
|
(solver.matrixL() * solver.matrixL().transpose() -
|
||||||
|
solver.shift() * Eigen::MatrixXd::Identity(A.rows(), A.cols())) *
|
||||||
|
solver.scalingS().asDiagonal().inverse()) *
|
||||||
|
solver.permutationP();
|
||||||
|
VERIFY_IS_APPROX(A.toDense(), M);
|
||||||
|
}
|
||||||
|
|
||||||
EIGEN_DECLARE_TEST(incomplete_cholesky) {
|
EIGEN_DECLARE_TEST(incomplete_cholesky) {
|
||||||
CALL_SUBTEST_1((test_incomplete_cholesky_T<double, int>()));
|
CALL_SUBTEST_1((test_incomplete_cholesky_T<double, int>()));
|
||||||
CALL_SUBTEST_2((test_incomplete_cholesky_T<std::complex<double>, int>()));
|
CALL_SUBTEST_2((test_incomplete_cholesky_T<std::complex<double>, int>()));
|
||||||
CALL_SUBTEST_3((test_incomplete_cholesky_T<double, long int>()));
|
CALL_SUBTEST_3((test_incomplete_cholesky_T<double, long int>()));
|
||||||
|
|
||||||
CALL_SUBTEST_1((bug1150<0>()));
|
CALL_SUBTEST_4((bug1150<0>()));
|
||||||
|
CALL_SUBTEST_4(test_non_spd());
|
||||||
}
|
}
|
||||||
|
Loading…
x
Reference in New Issue
Block a user