aboutsummaryrefslogtreecommitdiffhomepage
path: root/Eigen/src/LU
diff options
context:
space:
mode:
authorGravatar Benoit Jacob <jacob.benoit.1@gmail.com>2010-06-11 08:04:06 -0400
committerGravatar Benoit Jacob <jacob.benoit.1@gmail.com>2010-06-11 08:04:06 -0400
commitd72d538747d6c46f37bbcba2f0484414f10cb059 (patch)
treef64c6eeb7c33af9763e3e382197b233329bd132c /Eigen/src/LU
parentbdd7c6c88a0b8cb931480e04e33a17aa08022e06 (diff)
parent00267e3a471a10e842f771de474f0dca6407a693 (diff)
merge my Dynamic -> -1 change
Diffstat (limited to 'Eigen/src/LU')
-rw-r--r--Eigen/src/LU/FullPivLU.h2
-rw-r--r--Eigen/src/LU/Inverse.h14
-rw-r--r--Eigen/src/LU/PartialPivLU.h26
3 files changed, 25 insertions, 17 deletions
diff --git a/Eigen/src/LU/FullPivLU.h b/Eigen/src/LU/FullPivLU.h
index 32ac99d16..114cd4c7a 100644
--- a/Eigen/src/LU/FullPivLU.h
+++ b/Eigen/src/LU/FullPivLU.h
@@ -69,7 +69,7 @@ template<typename _MatrixType> class FullPivLU
typedef typename MatrixType::Scalar Scalar;
typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
typedef typename ei_traits<MatrixType>::StorageKind StorageKind;
- typedef typename ei_index<StorageKind>::type Index;
+ typedef typename MatrixType::Index Index;
typedef typename ei_plain_row_type<MatrixType, Index>::type IntRowVectorType;
typedef typename ei_plain_col_type<MatrixType, Index>::type IntColVectorType;
typedef PermutationMatrix<ColsAtCompileTime, MaxColsAtCompileTime> PermutationQType;
diff --git a/Eigen/src/LU/Inverse.h b/Eigen/src/LU/Inverse.h
index ed1724dda..5ccc59a32 100644
--- a/Eigen/src/LU/Inverse.h
+++ b/Eigen/src/LU/Inverse.h
@@ -284,7 +284,6 @@ struct ei_inverse_impl : public ReturnByValue<ei_inverse_impl<MatrixType> >
typedef typename MatrixType::Index Index;
typedef typename ei_eval<MatrixType>::type MatrixTypeNested;
typedef typename ei_cleantype<MatrixTypeNested>::type MatrixTypeNestedCleaned;
-
const MatrixTypeNested m_matrix;
ei_inverse_impl(const MatrixType& matrix)
@@ -296,6 +295,10 @@ struct ei_inverse_impl : public ReturnByValue<ei_inverse_impl<MatrixType> >
template<typename Dest> inline void evalTo(Dest& dst) const
{
+ const int Size = EIGEN_ENUM_MIN(MatrixType::ColsAtCompileTime,Dest::ColsAtCompileTime);
+ ei_assert(( (Size<=1) || (Size>4) || (ei_extract_data(m_matrix)!=ei_extract_data(dst)))
+ && "Aliasing problem detected in inverse(), you need to do inverse().eval() here.");
+
ei_compute_inverse<MatrixTypeNestedCleaned, Dest>::run(m_matrix, dst);
}
};
@@ -354,7 +357,14 @@ inline void MatrixBase<Derived>::computeInverseAndDetWithCheck(
{
// i'd love to put some static assertions there, but SFINAE means that they have no effect...
ei_assert(rows() == cols());
- ei_compute_inverse_and_det_with_check<PlainObject, ResultType>::run
+ // for 2x2, it's worth giving a chance to avoid evaluating.
+ // for larger sizes, evaluating has negligible cost and limits code size.
+ typedef typename ei_meta_if<
+ RowsAtCompileTime == 2,
+ typename ei_cleantype<typename ei_nested<Derived, 2>::type>::type,
+ PlainObject
+ >::ret MatrixType;
+ ei_compute_inverse_and_det_with_check<MatrixType, ResultType>::run
(derived(), absDeterminantThreshold, inverse, determinant, invertible);
}
diff --git a/Eigen/src/LU/PartialPivLU.h b/Eigen/src/LU/PartialPivLU.h
index 39c348e5e..7e7516a5c 100644
--- a/Eigen/src/LU/PartialPivLU.h
+++ b/Eigen/src/LU/PartialPivLU.h
@@ -72,9 +72,9 @@ template<typename _MatrixType> class PartialPivLU
typedef typename MatrixType::Scalar Scalar;
typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
typedef typename ei_traits<MatrixType>::StorageKind StorageKind;
- typedef typename ei_index<StorageKind>::type Index;
- typedef typename ei_plain_col_type<MatrixType, Index>::type PermutationVectorType;
+ typedef typename MatrixType::Index Index;
typedef PermutationMatrix<RowsAtCompileTime, MaxRowsAtCompileTime> PermutationType;
+ typedef Transpositions<RowsAtCompileTime, MaxRowsAtCompileTime> TranspositionType;
/**
@@ -186,7 +186,7 @@ template<typename _MatrixType> class PartialPivLU
protected:
MatrixType m_lu;
PermutationType m_p;
- PermutationVectorType m_rowsTranspositions;
+ TranspositionType m_rowsTranspositions;
Index m_det_p;
bool m_isInitialized;
};
@@ -339,9 +339,9 @@ struct ei_partial_lu_impl
Index tsize = size - k - bs; // trailing size
// partition the matrix:
- // A00 | A01 | A02
- // lu = A10 | A11 | A12
- // A20 | A21 | A22
+ // A00 | A01 | A02
+ // lu = A_0 | A_1 | A_2 = A10 | A11 | A12
+ // A20 | A21 | A22
BlockType A_0(lu,0,0,rows,k);
BlockType A_2(lu,0,k+bs,rows,tsize);
BlockType A11(lu,k,k,bs,bs);
@@ -350,8 +350,8 @@ struct ei_partial_lu_impl
BlockType A22(lu,k+bs,k+bs,trows,tsize);
Index nb_transpositions_in_panel;
- // recursively calls the blocked LU algorithm with a very small
- // blocking size:
+ // recursively call the blocked LU algorithm on [A11^T A21^T]^T
+ // with a very small blocking size:
if(!blocked_lu(trows+bs, bs, &lu.coeffRef(k,k), luStride,
row_transpositions+k, nb_transpositions_in_panel, 16))
{
@@ -364,7 +364,7 @@ struct ei_partial_lu_impl
}
nb_transpositions += nb_transpositions_in_panel;
- // update permutations and apply them to A10
+ // update permutations and apply them to A_0
for(Index i=k; i<k+bs; ++i)
{
Index piv = (row_transpositions[i] += k);
@@ -389,8 +389,8 @@ struct ei_partial_lu_impl
/** \internal performs the LU decomposition with partial pivoting in-place.
*/
-template<typename MatrixType, typename IntVector>
-void ei_partial_lu_inplace(MatrixType& lu, IntVector& row_transpositions, typename MatrixType::Index& nb_transpositions)
+template<typename MatrixType, typename TranspositionType>
+void ei_partial_lu_inplace(MatrixType& lu, TranspositionType& row_transpositions, typename MatrixType::Index& nb_transpositions)
{
ei_assert(lu.cols() == row_transpositions.size());
ei_assert((&row_transpositions.coeffRef(1)-&row_transpositions.coeffRef(0)) == 1);
@@ -414,9 +414,7 @@ PartialPivLU<MatrixType>& PartialPivLU<MatrixType>::compute(const MatrixType& ma
ei_partial_lu_inplace(m_lu, m_rowsTranspositions, nb_transpositions);
m_det_p = (nb_transpositions%2) ? -1 : 1;
- m_p.setIdentity(size);
- for(Index k = size-1; k >= 0; --k)
- m_p.applyTranspositionOnTheRight(k, m_rowsTranspositions.coeff(k));
+ m_p = m_rowsTranspositions;
m_isInitialized = true;
return *this;