aboutsummaryrefslogtreecommitdiffhomepage
path: root/Eigen/src/Cholesky
diff options
context:
space:
mode:
authorGravatar Rasmus Munk Larsen <rmlarsen@google.com>2016-04-04 14:20:01 -0700
committerGravatar Rasmus Munk Larsen <rmlarsen@google.com>2016-04-04 14:20:01 -0700
commit86e0ed81f8db5a0c9562b62a67a9ba60ec58dec0 (patch)
tree10d4e195015009b69a76bff39bc824f0e1e8fb09 /Eigen/src/Cholesky
parent30242b75653fa4128181dba364f540184beff5ac (diff)
Addresses comments on Eigen pull request PR-174.
* Get rid of code-duplication for real vs. complex matrices. * Fix flipped arguments to select. * Make the condition estimation functions free functions. * Use Vector::Unit() to generate canonical unit vectors. * Misc. cleanup.
Diffstat (limited to 'Eigen/src/Cholesky')
-rw-r--r--Eigen/src/Cholesky/LDLT.h12
-rw-r--r--Eigen/src/Cholesky/LLT.h12
2 files changed, 18 insertions, 6 deletions
diff --git a/Eigen/src/Cholesky/LDLT.h b/Eigen/src/Cholesky/LDLT.h
index 9753c84d8..80f18977c 100644
--- a/Eigen/src/Cholesky/LDLT.h
+++ b/Eigen/src/Cholesky/LDLT.h
@@ -198,7 +198,7 @@ template<typename _MatrixType, int _UpLo> class LDLT
RealScalar rcond() const
{
eigen_assert(m_isInitialized && "LDLT is not initialized.");
- return ConditionEstimator<LDLT<MatrixType, UpLo>, true >::rcond(m_l1_norm, *this);
+ return ReciprocalConditionNumberEstimate(m_l1_norm, *this);
}
template <typename Derived>
@@ -216,6 +216,12 @@ template<typename _MatrixType, int _UpLo> class LDLT
MatrixType reconstructedMatrix() const;
+ /** \returns the decomposition itself to allow generic code to do
+ * ldlt.transpose().solve(rhs).
+ */
+ const LDLT<MatrixType, UpLo>& transpose() const { return *this; };
+ const LDLT<MatrixType, UpLo>& adjoint() const { return *this; };
+
inline Index rows() const { return m_matrix.rows(); }
inline Index cols() const { return m_matrix.cols(); }
@@ -454,14 +460,14 @@ LDLT<MatrixType,_UpLo>& LDLT<MatrixType,_UpLo>::compute(const EigenBase<InputTyp
if (_UpLo == Lower) {
for (int col = 0; col < size; ++col) {
const RealScalar abs_col_sum = m_matrix.col(col).tail(size - col).cwiseAbs().sum() +
- m_matrix.row(col).tail(col).cwiseAbs().sum();
+ m_matrix.row(col).head(col).cwiseAbs().sum();
if (abs_col_sum > m_l1_norm) {
m_l1_norm = abs_col_sum;
}
}
} else {
for (int col = 0; col < a.cols(); ++col) {
- const RealScalar abs_col_sum = m_matrix.col(col).tail(col).cwiseAbs().sum() +
+ const RealScalar abs_col_sum = m_matrix.col(col).head(col).cwiseAbs().sum() +
m_matrix.row(col).tail(size - col).cwiseAbs().sum();
if (abs_col_sum > m_l1_norm) {
m_l1_norm = abs_col_sum;
diff --git a/Eigen/src/Cholesky/LLT.h b/Eigen/src/Cholesky/LLT.h
index b55c5bebf..94da1d52d 100644
--- a/Eigen/src/Cholesky/LLT.h
+++ b/Eigen/src/Cholesky/LLT.h
@@ -142,7 +142,7 @@ template<typename _MatrixType, int _UpLo> class LLT
{
eigen_assert(m_isInitialized && "LLT is not initialized.");
eigen_assert(m_info == Success && "LLT failed because matrix appears to be negative");
- return ConditionEstimator<LLT<MatrixType, UpLo>, true >::rcond(m_l1_norm, *this);
+ return ReciprocalConditionNumberEstimate(m_l1_norm, *this);
}
/** \returns the LLT decomposition matrix
@@ -169,6 +169,12 @@ template<typename _MatrixType, int _UpLo> class LLT
return m_info;
}
+ /** \returns the decomposition itself to allow generic code to do
+ * llt.transpose().solve(rhs).
+ */
+ const LLT<MatrixType, UpLo>& transpose() const { return *this; };
+ const LLT<MatrixType, UpLo>& adjoint() const { return *this; };
+
inline Index rows() const { return m_matrix.rows(); }
inline Index cols() const { return m_matrix.cols(); }
@@ -409,14 +415,14 @@ LLT<MatrixType,_UpLo>& LLT<MatrixType,_UpLo>::compute(const EigenBase<InputType>
if (_UpLo == Lower) {
for (int col = 0; col < size; ++col) {
const RealScalar abs_col_sum = m_matrix.col(col).tail(size - col).cwiseAbs().sum() +
- m_matrix.row(col).tail(col).cwiseAbs().sum();
+ m_matrix.row(col).head(col).cwiseAbs().sum();
if (abs_col_sum > m_l1_norm) {
m_l1_norm = abs_col_sum;
}
}
} else {
for (int col = 0; col < a.cols(); ++col) {
- const RealScalar abs_col_sum = m_matrix.col(col).tail(col).cwiseAbs().sum() +
+ const RealScalar abs_col_sum = m_matrix.col(col).head(col).cwiseAbs().sum() +
m_matrix.row(col).tail(size - col).cwiseAbs().sum();
if (abs_col_sum > m_l1_norm) {
m_l1_norm = abs_col_sum;