mirror of
https://gitlab.com/libeigen/eigen.git
synced 2025-06-04 18:54:00 +08:00
* make LU::kernel() and LU::image() also use ReturnByValue
* make them return zero vector in the degenerate case, instead of asserting (let's stick to the principle that we only assert on memory errors)
This commit is contained in:
parent
0ad3494bd3
commit
4f9e270343
@ -51,9 +51,9 @@ struct ei_nested<ReturnByValue<Derived>, n, PlainMatrixType>
|
|||||||
template<typename Derived>
|
template<typename Derived>
|
||||||
class ReturnByValue : public MatrixBase<ReturnByValue<Derived> >
|
class ReturnByValue : public MatrixBase<ReturnByValue<Derived> >
|
||||||
{
|
{
|
||||||
typedef typename ei_traits<Derived>::ReturnMatrixType ReturnMatrixType;
|
|
||||||
public:
|
public:
|
||||||
EIGEN_GENERIC_PUBLIC_INTERFACE(ReturnByValue)
|
EIGEN_GENERIC_PUBLIC_INTERFACE(ReturnByValue)
|
||||||
|
typedef typename ei_traits<Derived>::ReturnMatrixType ReturnMatrixType;
|
||||||
template<typename Dest>
|
template<typename Dest>
|
||||||
inline void evalTo(Dest& dst) const
|
inline void evalTo(Dest& dst) const
|
||||||
{ static_cast<const Derived* const>(this)->evalTo(dst); }
|
{ static_cast<const Derived* const>(this)->evalTo(dst); }
|
||||||
|
@ -26,88 +26,8 @@
|
|||||||
#define EIGEN_LU_H
|
#define EIGEN_LU_H
|
||||||
|
|
||||||
template<typename MatrixType, typename Rhs> struct ei_lu_solve_impl;
|
template<typename MatrixType, typename Rhs> struct ei_lu_solve_impl;
|
||||||
|
template<typename MatrixType> struct ei_lu_kernel_impl;
|
||||||
template<typename MatrixType,typename Rhs>
|
template<typename MatrixType> struct ei_lu_image_impl;
|
||||||
struct ei_traits<ei_lu_solve_impl<MatrixType,Rhs> >
|
|
||||||
{
|
|
||||||
typedef Matrix<typename Rhs::Scalar,
|
|
||||||
MatrixType::ColsAtCompileTime,
|
|
||||||
Rhs::ColsAtCompileTime,
|
|
||||||
Rhs::PlainMatrixType::Options,
|
|
||||||
MatrixType::MaxColsAtCompileTime,
|
|
||||||
Rhs::MaxColsAtCompileTime> ReturnMatrixType;
|
|
||||||
};
|
|
||||||
|
|
||||||
template<typename MatrixType, typename Rhs>
|
|
||||||
struct ei_lu_solve_impl : public ReturnByValue<ei_lu_solve_impl<MatrixType, Rhs> >
|
|
||||||
{
|
|
||||||
typedef typename ei_cleantype<typename Rhs::Nested>::type RhsNested;
|
|
||||||
typedef LU<MatrixType> LUType;
|
|
||||||
|
|
||||||
ei_lu_solve_impl(const LUType& lu, const Rhs& rhs)
|
|
||||||
: m_lu(lu), m_rhs(rhs)
|
|
||||||
{}
|
|
||||||
|
|
||||||
inline int rows() const { return m_lu.matrixLU().cols(); }
|
|
||||||
inline int cols() const { return m_rhs.cols(); }
|
|
||||||
|
|
||||||
template<typename Dest> void evalTo(Dest& dst) const
|
|
||||||
{
|
|
||||||
dst.resize(m_lu.matrixLU().cols(), m_rhs.cols());
|
|
||||||
if(m_lu.rank()==0)
|
|
||||||
{
|
|
||||||
dst.setZero();
|
|
||||||
return;
|
|
||||||
}
|
|
||||||
|
|
||||||
/* The decomposition PAQ = LU can be rewritten as A = P^{-1} L U Q^{-1}.
|
|
||||||
* So we proceed as follows:
|
|
||||||
* Step 1: compute c = P * rhs.
|
|
||||||
* Step 2: replace c by the solution x to Lx = c. Exists because L is invertible.
|
|
||||||
* Step 3: replace c by the solution x to Ux = c. May or may not exist.
|
|
||||||
* Step 4: result = Q * c;
|
|
||||||
*/
|
|
||||||
|
|
||||||
const int rows = m_lu.matrixLU().rows(),
|
|
||||||
cols = m_lu.matrixLU().cols(),
|
|
||||||
rank = m_lu.rank();
|
|
||||||
ei_assert(m_rhs.rows() == rows);
|
|
||||||
const int smalldim = std::min(rows, cols);
|
|
||||||
|
|
||||||
typename Rhs::PlainMatrixType c(m_rhs.rows(), m_rhs.cols());
|
|
||||||
|
|
||||||
// Step 1
|
|
||||||
for(int i = 0; i < rows; ++i)
|
|
||||||
c.row(m_lu.permutationP().coeff(i)) = m_rhs.row(i);
|
|
||||||
|
|
||||||
// Step 2
|
|
||||||
m_lu.matrixLU()
|
|
||||||
.corner(Eigen::TopLeft,smalldim,smalldim)
|
|
||||||
.template triangularView<UnitLowerTriangular>()
|
|
||||||
.solveInPlace(c.corner(Eigen::TopLeft, smalldim, c.cols()));
|
|
||||||
if(rows>cols)
|
|
||||||
{
|
|
||||||
c.corner(Eigen::BottomLeft, rows-cols, c.cols())
|
|
||||||
-= m_lu.matrixLU().corner(Eigen::BottomLeft, rows-cols, cols)
|
|
||||||
* c.corner(Eigen::TopLeft, cols, c.cols());
|
|
||||||
}
|
|
||||||
|
|
||||||
// Step 3
|
|
||||||
m_lu.matrixLU()
|
|
||||||
.corner(TopLeft, rank, rank)
|
|
||||||
.template triangularView<UpperTriangular>()
|
|
||||||
.solveInPlace(c.corner(TopLeft, rank, c.cols()));
|
|
||||||
|
|
||||||
// Step 4
|
|
||||||
for(int i = 0; i < rank; ++i)
|
|
||||||
dst.row(m_lu.permutationQ().coeff(i)) = c.row(i);
|
|
||||||
for(int i = rank; i < m_lu.matrixLU().cols(); ++i)
|
|
||||||
dst.row(m_lu.permutationQ().coeff(i)).setZero();
|
|
||||||
}
|
|
||||||
|
|
||||||
const LUType& m_lu;
|
|
||||||
const typename Rhs::Nested m_rhs;
|
|
||||||
};
|
|
||||||
|
|
||||||
/** \ingroup LU_Module
|
/** \ingroup LU_Module
|
||||||
*
|
*
|
||||||
@ -156,25 +76,6 @@ template<typename MatrixType> class LU
|
|||||||
MatrixType::MaxRowsAtCompileTime)
|
MatrixType::MaxRowsAtCompileTime)
|
||||||
};
|
};
|
||||||
|
|
||||||
typedef Matrix<typename MatrixType::Scalar,
|
|
||||||
MatrixType::ColsAtCompileTime, // the number of rows in the "kernel matrix" is the number of cols of the original matrix
|
|
||||||
// so that the product "matrix * kernel = zero" makes sense
|
|
||||||
Dynamic, // we don't know at compile-time the dimension of the kernel
|
|
||||||
MatrixType::Options,
|
|
||||||
MatrixType::MaxColsAtCompileTime, // see explanation for 2nd template parameter
|
|
||||||
MatrixType::MaxColsAtCompileTime // the kernel is a subspace of the domain space, whose dimension is the number
|
|
||||||
// of columns of the original matrix
|
|
||||||
> KernelResultType;
|
|
||||||
|
|
||||||
typedef Matrix<typename MatrixType::Scalar,
|
|
||||||
MatrixType::RowsAtCompileTime, // the image is a subspace of the destination space, whose
|
|
||||||
// dimension is the number of rows of the original matrix
|
|
||||||
Dynamic, // we don't know at compile time the dimension of the image (the rank)
|
|
||||||
MatrixType::Options,
|
|
||||||
MatrixType::MaxRowsAtCompileTime, // the image matrix will consist of columns from the original matrix,
|
|
||||||
MatrixType::MaxColsAtCompileTime // so it has the same number of rows and at most as many columns.
|
|
||||||
> ImageResultType;
|
|
||||||
|
|
||||||
/**
|
/**
|
||||||
* \brief Default Constructor.
|
* \brief Default Constructor.
|
||||||
*
|
*
|
||||||
@ -211,6 +112,14 @@ template<typename MatrixType> class LU
|
|||||||
return m_lu;
|
return m_lu;
|
||||||
}
|
}
|
||||||
|
|
||||||
|
/** \returns a pointer to the matrix of which *this is the LU decomposition.
|
||||||
|
*/
|
||||||
|
inline const MatrixType* originalMatrix() const
|
||||||
|
{
|
||||||
|
ei_assert(m_originalMatrix != 0 && "LU is not initialized.");
|
||||||
|
return m_originalMatrix;
|
||||||
|
}
|
||||||
|
|
||||||
/** \returns a vector of integers, whose size is the number of rows of the matrix being decomposed,
|
/** \returns a vector of integers, whose size is the number of rows of the matrix being decomposed,
|
||||||
* representing the P permutation i.e. the permutation of the rows. For its precise meaning,
|
* representing the P permutation i.e. the permutation of the rows. For its precise meaning,
|
||||||
* see the examples given in the documentation of class LU.
|
* see the examples given in the documentation of class LU.
|
||||||
@ -235,69 +144,37 @@ template<typename MatrixType> class LU
|
|||||||
return m_q;
|
return m_q;
|
||||||
}
|
}
|
||||||
|
|
||||||
/** Computes a basis of the kernel of the matrix, also called the null-space of the matrix.
|
|
||||||
*
|
|
||||||
* \note This method is only allowed on non-invertible matrices, as determined by
|
|
||||||
* isInvertible(). Calling it on an invertible matrix will make an assertion fail.
|
|
||||||
*
|
|
||||||
* \param result a pointer to the matrix in which to store the kernel. The columns of this
|
|
||||||
* matrix will be set to form a basis of the kernel (it will be resized
|
|
||||||
* if necessary).
|
|
||||||
*
|
|
||||||
* Example: \include LU_computeKernel.cpp
|
|
||||||
* Output: \verbinclude LU_computeKernel.out
|
|
||||||
*
|
|
||||||
* \sa kernel(), computeImage(), image()
|
|
||||||
*/
|
|
||||||
template<typename KernelMatrixType>
|
|
||||||
void computeKernel(KernelMatrixType *result) const;
|
|
||||||
|
|
||||||
/** Computes a basis of the image of the matrix, also called the column-space or range of he matrix.
|
|
||||||
*
|
|
||||||
* \note Calling this method on the zero matrix will make an assertion fail.
|
|
||||||
*
|
|
||||||
* \param result a pointer to the matrix in which to store the image. The columns of this
|
|
||||||
* matrix will be set to form a basis of the image (it will be resized
|
|
||||||
* if necessary).
|
|
||||||
*
|
|
||||||
* Example: \include LU_computeImage.cpp
|
|
||||||
* Output: \verbinclude LU_computeImage.out
|
|
||||||
*
|
|
||||||
* \sa image(), computeKernel(), kernel()
|
|
||||||
*/
|
|
||||||
template<typename ImageMatrixType>
|
|
||||||
void computeImage(ImageMatrixType *result) const;
|
|
||||||
|
|
||||||
/** \returns the kernel of the matrix, also called its null-space. The columns of the returned matrix
|
/** \returns the kernel of the matrix, also called its null-space. The columns of the returned matrix
|
||||||
* will form a basis of the kernel.
|
* will form a basis of the kernel.
|
||||||
*
|
*
|
||||||
* \note: this method is only allowed on non-invertible matrices, as determined by
|
* \note If the kernel has dimension zero, then the returned matrix is a column-vector filled with zeros.
|
||||||
* isInvertible(). Calling it on an invertible matrix will make an assertion fail.
|
|
||||||
*
|
|
||||||
* \note: this method returns a matrix by value, which induces some inefficiency.
|
|
||||||
* If you prefer to avoid this overhead, use computeKernel() instead.
|
|
||||||
*
|
*
|
||||||
* Example: \include LU_kernel.cpp
|
* Example: \include LU_kernel.cpp
|
||||||
* Output: \verbinclude LU_kernel.out
|
* Output: \verbinclude LU_kernel.out
|
||||||
*
|
*
|
||||||
* \sa computeKernel(), image()
|
* \sa image()
|
||||||
*/
|
*/
|
||||||
const KernelResultType kernel() const;
|
inline const ei_lu_kernel_impl<MatrixType> kernel() const
|
||||||
|
{
|
||||||
|
ei_assert(m_originalMatrix != 0 && "LU is not initialized.");
|
||||||
|
return ei_lu_kernel_impl<MatrixType>(*this);
|
||||||
|
}
|
||||||
|
|
||||||
/** \returns the image of the matrix, also called its column-space. The columns of the returned matrix
|
/** \returns the image of the matrix, also called its column-space. The columns of the returned matrix
|
||||||
* will form a basis of the kernel.
|
* will form a basis of the kernel.
|
||||||
*
|
*
|
||||||
* \note: Calling this method on the zero matrix will make an assertion fail.
|
* \note If the image has dimension zero, then the returned matrix is a column-vector filled with zeros.
|
||||||
*
|
|
||||||
* \note: this method returns a matrix by value, which induces some inefficiency.
|
|
||||||
* If you prefer to avoid this overhead, use computeImage() instead.
|
|
||||||
*
|
*
|
||||||
* Example: \include LU_image.cpp
|
* Example: \include LU_image.cpp
|
||||||
* Output: \verbinclude LU_image.out
|
* Output: \verbinclude LU_image.out
|
||||||
*
|
*
|
||||||
* \sa computeImage(), kernel()
|
* \sa kernel()
|
||||||
*/
|
*/
|
||||||
const ImageResultType image() const;
|
inline const ei_lu_image_impl<MatrixType> image() const
|
||||||
|
{
|
||||||
|
ei_assert(m_originalMatrix != 0 && "LU is not initialized.");
|
||||||
|
return ei_lu_image_impl<MatrixType>(*this);
|
||||||
|
}
|
||||||
|
|
||||||
/** This method returns a solution x to the equation Ax=b, where A is the matrix of which
|
/** This method returns a solution x to the equation Ax=b, where A is the matrix of which
|
||||||
* *this is the LU decomposition.
|
* *this is the LU decomposition.
|
||||||
@ -316,7 +193,7 @@ template<typename MatrixType> class LU
|
|||||||
* Example: \include LU_solve.cpp
|
* Example: \include LU_solve.cpp
|
||||||
* Output: \verbinclude LU_solve.out
|
* Output: \verbinclude LU_solve.out
|
||||||
*
|
*
|
||||||
* \sa TriangularView::solve(), kernel(), computeKernel(), inverse(), computeInverse()
|
* \sa TriangularView::solve(), kernel(), inverse(), computeInverse()
|
||||||
*/
|
*/
|
||||||
template<typename Rhs>
|
template<typename Rhs>
|
||||||
inline const ei_lu_solve_impl<MatrixType, Rhs>
|
inline const ei_lu_solve_impl<MatrixType, Rhs>
|
||||||
@ -547,13 +424,50 @@ typename ei_traits<MatrixType>::Scalar LU<MatrixType>::determinant() const
|
|||||||
return Scalar(m_det_pq) * m_lu.diagonal().prod();
|
return Scalar(m_det_pq) * m_lu.diagonal().prod();
|
||||||
}
|
}
|
||||||
|
|
||||||
|
/********* Implementation of kernel() **************************************************/
|
||||||
|
|
||||||
template<typename MatrixType>
|
template<typename MatrixType>
|
||||||
template<typename KernelMatrixType>
|
struct ei_traits<ei_lu_kernel_impl<MatrixType> >
|
||||||
void LU<MatrixType>::computeKernel(KernelMatrixType *result) const
|
|
||||||
{
|
{
|
||||||
ei_assert(m_originalMatrix != 0 && "LU is not initialized.");
|
typedef Matrix<
|
||||||
const int dimker = dimensionOfKernel(), cols = m_lu.cols();
|
typename MatrixType::Scalar,
|
||||||
result->resize(cols, dimker);
|
MatrixType::ColsAtCompileTime, // the number of rows in the "kernel matrix"
|
||||||
|
// is the number of cols of the original matrix
|
||||||
|
// so that the product "matrix * kernel = zero" makes sense
|
||||||
|
Dynamic, // we don't know at compile-time the dimension of the kernel
|
||||||
|
MatrixType::Options,
|
||||||
|
MatrixType::MaxColsAtCompileTime, // see explanation for 2nd template parameter
|
||||||
|
MatrixType::MaxColsAtCompileTime // the kernel is a subspace of the domain space,
|
||||||
|
// whose dimension is the number of columns of the original matrix
|
||||||
|
> ReturnMatrixType;
|
||||||
|
};
|
||||||
|
|
||||||
|
template<typename MatrixType>
|
||||||
|
struct ei_lu_kernel_impl : public ReturnByValue<ei_lu_kernel_impl<MatrixType> >
|
||||||
|
{
|
||||||
|
typedef LU<MatrixType> LUType;
|
||||||
|
const LUType& m_lu;
|
||||||
|
int m_dimKer;
|
||||||
|
|
||||||
|
ei_lu_kernel_impl(const LUType& lu) : m_lu(lu), m_dimKer(lu.dimensionOfKernel()) {}
|
||||||
|
|
||||||
|
inline int rows() const { return m_lu.matrixLU().cols(); }
|
||||||
|
inline int cols() const { return m_dimKer; }
|
||||||
|
|
||||||
|
template<typename Dest> void evalTo(Dest& dst) const
|
||||||
|
{
|
||||||
|
typedef typename MatrixType::Scalar Scalar;
|
||||||
|
const int rank = m_lu.rank(),
|
||||||
|
cols = m_lu.matrixLU().cols();
|
||||||
|
if(m_dimKer == 0)
|
||||||
|
{
|
||||||
|
// The Kernel is just {0}, so it doesn't have a basis properly speaking, but let's
|
||||||
|
// avoid crashing/asserting as that depends on floating point calculations. Let's
|
||||||
|
// just return a single column vector filled with zeros.
|
||||||
|
dst.resize(cols,1);
|
||||||
|
dst.setZero();
|
||||||
|
return;
|
||||||
|
}
|
||||||
|
|
||||||
/* Let us use the following lemma:
|
/* Let us use the following lemma:
|
||||||
*
|
*
|
||||||
@ -567,51 +481,156 @@ void LU<MatrixType>::computeKernel(KernelMatrixType *result) const
|
|||||||
*
|
*
|
||||||
* U is upper triangular, with eigenvalues sorted so that any zeros appear at the end.
|
* U is upper triangular, with eigenvalues sorted so that any zeros appear at the end.
|
||||||
* Thus, the diagonal of U ends with exactly
|
* Thus, the diagonal of U ends with exactly
|
||||||
* m_dimKer zero's. Let us use that to construct m_dimKer linearly
|
* m_dimKer zero's. Let us use that to construct dimKer linearly
|
||||||
* independent vectors in Ker U.
|
* independent vectors in Ker U.
|
||||||
*/
|
*/
|
||||||
|
|
||||||
Matrix<Scalar, Dynamic, Dynamic, MatrixType::Options,
|
dst.resize(cols, m_dimKer);
|
||||||
MatrixType::MaxColsAtCompileTime, MatrixType::MaxColsAtCompileTime>
|
|
||||||
y(-m_lu.corner(TopRight, m_rank, dimker));
|
|
||||||
|
|
||||||
m_lu.corner(TopLeft, m_rank, m_rank)
|
Matrix<typename MatrixType::Scalar, Dynamic, Dynamic, MatrixType::Options,
|
||||||
|
MatrixType::MaxColsAtCompileTime, MatrixType::MaxColsAtCompileTime>
|
||||||
|
y(-m_lu.matrixLU().corner(TopRight, rank, m_dimKer));
|
||||||
|
|
||||||
|
m_lu.matrixLU()
|
||||||
|
.corner(TopLeft, rank, rank)
|
||||||
.template triangularView<UpperTriangular>().solveInPlace(y);
|
.template triangularView<UpperTriangular>().solveInPlace(y);
|
||||||
|
|
||||||
for(int i = 0; i < m_rank; ++i) result->row(m_q.coeff(i)) = y.row(i);
|
for(int i = 0; i < rank; ++i) dst.row(m_lu.permutationQ().coeff(i)) = y.row(i);
|
||||||
for(int i = m_rank; i < cols; ++i) result->row(m_q.coeff(i)).setZero();
|
for(int i = rank; i < cols; ++i) dst.row(m_lu.permutationQ().coeff(i)).setZero();
|
||||||
for(int k = 0; k < dimker; ++k) result->coeffRef(m_q.coeff(m_rank+k), k) = Scalar(1);
|
for(int k = 0; k < m_dimKer; ++k) dst.coeffRef(m_lu.permutationQ().coeff(rank+k), k) = Scalar(1);
|
||||||
}
|
}
|
||||||
|
};
|
||||||
|
|
||||||
|
/***** Implementation of image() *****************************************************/
|
||||||
|
|
||||||
template<typename MatrixType>
|
template<typename MatrixType>
|
||||||
const typename LU<MatrixType>::KernelResultType
|
struct ei_traits<ei_lu_image_impl<MatrixType> >
|
||||||
LU<MatrixType>::kernel() const
|
|
||||||
{
|
{
|
||||||
ei_assert(m_originalMatrix != 0 && "LU is not initialized.");
|
typedef Matrix<
|
||||||
KernelResultType result(m_lu.cols(), dimensionOfKernel());
|
typename MatrixType::Scalar,
|
||||||
computeKernel(&result);
|
MatrixType::RowsAtCompileTime, // the image is a subspace of the destination space, whose
|
||||||
return result;
|
// dimension is the number of rows of the original matrix
|
||||||
}
|
Dynamic, // we don't know at compile time the dimension of the image (the rank)
|
||||||
|
MatrixType::Options,
|
||||||
|
MatrixType::MaxRowsAtCompileTime, // the image matrix will consist of columns from the original matrix,
|
||||||
|
MatrixType::MaxColsAtCompileTime // so it has the same number of rows and at most as many columns.
|
||||||
|
> ReturnMatrixType;
|
||||||
|
};
|
||||||
|
|
||||||
template<typename MatrixType>
|
template<typename MatrixType>
|
||||||
template<typename ImageMatrixType>
|
struct ei_lu_image_impl : public ReturnByValue<ei_lu_image_impl<MatrixType> >
|
||||||
void LU<MatrixType>::computeImage(ImageMatrixType *result) const
|
|
||||||
{
|
{
|
||||||
ei_assert(m_originalMatrix != 0 && "LU is not initialized.");
|
typedef LU<MatrixType> LUType;
|
||||||
result->resize(m_originalMatrix->rows(), m_rank);
|
const LUType& m_lu;
|
||||||
for(int i = 0; i < m_rank; ++i)
|
|
||||||
result->col(i) = m_originalMatrix->col(m_q.coeff(i));
|
ei_lu_image_impl(const LUType& lu) : m_lu(lu) {}
|
||||||
|
|
||||||
|
inline int rows() const { return m_lu.matrixLU().cols(); }
|
||||||
|
inline int cols() const { return m_lu.rank(); }
|
||||||
|
|
||||||
|
template<typename Dest> void evalTo(Dest& dst) const
|
||||||
|
{
|
||||||
|
int rank = m_lu.rank();
|
||||||
|
if(rank == 0)
|
||||||
|
{
|
||||||
|
// The Image is just {0}, so it doesn't have a basis properly speaking, but let's
|
||||||
|
// avoid crashing/asserting as that depends on floating point calculations. Let's
|
||||||
|
// just return a single column vector filled with zeros.
|
||||||
|
dst.resize(m_lu.originalMatrix()->rows(), 1);
|
||||||
|
dst.setZero();
|
||||||
|
return;
|
||||||
}
|
}
|
||||||
|
|
||||||
template<typename MatrixType>
|
dst.resize(m_lu.originalMatrix()->rows(), rank);
|
||||||
const typename LU<MatrixType>::ImageResultType
|
for(int i = 0; i < rank; ++i)
|
||||||
LU<MatrixType>::image() const
|
dst.col(i) = m_lu.originalMatrix()->col(m_lu.permutationQ().coeff(i));
|
||||||
{
|
|
||||||
ei_assert(m_originalMatrix != 0 && "LU is not initialized.");
|
|
||||||
ImageResultType result(m_originalMatrix->rows(), m_rank);
|
|
||||||
computeImage(&result);
|
|
||||||
return result;
|
|
||||||
}
|
}
|
||||||
|
};
|
||||||
|
|
||||||
|
/***** Implementation of solve() *****************************************************/
|
||||||
|
|
||||||
|
template<typename MatrixType,typename Rhs>
|
||||||
|
struct ei_traits<ei_lu_solve_impl<MatrixType,Rhs> >
|
||||||
|
{
|
||||||
|
typedef Matrix<typename Rhs::Scalar,
|
||||||
|
MatrixType::ColsAtCompileTime,
|
||||||
|
Rhs::ColsAtCompileTime,
|
||||||
|
Rhs::PlainMatrixType::Options,
|
||||||
|
MatrixType::MaxColsAtCompileTime,
|
||||||
|
Rhs::MaxColsAtCompileTime> ReturnMatrixType;
|
||||||
|
};
|
||||||
|
|
||||||
|
template<typename MatrixType, typename Rhs>
|
||||||
|
struct ei_lu_solve_impl : public ReturnByValue<ei_lu_solve_impl<MatrixType, Rhs> >
|
||||||
|
{
|
||||||
|
typedef typename ei_cleantype<typename Rhs::Nested>::type RhsNested;
|
||||||
|
typedef LU<MatrixType> LUType;
|
||||||
|
const LUType& m_lu;
|
||||||
|
const typename Rhs::Nested m_rhs;
|
||||||
|
|
||||||
|
ei_lu_solve_impl(const LUType& lu, const Rhs& rhs)
|
||||||
|
: m_lu(lu), m_rhs(rhs)
|
||||||
|
{}
|
||||||
|
|
||||||
|
inline int rows() const { return m_lu.matrixLU().cols(); }
|
||||||
|
inline int cols() const { return m_rhs.cols(); }
|
||||||
|
|
||||||
|
template<typename Dest> void evalTo(Dest& dst) const
|
||||||
|
{
|
||||||
|
dst.resize(m_lu.matrixLU().cols(), m_rhs.cols());
|
||||||
|
if(m_lu.rank()==0)
|
||||||
|
{
|
||||||
|
dst.setZero();
|
||||||
|
return;
|
||||||
|
}
|
||||||
|
|
||||||
|
/* The decomposition PAQ = LU can be rewritten as A = P^{-1} L U Q^{-1}.
|
||||||
|
* So we proceed as follows:
|
||||||
|
* Step 1: compute c = P * rhs.
|
||||||
|
* Step 2: replace c by the solution x to Lx = c. Exists because L is invertible.
|
||||||
|
* Step 3: replace c by the solution x to Ux = c. May or may not exist.
|
||||||
|
* Step 4: result = Q * c;
|
||||||
|
*/
|
||||||
|
|
||||||
|
const int rows = m_lu.matrixLU().rows(),
|
||||||
|
cols = m_lu.matrixLU().cols(),
|
||||||
|
rank = m_lu.rank();
|
||||||
|
ei_assert(m_rhs.rows() == rows);
|
||||||
|
const int smalldim = std::min(rows, cols);
|
||||||
|
|
||||||
|
typename Rhs::PlainMatrixType c(m_rhs.rows(), m_rhs.cols());
|
||||||
|
|
||||||
|
// Step 1
|
||||||
|
for(int i = 0; i < rows; ++i)
|
||||||
|
c.row(m_lu.permutationP().coeff(i)) = m_rhs.row(i);
|
||||||
|
|
||||||
|
// Step 2
|
||||||
|
m_lu.matrixLU()
|
||||||
|
.corner(Eigen::TopLeft,smalldim,smalldim)
|
||||||
|
.template triangularView<UnitLowerTriangular>()
|
||||||
|
.solveInPlace(c.corner(Eigen::TopLeft, smalldim, c.cols()));
|
||||||
|
if(rows>cols)
|
||||||
|
{
|
||||||
|
c.corner(Eigen::BottomLeft, rows-cols, c.cols())
|
||||||
|
-= m_lu.matrixLU().corner(Eigen::BottomLeft, rows-cols, cols)
|
||||||
|
* c.corner(Eigen::TopLeft, cols, c.cols());
|
||||||
|
}
|
||||||
|
|
||||||
|
// Step 3
|
||||||
|
m_lu.matrixLU()
|
||||||
|
.corner(TopLeft, rank, rank)
|
||||||
|
.template triangularView<UpperTriangular>()
|
||||||
|
.solveInPlace(c.corner(TopLeft, rank, c.cols()));
|
||||||
|
|
||||||
|
// Step 4
|
||||||
|
for(int i = 0; i < rank; ++i)
|
||||||
|
dst.row(m_lu.permutationQ().coeff(i)) = c.row(i);
|
||||||
|
for(int i = rank; i < m_lu.matrixLU().cols(); ++i)
|
||||||
|
dst.row(m_lu.permutationQ().coeff(i)).setZero();
|
||||||
|
}
|
||||||
|
};
|
||||||
|
|
||||||
|
/******* MatrixBase methods *****************************************************************/
|
||||||
|
|
||||||
/** \lu_module
|
/** \lu_module
|
||||||
*
|
*
|
||||||
|
@ -1,7 +1,7 @@
|
|||||||
// This file is part of Eigen, a lightweight C++ template library
|
// This file is part of Eigen, a lightweight C++ template library
|
||||||
// for linear algebra.
|
// for linear algebra.
|
||||||
//
|
//
|
||||||
// Copyright (C) 2008 Benoit Jacob <jacob.benoit.1@gmail.com>
|
// Copyright (C) 2008-2009 Benoit Jacob <jacob.benoit.1@gmail.com>
|
||||||
//
|
//
|
||||||
// Eigen is free software; you can redistribute it and/or
|
// Eigen is free software; you can redistribute it and/or
|
||||||
// modify it under the terms of the GNU Lesser General Public
|
// modify it under the terms of the GNU Lesser General Public
|
||||||
@ -37,9 +37,10 @@ template<typename MatrixType> void lu_non_invertible()
|
|||||||
createRandomMatrixOfRank(rank, rows, cols, m1);
|
createRandomMatrixOfRank(rank, rows, cols, m1);
|
||||||
|
|
||||||
LU<MatrixType> lu(m1);
|
LU<MatrixType> lu(m1);
|
||||||
typename LU<MatrixType>::KernelResultType m1kernel = lu.kernel();
|
typename ei_lu_kernel_impl<MatrixType>::ReturnMatrixType m1kernel = lu.kernel();
|
||||||
typename LU<MatrixType>::ImageResultType m1image = lu.image();
|
typename ei_lu_image_impl <MatrixType>::ReturnMatrixType m1image = lu.image();
|
||||||
|
|
||||||
|
std::cout << "rows:" << rows << " cols:" << cols << " | " << rank << " ----- " << lu.rank() << std::endl;
|
||||||
VERIFY(rank == lu.rank());
|
VERIFY(rank == lu.rank());
|
||||||
VERIFY(cols - lu.rank() == lu.dimensionOfKernel());
|
VERIFY(cols - lu.rank() == lu.dimensionOfKernel());
|
||||||
VERIFY(!lu.isInjective());
|
VERIFY(!lu.isInjective());
|
||||||
@ -101,8 +102,6 @@ template<typename MatrixType> void lu_verify_assert()
|
|||||||
VERIFY_RAISES_ASSERT(lu.matrixLU())
|
VERIFY_RAISES_ASSERT(lu.matrixLU())
|
||||||
VERIFY_RAISES_ASSERT(lu.permutationP())
|
VERIFY_RAISES_ASSERT(lu.permutationP())
|
||||||
VERIFY_RAISES_ASSERT(lu.permutationQ())
|
VERIFY_RAISES_ASSERT(lu.permutationQ())
|
||||||
VERIFY_RAISES_ASSERT(lu.computeKernel(&tmp))
|
|
||||||
VERIFY_RAISES_ASSERT(lu.computeImage(&tmp))
|
|
||||||
VERIFY_RAISES_ASSERT(lu.kernel())
|
VERIFY_RAISES_ASSERT(lu.kernel())
|
||||||
VERIFY_RAISES_ASSERT(lu.image())
|
VERIFY_RAISES_ASSERT(lu.image())
|
||||||
VERIFY_RAISES_ASSERT(lu.solve(tmp))
|
VERIFY_RAISES_ASSERT(lu.solve(tmp))
|
||||||
|
Loading…
x
Reference in New Issue
Block a user