aboutsummaryrefslogtreecommitdiffhomepage
path: root/Eigen/src/LU
diff options
context:
space:
mode:
authorGravatar Benoit Jacob <jacob.benoit.1@gmail.com>2009-11-09 07:51:31 -0500
committerGravatar Benoit Jacob <jacob.benoit.1@gmail.com>2009-11-09 07:51:31 -0500
commit9a0900e33e9ca4bc174cfccc26dbf4bdc38f7627 (patch)
tree5d24519203764f0c6b6dc479946e56fd92025cb2 /Eigen/src/LU
parente4e58e8337e82ba76f6bf4fe7000acac9337056c (diff)
last round of changes, mainly to return derived types instead of base types, and fix various compilation issues
Diffstat (limited to 'Eigen/src/LU')
-rw-r--r--Eigen/src/LU/FullPivLU.h28
-rw-r--r--Eigen/src/LU/PartialPivLU.h12
2 files changed, 20 insertions, 20 deletions
diff --git a/Eigen/src/LU/FullPivLU.h b/Eigen/src/LU/FullPivLU.h
index c38789bd9..d8975b0b6 100644
--- a/Eigen/src/LU/FullPivLU.h
+++ b/Eigen/src/LU/FullPivLU.h
@@ -158,10 +158,10 @@ template<typename _MatrixType> class FullPivLU
*
* \sa image()
*/
- inline const ei_kernel_return_value<FullPivLU> kernel() const
+ inline const ei_kernel_retval<FullPivLU> kernel() const
{
ei_assert(m_isInitialized && "LU is not initialized.");
- return ei_kernel_return_value<FullPivLU>(*this);
+ return ei_kernel_retval<FullPivLU>(*this);
}
/** \returns the image of the matrix, also called its column-space. The columns of the returned matrix
@@ -183,11 +183,11 @@ template<typename _MatrixType> class FullPivLU
*
* \sa kernel()
*/
- inline const ei_image_return_value<FullPivLU>
+ inline const ei_image_retval<FullPivLU>
image(const MatrixType& originalMatrix) const
{
ei_assert(m_isInitialized && "LU is not initialized.");
- return ei_image_return_value<FullPivLU>(*this, originalMatrix);
+ return ei_image_retval<FullPivLU>(*this, originalMatrix);
}
/** \return a solution x to the equation Ax=b, where A is the matrix of which
@@ -210,11 +210,11 @@ template<typename _MatrixType> class FullPivLU
* \sa TriangularView::solve(), kernel(), inverse()
*/
template<typename Rhs>
- inline const ei_solve_return_value<FullPivLU, Rhs>
+ inline const ei_solve_retval<FullPivLU, Rhs>
solve(const MatrixBase<Rhs>& b) const
{
ei_assert(m_isInitialized && "LU is not initialized.");
- return ei_solve_return_value<FullPivLU, Rhs>(*this, b.derived());
+ return ei_solve_retval<FullPivLU, Rhs>(*this, b.derived());
}
/** \returns the determinant of the matrix of which
@@ -355,11 +355,11 @@ template<typename _MatrixType> class FullPivLU
*
* \sa MatrixBase::inverse()
*/
- inline const ei_solve_return_value<FullPivLU,NestByValue<typename MatrixType::IdentityReturnType> > inverse() const
+ inline const ei_solve_retval<FullPivLU,NestByValue<typename MatrixType::IdentityReturnType> > inverse() const
{
ei_assert(m_isInitialized && "LU is not initialized.");
ei_assert(m_lu.rows() == m_lu.cols() && "You can't take the inverse of a non-square matrix!");
- return ei_solve_return_value<FullPivLU,NestByValue<typename MatrixType::IdentityReturnType> >
+ return ei_solve_retval<FullPivLU,NestByValue<typename MatrixType::IdentityReturnType> >
(*this, MatrixType::Identity(m_lu.rows(), m_lu.cols()).nestByValue());
}
@@ -486,8 +486,8 @@ typename ei_traits<MatrixType>::Scalar FullPivLU<MatrixType>::determinant() cons
/********* Implementation of kernel() **************************************************/
template<typename _MatrixType>
-struct ei_kernel_impl<FullPivLU<_MatrixType> >
- : ei_kernel_return_value<FullPivLU<_MatrixType> >
+struct ei_kernel_retval<FullPivLU<_MatrixType> >
+ : ei_kernel_retval_base<FullPivLU<_MatrixType> >
{
EIGEN_MAKE_KERNEL_HELPERS(FullPivLU<_MatrixType>)
@@ -571,8 +571,8 @@ struct ei_kernel_impl<FullPivLU<_MatrixType> >
/***** Implementation of image() *****************************************************/
template<typename _MatrixType>
-struct ei_image_impl<FullPivLU<_MatrixType> >
- : ei_image_return_value<FullPivLU<_MatrixType> >
+struct ei_image_retval<FullPivLU<_MatrixType> >
+ : ei_image_retval_base<FullPivLU<_MatrixType> >
{
EIGEN_MAKE_IMAGE_HELPERS(FullPivLU<_MatrixType>)
@@ -608,8 +608,8 @@ struct ei_image_impl<FullPivLU<_MatrixType> >
/***** Implementation of solve() *****************************************************/
template<typename _MatrixType, typename Rhs>
-struct ei_solve_impl<FullPivLU<_MatrixType>, Rhs>
- : ei_solve_return_value<FullPivLU<_MatrixType>, Rhs>
+struct ei_solve_retval<FullPivLU<_MatrixType>, Rhs>
+ : ei_solve_retval_base<FullPivLU<_MatrixType>, Rhs>
{
EIGEN_MAKE_SOLVE_HELPERS(FullPivLU<_MatrixType>,Rhs)
diff --git a/Eigen/src/LU/PartialPivLU.h b/Eigen/src/LU/PartialPivLU.h
index 4f6cda68f..03c91456a 100644
--- a/Eigen/src/LU/PartialPivLU.h
+++ b/Eigen/src/LU/PartialPivLU.h
@@ -133,11 +133,11 @@ template<typename _MatrixType> class PartialPivLU
* \sa TriangularView::solve(), inverse(), computeInverse()
*/
template<typename Rhs>
- inline const ei_solve_return_value<PartialPivLU, Rhs>
+ inline const ei_solve_retval<PartialPivLU, Rhs>
solve(const MatrixBase<Rhs>& b) const
{
ei_assert(m_isInitialized && "PartialPivLU is not initialized.");
- return ei_solve_return_value<PartialPivLU, Rhs>(*this, b.derived());
+ return ei_solve_retval<PartialPivLU, Rhs>(*this, b.derived());
}
/** \returns the inverse of the matrix of which *this is the LU decomposition.
@@ -147,10 +147,10 @@ template<typename _MatrixType> class PartialPivLU
*
* \sa MatrixBase::inverse(), LU::inverse()
*/
- inline const ei_solve_return_value<PartialPivLU,NestByValue<typename MatrixType::IdentityReturnType> > inverse() const
+ inline const ei_solve_retval<PartialPivLU,NestByValue<typename MatrixType::IdentityReturnType> > inverse() const
{
ei_assert(m_isInitialized && "PartialPivLU is not initialized.");
- return ei_solve_return_value<PartialPivLU,NestByValue<typename MatrixType::IdentityReturnType> >
+ return ei_solve_retval<PartialPivLU,NestByValue<typename MatrixType::IdentityReturnType> >
(*this, MatrixType::Identity(m_lu.rows(), m_lu.cols()).nestByValue());
}
@@ -408,8 +408,8 @@ typename ei_traits<MatrixType>::Scalar PartialPivLU<MatrixType>::determinant() c
/***** Implementation of solve() *****************************************************/
template<typename _MatrixType, typename Rhs>
-struct ei_solve_impl<PartialPivLU<_MatrixType>, Rhs>
- : ei_solve_return_value<PartialPivLU<_MatrixType>, Rhs>
+struct ei_solve_retval<PartialPivLU<_MatrixType>, Rhs>
+ : ei_solve_retval_base<PartialPivLU<_MatrixType>, Rhs>
{
EIGEN_MAKE_SOLVE_HELPERS(PartialPivLU<_MatrixType>,Rhs)