mirror of
https://gitlab.com/libeigen/eigen.git
synced 2025-04-16 14:49:39 +08:00
TriangularMatrix: extend to rectangular matrices
This commit is contained in:
parent
2275f98d7b
commit
a20a744adc
@ -26,22 +26,11 @@
|
||||
#ifndef EIGEN_TRIANGULARMATRIX_H
|
||||
#define EIGEN_TRIANGULARMATRIX_H
|
||||
|
||||
/** \nonstableyet
|
||||
/** \internal
|
||||
*
|
||||
* \class TriangularBase
|
||||
*
|
||||
* \brief Expression of a triangular matrix extracted from a given matrix
|
||||
*
|
||||
* \param MatrixType the type of the object in which we are taking the triangular part
|
||||
* \param Mode the kind of triangular matrix expression to construct. Can be UpperTriangular,
|
||||
* LowerTriangular, UpperSelfadjoint, or LowerSelfadjoint. This is in fact a bit field;
|
||||
* it must have either UpperBit or LowerBit, and additionnaly it may have either
|
||||
* TraingularBit or SelfadjointBit.
|
||||
*
|
||||
* This class represents an expression of the upper or lower triangular part of
|
||||
* a square matrix, possibly with a further assumption on the diagonal. It is the return type
|
||||
* of MatrixBase::part() and most of the time this is the only way it is used.
|
||||
*
|
||||
* \sa MatrixBase::part()
|
||||
* \brief Base class for triangular part in a matrix
|
||||
*/
|
||||
template<typename Derived> class TriangularBase : public AnyMatrixBase<Derived>
|
||||
{
|
||||
@ -115,19 +104,21 @@ template<typename Derived> class TriangularBase : public AnyMatrixBase<Derived>
|
||||
|
||||
};
|
||||
|
||||
|
||||
/** \class TriangularView
|
||||
* \nonstableyet
|
||||
*
|
||||
* \brief Expression of a triangular part of a dense matrix
|
||||
* \brief Base class for triangular part in a matrix
|
||||
*
|
||||
* \param MatrixType the type of the dense matrix storing the coefficients
|
||||
* \param MatrixType the type of the object in which we are taking the triangular part
|
||||
* \param Mode the kind of triangular matrix expression to construct. Can be UpperTriangular,
|
||||
* LowerTriangular, UpperSelfadjoint, or LowerSelfadjoint. This is in fact a bit field;
|
||||
* it must have either UpperBit or LowerBit, and additionnaly it may have either
|
||||
* TraingularBit or SelfadjointBit.
|
||||
*
|
||||
* This class is an expression of a triangular part of a matrix with given dense
|
||||
* storage of the coefficients. It is the return type of MatrixBase::triangularPart()
|
||||
* and most of the time this is the only way that it is used.
|
||||
* This class represents a triangular part of a matrix, not necessarily square. Strictly speaking, for rectangular
|
||||
* matrices one should speak ok "trapezoid" parts. This class is the return type
|
||||
* of MatrixBase::triangularView() and most of the time this is the only way it is used.
|
||||
*
|
||||
* \sa class TriangularBase, MatrixBase::triangularPart(), class DiagonalWrapper
|
||||
* \sa MatrixBase::triangularView()
|
||||
*/
|
||||
template<typename MatrixType, unsigned int _Mode>
|
||||
struct ei_traits<TriangularView<MatrixType, _Mode> > : ei_traits<MatrixType>
|
||||
@ -155,7 +146,7 @@ template<typename _MatrixType, unsigned int _Mode> class TriangularView
|
||||
typedef TriangularBase<TriangularView> Base;
|
||||
typedef typename ei_traits<TriangularView>::Scalar Scalar;
|
||||
typedef _MatrixType MatrixType;
|
||||
typedef typename MatrixType::PlainMatrixType PlainMatrixType;
|
||||
typedef typename MatrixType::PlainMatrixType DenseMatrixType;
|
||||
typedef typename MatrixType::Nested MatrixTypeNested;
|
||||
typedef typename ei_cleantype<MatrixTypeNested>::type _MatrixTypeNested;
|
||||
|
||||
@ -244,9 +235,9 @@ template<typename _MatrixType, unsigned int _Mode> class TriangularView
|
||||
inline const TriangularView<NestByValue<Transpose<MatrixType> >,TransposeMode> transpose() const
|
||||
{ return m_matrix.transpose().nestByValue(); }
|
||||
|
||||
PlainMatrixType toDense() const
|
||||
DenseMatrixType toDenseMatrix() const
|
||||
{
|
||||
PlainMatrixType res(rows(), cols());
|
||||
DenseMatrixType res(rows(), cols());
|
||||
res = *this;
|
||||
return res;
|
||||
}
|
||||
@ -351,6 +342,7 @@ struct ei_triangular_assignment_selector
|
||||
}
|
||||
}
|
||||
};
|
||||
|
||||
template<typename Derived1, typename Derived2, unsigned int Mode, bool ClearOpposite>
|
||||
struct ei_triangular_assignment_selector<Derived1, Derived2, Mode, 1, ClearOpposite>
|
||||
{
|
||||
@ -365,6 +357,7 @@ struct ei_triangular_assignment_selector<Derived1, Derived2, Mode, 1, ClearOppos
|
||||
dst.copyCoeff(0, 0, src);
|
||||
}
|
||||
};
|
||||
|
||||
// prevent buggy user code from causing an infinite recursion
|
||||
template<typename Derived1, typename Derived2, unsigned int Mode, bool ClearOpposite>
|
||||
struct ei_triangular_assignment_selector<Derived1, Derived2, Mode, 0, ClearOpposite>
|
||||
@ -379,14 +372,16 @@ struct ei_triangular_assignment_selector<Derived1, Derived2, UpperTriangular, Dy
|
||||
{
|
||||
for(int j = 0; j < dst.cols(); ++j)
|
||||
{
|
||||
for(int i = 0; i <= j; ++i)
|
||||
int maxi = std::min(j, dst.rows()-1);
|
||||
for(int i = 0; i <= maxi; ++i)
|
||||
dst.copyCoeff(i, j, src);
|
||||
if (ClearOpposite)
|
||||
for(int i = j+1; i < dst.rows(); ++i)
|
||||
for(int i = maxi+1; i < dst.rows(); ++i)
|
||||
dst.coeffRef(i, j) = 0;
|
||||
}
|
||||
}
|
||||
};
|
||||
|
||||
template<typename Derived1, typename Derived2, bool ClearOpposite>
|
||||
struct ei_triangular_assignment_selector<Derived1, Derived2, LowerTriangular, Dynamic, ClearOpposite>
|
||||
{
|
||||
@ -396,8 +391,9 @@ struct ei_triangular_assignment_selector<Derived1, Derived2, LowerTriangular, Dy
|
||||
{
|
||||
for(int i = j; i < dst.rows(); ++i)
|
||||
dst.copyCoeff(i, j, src);
|
||||
int maxi = std::min(j, dst.rows());
|
||||
if (ClearOpposite)
|
||||
for(int i = 0; i < j; ++i)
|
||||
for(int i = 0; i < maxi; ++i)
|
||||
dst.coeffRef(i, j) = 0;
|
||||
}
|
||||
}
|
||||
@ -410,14 +406,16 @@ struct ei_triangular_assignment_selector<Derived1, Derived2, StrictlyUpperTriang
|
||||
{
|
||||
for(int j = 0; j < dst.cols(); ++j)
|
||||
{
|
||||
for(int i = 0; i < j; ++i)
|
||||
int maxi = std::min(j, dst.rows());
|
||||
for(int i = 0; i < maxi; ++i)
|
||||
dst.copyCoeff(i, j, src);
|
||||
if (ClearOpposite)
|
||||
for(int i = j; i < dst.rows(); ++i)
|
||||
for(int i = maxi; i < dst.rows(); ++i)
|
||||
dst.coeffRef(i, j) = 0;
|
||||
}
|
||||
}
|
||||
};
|
||||
|
||||
template<typename Derived1, typename Derived2, bool ClearOpposite>
|
||||
struct ei_triangular_assignment_selector<Derived1, Derived2, StrictlyLowerTriangular, Dynamic, ClearOpposite>
|
||||
{
|
||||
@ -427,8 +425,9 @@ struct ei_triangular_assignment_selector<Derived1, Derived2, StrictlyLowerTriang
|
||||
{
|
||||
for(int i = j+1; i < dst.rows(); ++i)
|
||||
dst.copyCoeff(i, j, src);
|
||||
int maxi = std::min(j, dst.rows());
|
||||
if (ClearOpposite)
|
||||
for(int i = 0; i <= j; ++i)
|
||||
for(int i = 0; i <= maxi; ++i)
|
||||
dst.coeffRef(i, j) = 0;
|
||||
}
|
||||
}
|
||||
@ -441,11 +440,12 @@ struct ei_triangular_assignment_selector<Derived1, Derived2, UnitUpperTriangular
|
||||
{
|
||||
for(int j = 0; j < dst.cols(); ++j)
|
||||
{
|
||||
for(int i = 0; i < j; ++i)
|
||||
int maxi = std::min(j, dst.rows());
|
||||
for(int i = 0; i < maxi; ++i)
|
||||
dst.copyCoeff(i, j, src);
|
||||
if (ClearOpposite)
|
||||
{
|
||||
for(int i = j+1; i < dst.rows(); ++i)
|
||||
for(int i = maxi+1; i < dst.rows(); ++i)
|
||||
dst.coeffRef(i, j) = 0;
|
||||
dst.coeffRef(j, j) = 1;
|
||||
}
|
||||
@ -459,11 +459,12 @@ struct ei_triangular_assignment_selector<Derived1, Derived2, UnitLowerTriangular
|
||||
{
|
||||
for(int j = 0; j < dst.cols(); ++j)
|
||||
{
|
||||
for(int i = j+1; i < dst.rows(); ++i)
|
||||
int maxi = std::min(j, dst.rows());
|
||||
for(int i = maxi+1; i < dst.rows(); ++i)
|
||||
dst.copyCoeff(i, j, src);
|
||||
if (ClearOpposite)
|
||||
{
|
||||
for(int i = 0; i < j; ++i)
|
||||
for(int i = 0; i < maxi; ++i)
|
||||
dst.coeffRef(i, j) = 0;
|
||||
dst.coeffRef(j, j) = 1;
|
||||
}
|
||||
@ -514,7 +515,7 @@ TriangularView<MatrixType, Mode>::operator=(const TriangularBase<OtherDerived>&
|
||||
ei_assert(Mode == OtherDerived::Mode);
|
||||
if(ei_traits<OtherDerived>::Flags & EvalBeforeAssigningBit)
|
||||
{
|
||||
typename OtherDerived::PlainMatrixType other_evaluated(other.rows(), other.cols());
|
||||
typename OtherDerived::DenseMatrixType other_evaluated(other.rows(), other.cols());
|
||||
other_evaluated.template triangularView<Mode>().lazyAssign(other.derived());
|
||||
lazyAssign(other_evaluated);
|
||||
}
|
||||
@ -633,17 +634,20 @@ const TriangularView<Derived, Mode> MatrixBase<Derived>::triangularView() const
|
||||
template<typename Derived>
|
||||
bool MatrixBase<Derived>::isUpperTriangular(RealScalar prec) const
|
||||
{
|
||||
if(cols() != rows()) return false;
|
||||
RealScalar maxAbsOnUpperTriangularPart = static_cast<RealScalar>(-1);
|
||||
for(int j = 0; j < cols(); ++j)
|
||||
for(int i = 0; i <= j; ++i)
|
||||
{
|
||||
int maxi = std::min(j, rows()-1);
|
||||
for(int i = 0; i <= maxi; ++i)
|
||||
{
|
||||
RealScalar absValue = ei_abs(coeff(i,j));
|
||||
if(absValue > maxAbsOnUpperTriangularPart) maxAbsOnUpperTriangularPart = absValue;
|
||||
}
|
||||
for(int j = 0; j < cols()-1; ++j)
|
||||
}
|
||||
RealScalar threshold = maxAbsOnUpperTriangularPart * prec;
|
||||
for(int j = 0; j < cols(); ++j)
|
||||
for(int i = j+1; i < rows(); ++i)
|
||||
if(!ei_isMuchSmallerThan(coeff(i, j), maxAbsOnUpperTriangularPart, prec)) return false;
|
||||
if(ei_abs(coeff(i, j)) > threshold) return false;
|
||||
return true;
|
||||
}
|
||||
|
||||
@ -655,7 +659,6 @@ bool MatrixBase<Derived>::isUpperTriangular(RealScalar prec) const
|
||||
template<typename Derived>
|
||||
bool MatrixBase<Derived>::isLowerTriangular(RealScalar prec) const
|
||||
{
|
||||
if(cols() != rows()) return false;
|
||||
RealScalar maxAbsOnLowerTriangularPart = static_cast<RealScalar>(-1);
|
||||
for(int j = 0; j < cols(); ++j)
|
||||
for(int i = j; i < rows(); ++i)
|
||||
@ -663,9 +666,13 @@ bool MatrixBase<Derived>::isLowerTriangular(RealScalar prec) const
|
||||
RealScalar absValue = ei_abs(coeff(i,j));
|
||||
if(absValue > maxAbsOnLowerTriangularPart) maxAbsOnLowerTriangularPart = absValue;
|
||||
}
|
||||
RealScalar threshold = maxAbsOnLowerTriangularPart * prec;
|
||||
for(int j = 1; j < cols(); ++j)
|
||||
for(int i = 0; i < j; ++i)
|
||||
if(!ei_isMuchSmallerThan(coeff(i, j), maxAbsOnLowerTriangularPart, prec)) return false;
|
||||
{
|
||||
int maxi = std::min(j, rows()-1);
|
||||
for(int i = 0; i < maxi; ++i)
|
||||
if(ei_abs(coeff(i, j)) > threshold) return false;
|
||||
}
|
||||
return true;
|
||||
}
|
||||
|
||||
|
@ -24,7 +24,9 @@
|
||||
|
||||
#include "main.h"
|
||||
|
||||
template<typename MatrixType> void triangular(const MatrixType& m)
|
||||
|
||||
|
||||
template<typename MatrixType> void triangular_square(const MatrixType& m)
|
||||
{
|
||||
typedef typename MatrixType::Scalar Scalar;
|
||||
typedef typename NumTraits<Scalar>::Real RealScalar;
|
||||
@ -51,8 +53,8 @@ template<typename MatrixType> void triangular(const MatrixType& m)
|
||||
v2 = VectorType::Random(rows),
|
||||
vzero = VectorType::Zero(rows);
|
||||
|
||||
MatrixType m1up = m1.template triangularView<Eigen::UpperTriangular>();
|
||||
MatrixType m2up = m2.template triangularView<Eigen::UpperTriangular>();
|
||||
MatrixType m1up = m1.template triangularView<UpperTriangular>();
|
||||
MatrixType m2up = m2.template triangularView<UpperTriangular>();
|
||||
|
||||
if (rows*cols>1)
|
||||
{
|
||||
@ -66,20 +68,20 @@ template<typename MatrixType> void triangular(const MatrixType& m)
|
||||
// test overloaded operator+=
|
||||
r1.setZero();
|
||||
r2.setZero();
|
||||
r1.template triangularView<Eigen::UpperTriangular>() += m1;
|
||||
r1.template triangularView<UpperTriangular>() += m1;
|
||||
r2 += m1up;
|
||||
VERIFY_IS_APPROX(r1,r2);
|
||||
|
||||
// test overloaded operator=
|
||||
m1.setZero();
|
||||
m1.template triangularView<Eigen::UpperTriangular>() = m2.transpose() + m2;
|
||||
m1.template triangularView<UpperTriangular>() = m2.transpose() + m2;
|
||||
m3 = m2.transpose() + m2;
|
||||
VERIFY_IS_APPROX(m3.template triangularView<Eigen::LowerTriangular>().transpose().toDense(), m1);
|
||||
VERIFY_IS_APPROX(m3.template triangularView<LowerTriangular>().transpose().toDenseMatrix(), m1);
|
||||
|
||||
// test overloaded operator=
|
||||
m1.setZero();
|
||||
m1.template triangularView<Eigen::LowerTriangular>() = m2.transpose() + m2;
|
||||
VERIFY_IS_APPROX(m3.template triangularView<Eigen::LowerTriangular>().toDense(), m1);
|
||||
m1.template triangularView<LowerTriangular>() = m2.transpose() + m2;
|
||||
VERIFY_IS_APPROX(m3.template triangularView<LowerTriangular>().toDenseMatrix(), m1);
|
||||
|
||||
m1 = MatrixType::Random(rows, cols);
|
||||
for (int i=0; i<rows; ++i)
|
||||
@ -87,49 +89,143 @@ template<typename MatrixType> void triangular(const MatrixType& m)
|
||||
|
||||
Transpose<MatrixType> trm4(m4);
|
||||
// test back and forward subsitution with a vector as the rhs
|
||||
m3 = m1.template triangularView<Eigen::UpperTriangular>();
|
||||
VERIFY(v2.isApprox(m3.adjoint() * (m1.adjoint().template triangularView<Eigen::LowerTriangular>().solve(v2)), largerEps));
|
||||
m3 = m1.template triangularView<Eigen::LowerTriangular>();
|
||||
VERIFY(v2.isApprox(m3.transpose() * (m1.transpose().template triangularView<Eigen::UpperTriangular>().solve(v2)), largerEps));
|
||||
m3 = m1.template triangularView<Eigen::UpperTriangular>();
|
||||
VERIFY(v2.isApprox(m3 * (m1.template triangularView<Eigen::UpperTriangular>().solve(v2)), largerEps));
|
||||
m3 = m1.template triangularView<Eigen::LowerTriangular>();
|
||||
VERIFY(v2.isApprox(m3.conjugate() * (m1.conjugate().template triangularView<Eigen::LowerTriangular>().solve(v2)), largerEps));
|
||||
m3 = m1.template triangularView<UpperTriangular>();
|
||||
VERIFY(v2.isApprox(m3.adjoint() * (m1.adjoint().template triangularView<LowerTriangular>().solve(v2)), largerEps));
|
||||
m3 = m1.template triangularView<LowerTriangular>();
|
||||
VERIFY(v2.isApprox(m3.transpose() * (m1.transpose().template triangularView<UpperTriangular>().solve(v2)), largerEps));
|
||||
m3 = m1.template triangularView<UpperTriangular>();
|
||||
VERIFY(v2.isApprox(m3 * (m1.template triangularView<UpperTriangular>().solve(v2)), largerEps));
|
||||
m3 = m1.template triangularView<LowerTriangular>();
|
||||
VERIFY(v2.isApprox(m3.conjugate() * (m1.conjugate().template triangularView<LowerTriangular>().solve(v2)), largerEps));
|
||||
|
||||
// test back and forward subsitution with a matrix as the rhs
|
||||
m3 = m1.template triangularView<Eigen::UpperTriangular>();
|
||||
VERIFY(m2.isApprox(m3.adjoint() * (m1.adjoint().template triangularView<Eigen::LowerTriangular>().solve(m2)), largerEps));
|
||||
m3 = m1.template triangularView<Eigen::LowerTriangular>();
|
||||
VERIFY(m2.isApprox(m3.transpose() * (m1.transpose().template triangularView<Eigen::UpperTriangular>().solve(m2)), largerEps));
|
||||
m3 = m1.template triangularView<Eigen::UpperTriangular>();
|
||||
VERIFY(m2.isApprox(m3 * (m1.template triangularView<Eigen::UpperTriangular>().solve(m2)), largerEps));
|
||||
m3 = m1.template triangularView<Eigen::LowerTriangular>();
|
||||
VERIFY(m2.isApprox(m3.conjugate() * (m1.conjugate().template triangularView<Eigen::LowerTriangular>().solve(m2)), largerEps));
|
||||
m3 = m1.template triangularView<UpperTriangular>();
|
||||
VERIFY(m2.isApprox(m3.adjoint() * (m1.adjoint().template triangularView<LowerTriangular>().solve(m2)), largerEps));
|
||||
m3 = m1.template triangularView<LowerTriangular>();
|
||||
VERIFY(m2.isApprox(m3.transpose() * (m1.transpose().template triangularView<UpperTriangular>().solve(m2)), largerEps));
|
||||
m3 = m1.template triangularView<UpperTriangular>();
|
||||
VERIFY(m2.isApprox(m3 * (m1.template triangularView<UpperTriangular>().solve(m2)), largerEps));
|
||||
m3 = m1.template triangularView<LowerTriangular>();
|
||||
VERIFY(m2.isApprox(m3.conjugate() * (m1.conjugate().template triangularView<LowerTriangular>().solve(m2)), largerEps));
|
||||
|
||||
// check M * inv(L) using in place API
|
||||
m4 = m3;
|
||||
m3.transpose().template triangularView<Eigen::UpperTriangular>().solveInPlace(trm4);
|
||||
m3.transpose().template triangularView<UpperTriangular>().solveInPlace(trm4);
|
||||
VERIFY(m4.cwise().abs().isIdentity(test_precision<RealScalar>()));
|
||||
|
||||
// check M * inv(U) using in place API
|
||||
m3 = m1.template triangularView<Eigen::UpperTriangular>();
|
||||
m3 = m1.template triangularView<UpperTriangular>();
|
||||
m4 = m3;
|
||||
m3.transpose().template triangularView<Eigen::LowerTriangular>().solveInPlace(trm4);
|
||||
m3.transpose().template triangularView<LowerTriangular>().solveInPlace(trm4);
|
||||
VERIFY(m4.cwise().abs().isIdentity(test_precision<RealScalar>()));
|
||||
|
||||
// check solve with unit diagonal
|
||||
m3 = m1.template triangularView<Eigen::UnitUpperTriangular>();
|
||||
VERIFY(m2.isApprox(m3 * (m1.template triangularView<Eigen::UnitUpperTriangular>().solve(m2)), largerEps));
|
||||
m3 = m1.template triangularView<UnitUpperTriangular>();
|
||||
VERIFY(m2.isApprox(m3 * (m1.template triangularView<UnitUpperTriangular>().solve(m2)), largerEps));
|
||||
|
||||
// VERIFY(( m1.template triangularView<Eigen::UpperTriangular>()
|
||||
// * m2.template triangularView<Eigen::UpperTriangular>()).isUpperTriangular());
|
||||
// VERIFY(( m1.template triangularView<UpperTriangular>()
|
||||
// * m2.template triangularView<UpperTriangular>()).isUpperTriangular());
|
||||
|
||||
// test swap
|
||||
m1.setOnes();
|
||||
m2.setZero();
|
||||
m2.template triangularView<Eigen::UpperTriangular>().swap(m1);
|
||||
m2.template triangularView<UpperTriangular>().swap(m1);
|
||||
m3.setZero();
|
||||
m3.template triangularView<Eigen::UpperTriangular>().setOnes();
|
||||
m3.template triangularView<UpperTriangular>().setOnes();
|
||||
VERIFY_IS_APPROX(m2,m3);
|
||||
|
||||
}
|
||||
|
||||
|
||||
template<typename MatrixType> void triangular_rect(const MatrixType& m)
|
||||
{
|
||||
typedef typename MatrixType::Scalar Scalar;
|
||||
typedef typename NumTraits<Scalar>::Real RealScalar;
|
||||
typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
|
||||
|
||||
int rows = m.rows();
|
||||
int cols = m.cols();
|
||||
|
||||
MatrixType m1 = MatrixType::Random(rows, cols),
|
||||
m2 = MatrixType::Random(rows, cols),
|
||||
m3(rows, cols),
|
||||
m4(rows, cols),
|
||||
r1(rows, cols),
|
||||
r2(rows, cols),
|
||||
mzero = MatrixType::Zero(rows, cols),
|
||||
mones = MatrixType::Ones(rows, cols),
|
||||
identity = Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime>
|
||||
::Identity(rows, rows),
|
||||
square = Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime>
|
||||
::Random(rows, rows);
|
||||
VectorType v1 = VectorType::Random(rows),
|
||||
v2 = VectorType::Random(rows),
|
||||
vzero = VectorType::Zero(rows);
|
||||
|
||||
MatrixType m1up = m1.template triangularView<UpperTriangular>();
|
||||
MatrixType m2up = m2.template triangularView<UpperTriangular>();
|
||||
|
||||
if (rows*cols>1)
|
||||
{
|
||||
VERIFY(m1up.isUpperTriangular());
|
||||
VERIFY(m2up.transpose().isLowerTriangular());
|
||||
VERIFY(!m2.isLowerTriangular());
|
||||
}
|
||||
|
||||
// VERIFY_IS_APPROX(m1up.transpose() * m2, m1.upper().transpose().lower() * m2);
|
||||
|
||||
// test overloaded operator+=
|
||||
r1.setZero();
|
||||
r2.setZero();
|
||||
r1.template triangularView<UpperTriangular>() += m1;
|
||||
r2 += m1up;
|
||||
VERIFY_IS_APPROX(r1,r2);
|
||||
|
||||
// test overloaded operator=
|
||||
m1.setZero();
|
||||
m1.template triangularView<UpperTriangular>() = 3 * m2;
|
||||
m3 = 3 * m2;
|
||||
VERIFY_IS_APPROX(m3.template triangularView<UpperTriangular>().toDenseMatrix(), m1);
|
||||
|
||||
m1.setZero();
|
||||
m1.template triangularView<LowerTriangular>() = 3 * m2;
|
||||
VERIFY_IS_APPROX(m3.template triangularView<LowerTriangular>().toDenseMatrix(), m1);
|
||||
|
||||
m1.setZero();
|
||||
m1.template triangularView<StrictlyUpperTriangular>() = 3 * m2;
|
||||
VERIFY_IS_APPROX(m3.template triangularView<StrictlyUpperTriangular>().toDenseMatrix(), m1);
|
||||
|
||||
m1.setZero();
|
||||
m1.template triangularView<StrictlyLowerTriangular>() = 3 * m2;
|
||||
VERIFY_IS_APPROX(m3.template triangularView<StrictlyLowerTriangular>().toDenseMatrix(), m1);
|
||||
|
||||
m1.setRandom();
|
||||
m2 = m1.template triangularView<UpperTriangular>();
|
||||
VERIFY(m2.isUpperTriangular());
|
||||
VERIFY(!m2.isLowerTriangular());
|
||||
m2 = m1.template triangularView<StrictlyUpperTriangular>();
|
||||
VERIFY(m2.isUpperTriangular());
|
||||
VERIFY(m2.diagonal().isMuchSmallerThan(RealScalar(1)));
|
||||
m2 = m1.template triangularView<UnitUpperTriangular>();
|
||||
VERIFY(m2.isUpperTriangular());
|
||||
m2.diagonal().cwise() -= Scalar(1);
|
||||
VERIFY(m2.diagonal().isMuchSmallerThan(RealScalar(1)));
|
||||
m2 = m1.template triangularView<LowerTriangular>();
|
||||
VERIFY(m2.isLowerTriangular());
|
||||
VERIFY(!m2.isUpperTriangular());
|
||||
m2 = m1.template triangularView<StrictlyLowerTriangular>();
|
||||
VERIFY(m2.isLowerTriangular());
|
||||
VERIFY(m2.diagonal().isMuchSmallerThan(RealScalar(1)));
|
||||
m2 = m1.template triangularView<UnitLowerTriangular>();
|
||||
VERIFY(m2.isLowerTriangular());
|
||||
m2.diagonal().cwise() -= Scalar(1);
|
||||
VERIFY(m2.diagonal().isMuchSmallerThan(RealScalar(1)));
|
||||
// test swap
|
||||
m1.setOnes();
|
||||
m2.setZero();
|
||||
m2.template triangularView<UpperTriangular>().swap(m1);
|
||||
m3.setZero();
|
||||
m3.template triangularView<UpperTriangular>().setOnes();
|
||||
VERIFY_IS_APPROX(m2,m3);
|
||||
|
||||
}
|
||||
@ -137,12 +233,19 @@ template<typename MatrixType> void triangular(const MatrixType& m)
|
||||
void test_triangular()
|
||||
{
|
||||
for(int i = 0; i < g_repeat ; i++) {
|
||||
CALL_SUBTEST_1( triangular(Matrix<float, 1, 1>()) );
|
||||
CALL_SUBTEST_2( triangular(Matrix<float, 2, 2>()) );
|
||||
CALL_SUBTEST_3( triangular(Matrix3d()) );
|
||||
CALL_SUBTEST_4( triangular(MatrixXcf(4, 4)) );
|
||||
CALL_SUBTEST_5( triangular(Matrix<std::complex<float>,8, 8>()) );
|
||||
CALL_SUBTEST_6( triangular(MatrixXcd(17,17)) );
|
||||
CALL_SUBTEST_7( triangular(Matrix<float,Dynamic,Dynamic,RowMajor>(5, 5)) );
|
||||
CALL_SUBTEST_1( triangular_square(Matrix<float, 1, 1>()) );
|
||||
CALL_SUBTEST_2( triangular_square(Matrix<float, 2, 2>()) );
|
||||
CALL_SUBTEST_3( triangular_square(Matrix3d()) );
|
||||
CALL_SUBTEST_4( triangular_square(MatrixXcf(4, 4)) );
|
||||
CALL_SUBTEST_5( triangular_square(Matrix<std::complex<float>,8, 8>()) );
|
||||
CALL_SUBTEST_6( triangular_square(MatrixXcd(17,17)) );
|
||||
CALL_SUBTEST_7( triangular_square(Matrix<float,Dynamic,Dynamic,RowMajor>(5, 5)) );
|
||||
|
||||
CALL_SUBTEST_8( triangular_rect(Matrix<float, 4, 5>()) );
|
||||
CALL_SUBTEST_9( triangular_rect(Matrix<double, 6, 2>()) );
|
||||
CALL_SUBTEST_4( triangular_rect(MatrixXcf(4, 10)) );
|
||||
CALL_SUBTEST_6( triangular_rect(MatrixXcd(11, 3)) );
|
||||
CALL_SUBTEST_7( triangular_rect(Matrix<float,Dynamic,Dynamic,RowMajor>(7, 6)) );
|
||||
|
||||
}
|
||||
}
|
||||
|
Loading…
x
Reference in New Issue
Block a user