aboutsummaryrefslogtreecommitdiffhomepage
path: root/Eigen/src
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
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')
-rw-r--r--Eigen/src/Cholesky/LDLT.h8
-rw-r--r--Eigen/src/Cholesky/LLT.h8
-rw-r--r--Eigen/src/Core/util/ForwardDeclarations.h12
-rw-r--r--Eigen/src/Geometry/arch/Geometry_SSE.h4
-rw-r--r--Eigen/src/LU/FullPivLU.h28
-rw-r--r--Eigen/src/LU/PartialPivLU.h12
-rw-r--r--Eigen/src/QR/ColPivHouseholderQR.h12
-rw-r--r--Eigen/src/QR/FullPivHouseholderQR.h12
-rw-r--r--Eigen/src/QR/HouseholderQR.h8
-rw-r--r--Eigen/src/SVD/SVD.h8
-rw-r--r--Eigen/src/misc/Image.h35
-rw-r--r--Eigen/src/misc/Kernel.h28
-rw-r--r--Eigen/src/misc/Solve.h30
13 files changed, 115 insertions, 90 deletions
diff --git a/Eigen/src/Cholesky/LDLT.h b/Eigen/src/Cholesky/LDLT.h
index a1faae49c..d0f292634 100644
--- a/Eigen/src/Cholesky/LDLT.h
+++ b/Eigen/src/Cholesky/LDLT.h
@@ -124,13 +124,13 @@ template<typename _MatrixType> class LDLT
* \sa solveInPlace(), MatrixBase::ldlt()
*/
template<typename Rhs>
- inline const ei_solve_return_value<LDLT, Rhs>
+ inline const ei_solve_retval<LDLT, Rhs>
solve(const MatrixBase<Rhs>& b) const
{
ei_assert(m_isInitialized && "LDLT is not initialized.");
ei_assert(m_matrix.rows()==b.rows()
&& "LDLT::solve(): invalid number of rows of the right hand side matrix b");
- return ei_solve_return_value<LDLT, Rhs>(*this, b.derived());
+ return ei_solve_retval<LDLT, Rhs>(*this, b.derived());
}
template<typename Derived>
@@ -265,8 +265,8 @@ LDLT<MatrixType>& LDLT<MatrixType>::compute(const MatrixType& a)
}
template<typename _MatrixType, typename Rhs>
-struct ei_solve_impl<LDLT<_MatrixType>, Rhs>
- : ei_solve_return_value<LDLT<_MatrixType>, Rhs>
+struct ei_solve_retval<LDLT<_MatrixType>, Rhs>
+ : ei_solve_retval_base<LDLT<_MatrixType>, Rhs>
{
EIGEN_MAKE_SOLVE_HELPERS(LDLT<_MatrixType>,Rhs)
diff --git a/Eigen/src/Cholesky/LLT.h b/Eigen/src/Cholesky/LLT.h
index 30c48578a..a1706e53e 100644
--- a/Eigen/src/Cholesky/LLT.h
+++ b/Eigen/src/Cholesky/LLT.h
@@ -109,13 +109,13 @@ template<typename _MatrixType, int _UpLo> class LLT
* \sa solveInPlace(), MatrixBase::llt()
*/
template<typename Rhs>
- inline const ei_solve_return_value<LLT, Rhs>
+ inline const ei_solve_retval<LLT, Rhs>
solve(const MatrixBase<Rhs>& b) const
{
ei_assert(m_isInitialized && "LLT is not initialized.");
ei_assert(m_matrix.rows()==b.rows()
&& "LLT::solve(): invalid number of rows of the right hand side matrix b");
- return ei_solve_return_value<LLT, Rhs>(*this, b.derived());
+ return ei_solve_retval<LLT, Rhs>(*this, b.derived());
}
template<typename Derived>
@@ -259,8 +259,8 @@ LLT<MatrixType,_UpLo>& LLT<MatrixType,_UpLo>::compute(const MatrixType& a)
}
template<typename _MatrixType, int UpLo, typename Rhs>
-struct ei_solve_impl<LLT<_MatrixType, UpLo>, Rhs>
- : ei_solve_return_value<LLT<_MatrixType, UpLo>, Rhs>
+struct ei_solve_retval<LLT<_MatrixType, UpLo>, Rhs>
+ : ei_solve_retval_base<LLT<_MatrixType, UpLo>, Rhs>
{
typedef LLT<_MatrixType,UpLo> LLTType;
EIGEN_MAKE_SOLVE_HELPERS(LLTType,Rhs)
diff --git a/Eigen/src/Core/util/ForwardDeclarations.h b/Eigen/src/Core/util/ForwardDeclarations.h
index 323848919..86df5395e 100644
--- a/Eigen/src/Core/util/ForwardDeclarations.h
+++ b/Eigen/src/Core/util/ForwardDeclarations.h
@@ -66,12 +66,12 @@ template<typename ExpressionType> class WithFormat;
template<typename MatrixType> struct CommaInitializer;
template<typename Derived> class ReturnByValue;
-template<typename DecompositionType, typename Rhs> struct ei_solve_return_value;
-template<typename DecompositionType, typename Rhs> struct ei_solve_impl;
-template<typename DecompositionType> struct ei_kernel_return_value;
-template<typename DecompositionType> struct ei_kernel_impl;
-template<typename DecompositionType> struct ei_image_return_value;
-template<typename DecompositionType> struct ei_image_impl;
+template<typename DecompositionType, typename Rhs> struct ei_solve_retval_base;
+template<typename DecompositionType, typename Rhs> struct ei_solve_retval;
+template<typename DecompositionType> struct ei_kernel_retval_base;
+template<typename DecompositionType> struct ei_kernel_retval;
+template<typename DecompositionType> struct ei_image_retval_base;
+template<typename DecompositionType> struct ei_image_retval;
template<typename _Scalar, int Rows=Dynamic, int Cols=Dynamic, int Supers=Dynamic, int Subs=Dynamic, int Options=0> class BandMatrix;
diff --git a/Eigen/src/Geometry/arch/Geometry_SSE.h b/Eigen/src/Geometry/arch/Geometry_SSE.h
index 1b8f6aead..a6ed10d82 100644
--- a/Eigen/src/Geometry/arch/Geometry_SSE.h
+++ b/Eigen/src/Geometry/arch/Geometry_SSE.h
@@ -32,8 +32,8 @@ template<class Derived, class OtherDerived> struct ei_quat_product<EiArch_SSE, D
{
const __m128 mask = _mm_castsi128_ps(_mm_setr_epi32(0,0,0,0x80000000));
Quaternion<float> res;
- __m128 a = _a.coeffs().packet<Aligned>(0);
- __m128 b = _b.coeffs().packet<Aligned>(0);
+ __m128 a = _a.coeffs().template packet<Aligned>(0);
+ __m128 b = _b.coeffs().template packet<Aligned>(0);
__m128 flip1 = _mm_xor_ps(_mm_mul_ps(ei_vec4f_swizzle1(a,1,2,0,2),
ei_vec4f_swizzle1(b,2,0,1,2)),mask);
__m128 flip2 = _mm_xor_ps(_mm_mul_ps(ei_vec4f_swizzle1(a,3,3,3,1),
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)
diff --git a/Eigen/src/QR/ColPivHouseholderQR.h b/Eigen/src/QR/ColPivHouseholderQR.h
index 4e2359e0b..d6be9a506 100644
--- a/Eigen/src/QR/ColPivHouseholderQR.h
+++ b/Eigen/src/QR/ColPivHouseholderQR.h
@@ -98,11 +98,11 @@ template<typename _MatrixType> class ColPivHouseholderQR
* Output: \verbinclude ColPivHouseholderQR_solve.out
*/
template<typename Rhs>
- inline const ei_solve_return_value<ColPivHouseholderQR, Rhs>
+ inline const ei_solve_retval<ColPivHouseholderQR, Rhs>
solve(const MatrixBase<Rhs>& b) const
{
ei_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
- return ei_solve_return_value<ColPivHouseholderQR, Rhs>(*this, b.derived());
+ return ei_solve_retval<ColPivHouseholderQR, Rhs>(*this, b.derived());
}
HouseholderSequenceType matrixQ(void) const;
@@ -215,11 +215,11 @@ template<typename _MatrixType> class ColPivHouseholderQR
* Use isInvertible() to first determine whether this matrix is invertible.
*/
inline const
- ei_solve_return_value<ColPivHouseholderQR, NestByValue<typename MatrixType::IdentityReturnType> >
+ ei_solve_retval<ColPivHouseholderQR, NestByValue<typename MatrixType::IdentityReturnType> >
inverse() const
{
ei_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
- return ei_solve_return_value<ColPivHouseholderQR,NestByValue<typename MatrixType::IdentityReturnType> >
+ return ei_solve_retval<ColPivHouseholderQR,NestByValue<typename MatrixType::IdentityReturnType> >
(*this, MatrixType::Identity(m_qr.rows(), m_qr.cols()).nestByValue());
}
@@ -325,8 +325,8 @@ ColPivHouseholderQR<MatrixType>& ColPivHouseholderQR<MatrixType>::compute(const
}
template<typename _MatrixType, typename Rhs>
-struct ei_solve_impl<ColPivHouseholderQR<_MatrixType>, Rhs>
- : ei_solve_return_value<ColPivHouseholderQR<_MatrixType>, Rhs>
+struct ei_solve_retval<ColPivHouseholderQR<_MatrixType>, Rhs>
+ : ei_solve_retval_base<ColPivHouseholderQR<_MatrixType>, Rhs>
{
EIGEN_MAKE_SOLVE_HELPERS(ColPivHouseholderQR<_MatrixType>,Rhs)
diff --git a/Eigen/src/QR/FullPivHouseholderQR.h b/Eigen/src/QR/FullPivHouseholderQR.h
index ba6866fc0..8be90a3c9 100644
--- a/Eigen/src/QR/FullPivHouseholderQR.h
+++ b/Eigen/src/QR/FullPivHouseholderQR.h
@@ -93,11 +93,11 @@ template<typename _MatrixType> class FullPivHouseholderQR
* Output: \verbinclude FullPivHouseholderQR_solve.out
*/
template<typename Rhs>
- inline const ei_solve_return_value<FullPivHouseholderQR, Rhs>
+ inline const ei_solve_retval<FullPivHouseholderQR, Rhs>
solve(const MatrixBase<Rhs>& b) const
{
ei_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
- return ei_solve_return_value<FullPivHouseholderQR, Rhs>(*this, b.derived());
+ return ei_solve_retval<FullPivHouseholderQR, Rhs>(*this, b.derived());
}
MatrixQType matrixQ(void) const;
@@ -215,11 +215,11 @@ template<typename _MatrixType> class FullPivHouseholderQR
* \note If this matrix is not invertible, the returned matrix has undefined coefficients.
* Use isInvertible() to first determine whether this matrix is invertible.
*/ inline const
- ei_solve_return_value<FullPivHouseholderQR, NestByValue<typename MatrixType::IdentityReturnType> >
+ ei_solve_retval<FullPivHouseholderQR, NestByValue<typename MatrixType::IdentityReturnType> >
inverse() const
{
ei_assert(m_isInitialized && "FullPivHouseholderQR is not initialized.");
- return ei_solve_return_value<FullPivHouseholderQR,NestByValue<typename MatrixType::IdentityReturnType> >
+ return ei_solve_retval<FullPivHouseholderQR,NestByValue<typename MatrixType::IdentityReturnType> >
(*this, MatrixType::Identity(m_qr.rows(), m_qr.cols()).nestByValue());
}
@@ -333,8 +333,8 @@ FullPivHouseholderQR<MatrixType>& FullPivHouseholderQR<MatrixType>::compute(cons
}
template<typename _MatrixType, typename Rhs>
-struct ei_solve_impl<FullPivHouseholderQR<_MatrixType>, Rhs>
- : ei_solve_return_value<FullPivHouseholderQR<_MatrixType>, Rhs>
+struct ei_solve_retval<FullPivHouseholderQR<_MatrixType>, Rhs>
+ : ei_solve_retval_base<FullPivHouseholderQR<_MatrixType>, Rhs>
{
EIGEN_MAKE_SOLVE_HELPERS(FullPivHouseholderQR<_MatrixType>,Rhs)
diff --git a/Eigen/src/QR/HouseholderQR.h b/Eigen/src/QR/HouseholderQR.h
index 8742fdf71..04b153cc4 100644
--- a/Eigen/src/QR/HouseholderQR.h
+++ b/Eigen/src/QR/HouseholderQR.h
@@ -98,11 +98,11 @@ template<typename _MatrixType> class HouseholderQR
* Output: \verbinclude HouseholderQR_solve.out
*/
template<typename Rhs>
- inline const ei_solve_return_value<HouseholderQR, Rhs>
+ inline const ei_solve_retval<HouseholderQR, Rhs>
solve(const MatrixBase<Rhs>& b) const
{
ei_assert(m_isInitialized && "HouseholderQR is not initialized.");
- return ei_solve_return_value<HouseholderQR, Rhs>(*this, b.derived());
+ return ei_solve_retval<HouseholderQR, Rhs>(*this, b.derived());
}
MatrixQType matrixQ() const;
@@ -210,8 +210,8 @@ HouseholderQR<MatrixType>& HouseholderQR<MatrixType>::compute(const MatrixType&
}
template<typename _MatrixType, typename Rhs>
-struct ei_solve_impl<HouseholderQR<_MatrixType>, Rhs>
- : ei_solve_return_value<HouseholderQR<_MatrixType>, Rhs>
+struct ei_solve_retval<HouseholderQR<_MatrixType>, Rhs>
+ : ei_solve_retval_base<HouseholderQR<_MatrixType>, Rhs>
{
EIGEN_MAKE_SOLVE_HELPERS(HouseholderQR<_MatrixType>,Rhs)
diff --git a/Eigen/src/SVD/SVD.h b/Eigen/src/SVD/SVD.h
index 2ae7c4859..3b424be04 100644
--- a/Eigen/src/SVD/SVD.h
+++ b/Eigen/src/SVD/SVD.h
@@ -90,11 +90,11 @@ template<typename _MatrixType> class SVD
* \sa MatrixBase::svd(),
*/
template<typename Rhs>
- inline const ei_solve_return_value<SVD, Rhs>
+ inline const ei_solve_retval<SVD, Rhs>
solve(const MatrixBase<Rhs>& b) const
{
ei_assert(m_isInitialized && "SVD is not initialized.");
- return ei_solve_return_value<SVD, Rhs>(*this, b.derived());
+ return ei_solve_retval<SVD, Rhs>(*this, b.derived());
}
const MatrixUType& matrixU() const
@@ -429,8 +429,8 @@ SVD<MatrixType>& SVD<MatrixType>::compute(const MatrixType& matrix)
}
template<typename _MatrixType, typename Rhs>
-struct ei_solve_impl<SVD<_MatrixType>, Rhs>
- : ei_solve_return_value<SVD<_MatrixType>, Rhs>
+struct ei_solve_retval<SVD<_MatrixType>, Rhs>
+ : ei_solve_retval_base<SVD<_MatrixType>, Rhs>
{
EIGEN_MAKE_SOLVE_HELPERS(SVD<_MatrixType>,Rhs)
diff --git a/Eigen/src/misc/Image.h b/Eigen/src/misc/Image.h
index 01e2b160d..9ed5d5f70 100644
--- a/Eigen/src/misc/Image.h
+++ b/Eigen/src/misc/Image.h
@@ -25,11 +25,11 @@
#ifndef EIGEN_MISC_IMAGE_H
#define EIGEN_MISC_IMAGE_H
-/** \class ei_image_return_value
+/** \class ei_image_retval_base
*
*/
template<typename DecompositionType>
-struct ei_traits<ei_image_return_value<DecompositionType> >
+struct ei_traits<ei_image_retval_base<DecompositionType> >
{
typedef typename DecompositionType::MatrixType MatrixType;
typedef Matrix<
@@ -43,17 +43,13 @@ struct ei_traits<ei_image_return_value<DecompositionType> >
> ReturnMatrixType;
};
-template<typename _DecompositionType> struct ei_image_return_value
- : public ReturnByValue<ei_image_return_value<_DecompositionType> >
+template<typename _DecompositionType> struct ei_image_retval_base
+ : public ReturnByValue<ei_image_retval_base<_DecompositionType> >
{
typedef _DecompositionType DecompositionType;
typedef typename DecompositionType::MatrixType MatrixType;
- const DecompositionType& m_dec;
- int m_rank, m_cols;
- const MatrixType& m_originalMatrix;
-
- ei_image_return_value(const DecompositionType& dec, const MatrixType& originalMatrix)
+ ei_image_retval_base(const DecompositionType& dec, const MatrixType& originalMatrix)
: m_dec(dec), m_rank(dec.rank()),
m_cols(m_rank == 0 ? 1 : m_rank),
m_originalMatrix(originalMatrix)
@@ -61,19 +57,32 @@ template<typename _DecompositionType> struct ei_image_return_value
inline int rows() const { return m_dec.rows(); }
inline int cols() const { return m_cols; }
+ inline int rank() const { return m_rank; }
+ inline const DecompositionType& dec() const { return m_dec; }
+ inline const MatrixType& originalMatrix() const { return m_originalMatrix; }
template<typename Dest> inline void evalTo(Dest& dst) const
{
- static_cast<const ei_image_impl<DecompositionType>*>(this)->evalTo(dst);
+ static_cast<const ei_image_retval<DecompositionType>*>(this)->evalTo(dst);
}
+
+ protected:
+ const DecompositionType& m_dec;
+ int m_rank, m_cols;
+ const MatrixType& m_originalMatrix;
};
#define EIGEN_MAKE_IMAGE_HELPERS(DecompositionType) \
typedef typename DecompositionType::MatrixType MatrixType; \
typedef typename MatrixType::Scalar Scalar; \
typedef typename MatrixType::RealScalar RealScalar; \
- inline const DecompositionType& dec() const { return this->m_dec; } \
- inline const MatrixType& originalMatrix() const { return this->m_originalMatrix; } \
- inline int rank() const { return this->m_rank; }
+ typedef ei_image_retval_base<DecompositionType> Base; \
+ using Base::dec; \
+ using Base::originalMatrix; \
+ using Base::rank; \
+ using Base::rows; \
+ using Base::cols; \
+ ei_image_retval(const DecompositionType& dec, const MatrixType& originalMatrix) \
+ : Base(dec, originalMatrix) {}
#endif // EIGEN_MISC_IMAGE_H
diff --git a/Eigen/src/misc/Kernel.h b/Eigen/src/misc/Kernel.h
index 74ef16abc..717eef450 100644
--- a/Eigen/src/misc/Kernel.h
+++ b/Eigen/src/misc/Kernel.h
@@ -25,11 +25,11 @@
#ifndef EIGEN_MISC_KERNEL_H
#define EIGEN_MISC_KERNEL_H
-/** \class ei_kernel_return_value
+/** \class ei_kernel_retval_base
*
*/
template<typename DecompositionType>
-struct ei_traits<ei_kernel_return_value<DecompositionType> >
+struct ei_traits<ei_kernel_retval_base<DecompositionType> >
{
typedef typename DecompositionType::MatrixType MatrixType;
typedef Matrix<
@@ -45,14 +45,12 @@ struct ei_traits<ei_kernel_return_value<DecompositionType> >
> ReturnMatrixType;
};
-template<typename _DecompositionType> struct ei_kernel_return_value
- : public ReturnByValue<ei_kernel_return_value<_DecompositionType> >
+template<typename _DecompositionType> struct ei_kernel_retval_base
+ : public ReturnByValue<ei_kernel_retval_base<_DecompositionType> >
{
typedef _DecompositionType DecompositionType;
- const DecompositionType& m_dec;
- int m_rank, m_cols;
- ei_kernel_return_value(const DecompositionType& dec)
+ ei_kernel_retval_base(const DecompositionType& dec)
: m_dec(dec),
m_rank(dec.rank()),
m_cols(m_rank==dec.cols() ? 1 : dec.cols() - m_rank)
@@ -60,18 +58,28 @@ template<typename _DecompositionType> struct ei_kernel_return_value
inline int rows() const { return m_dec.cols(); }
inline int cols() const { return m_cols; }
+ inline int rank() const { return m_rank; }
+ inline const DecompositionType& dec() const { return m_dec; }
template<typename Dest> inline void evalTo(Dest& dst) const
{
- static_cast<const ei_kernel_impl<DecompositionType>*>(this)->evalTo(dst);
+ static_cast<const ei_kernel_retval<DecompositionType>*>(this)->evalTo(dst);
}
+
+ protected:
+ const DecompositionType& m_dec;
+ int m_rank, m_cols;
};
#define EIGEN_MAKE_KERNEL_HELPERS(DecompositionType) \
typedef typename DecompositionType::MatrixType MatrixType; \
typedef typename MatrixType::Scalar Scalar; \
typedef typename MatrixType::RealScalar RealScalar; \
- inline const DecompositionType& dec() const { return this->m_dec; } \
- inline int rank() const { return this->m_rank; }
+ typedef ei_kernel_retval_base<DecompositionType> Base; \
+ using Base::dec; \
+ using Base::rank; \
+ using Base::rows; \
+ using Base::cols; \
+ ei_kernel_retval(const DecompositionType& dec) : Base(dec) {}
#endif // EIGEN_MISC_KERNEL_H
diff --git a/Eigen/src/misc/Solve.h b/Eigen/src/misc/Solve.h
index bf8f15773..d93869121 100644
--- a/Eigen/src/misc/Solve.h
+++ b/Eigen/src/misc/Solve.h
@@ -25,11 +25,11 @@
#ifndef EIGEN_MISC_SOLVE_H
#define EIGEN_MISC_SOLVE_H
-/** \class ei_solve_return_value
+/** \class ei_solve_retval_base
*
*/
template<typename DecompositionType, typename Rhs>
-struct ei_traits<ei_solve_return_value<DecompositionType, Rhs> >
+struct ei_traits<ei_solve_retval_base<DecompositionType, Rhs> >
{
typedef typename DecompositionType::MatrixType MatrixType;
typedef Matrix<typename Rhs::Scalar,
@@ -40,33 +40,41 @@ struct ei_traits<ei_solve_return_value<DecompositionType, Rhs> >
Rhs::MaxColsAtCompileTime> ReturnMatrixType;
};
-template<typename _DecompositionType, typename Rhs> struct ei_solve_return_value
- : public ReturnByValue<ei_solve_return_value<_DecompositionType, Rhs> >
+template<typename _DecompositionType, typename Rhs> struct ei_solve_retval_base
+ : public ReturnByValue<ei_solve_retval_base<_DecompositionType, Rhs> >
{
typedef typename ei_cleantype<typename Rhs::Nested>::type RhsNestedCleaned;
typedef _DecompositionType DecompositionType;
- const DecompositionType& m_dec;
- const typename Rhs::Nested m_rhs;
- ei_solve_return_value(const DecompositionType& dec, const Rhs& rhs)
+ ei_solve_retval_base(const DecompositionType& dec, const Rhs& rhs)
: m_dec(dec), m_rhs(rhs)
{}
inline int rows() const { return m_dec.cols(); }
inline int cols() const { return m_rhs.cols(); }
+ inline const DecompositionType& dec() const { return m_dec; }
+ inline const RhsNestedCleaned& rhs() const { return m_rhs; }
template<typename Dest> inline void evalTo(Dest& dst) const
{
- static_cast<const ei_solve_impl<DecompositionType,Rhs>*>(this)->evalTo(dst);
+ static_cast<const ei_solve_retval<DecompositionType,Rhs>*>(this)->evalTo(dst);
}
+
+ protected:
+ const DecompositionType& m_dec;
+ const typename Rhs::Nested m_rhs;
};
#define EIGEN_MAKE_SOLVE_HELPERS(DecompositionType,Rhs) \
typedef typename DecompositionType::MatrixType MatrixType; \
typedef typename MatrixType::Scalar Scalar; \
typedef typename MatrixType::RealScalar RealScalar; \
- typedef typename ei_cleantype<typename Rhs::Nested>::type RhsNestedCleaned; \
- inline const DecompositionType& dec() const { return this->m_dec; } \
- inline const RhsNestedCleaned& rhs() const { return this->m_rhs; }
+ typedef ei_solve_retval_base<DecompositionType,Rhs> Base; \
+ using Base::dec; \
+ using Base::rhs; \
+ using Base::rows; \
+ using Base::cols; \
+ ei_solve_retval(const DecompositionType& dec, const Rhs& rhs) \
+ : Base(dec, rhs) {}
#endif // EIGEN_MISC_SOLVE_H