aboutsummaryrefslogtreecommitdiffhomepage
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2009-07-27 13:17:39 +0200
committerGravatar Gael Guennebaud <g.gael@free.fr>2009-07-27 13:17:39 +0200
commit0590c18555bd5d195e29ee6a131285cf0f80f9d1 (patch)
treea7c26663e1f1c7b72e666889cb1bf7d4c0df01b7
parentb5e40642898e3af17e67a7b78fa3b63191e11bb7 (diff)
various compilation and bug fixes in selfadjoint stuff
-rw-r--r--Eigen/src/Core/SelfAdjointView.h13
-rw-r--r--Eigen/src/Core/products/SelfadjointMatrixVector.h3
-rw-r--r--Eigen/src/Core/products/SelfadjointProduct.h2
-rw-r--r--Eigen/src/Core/products/SelfadjointRank2Update.h10
-rw-r--r--Eigen/src/QR/Tridiagonalization.h18
-rw-r--r--test/eigensolver_selfadjoint.cpp4
-rw-r--r--test/product_extra.cpp15
-rw-r--r--test/product_selfadjoint.cpp27
-rw-r--r--test/product_symm.cpp84
-rw-r--r--test/product_syrk.cpp12
10 files changed, 93 insertions, 95 deletions
diff --git a/Eigen/src/Core/SelfAdjointView.h b/Eigen/src/Core/SelfAdjointView.h
index 4787bcbd8..3faebfc5d 100644
--- a/Eigen/src/Core/SelfAdjointView.h
+++ b/Eigen/src/Core/SelfAdjointView.h
@@ -125,20 +125,24 @@ template<typename MatrixType, unsigned int UpLo> class SelfAdjointView
* The vectors \a u and \c v \b must be column vectors, however they can be
* a adjoint expression without any overhead. Only the meaningful triangular
* part of the matrix is updated, the rest is left unchanged.
+ *
+ * \sa rankUpdate(const MatrixBase<DerivedU>&, Scalar)
*/
template<typename DerivedU, typename DerivedV>
- SelfAdjointView& rank2update(const MatrixBase<DerivedU>& u, const MatrixBase<DerivedV>& v, Scalar alpha = Scalar(1));
+ SelfAdjointView& rankUpdate(const MatrixBase<DerivedU>& u, const MatrixBase<DerivedV>& v, Scalar alpha = Scalar(1));
/** Perform a symmetric rank K update of the selfadjoint matrix \c *this:
* \f$ this = this + \alpha ( u u^* ) \f$ where \a u is a vector or matrix.
- *
+ *
* \returns a reference to \c *this
*
* Note that to perform \f$ this = this + \alpha ( u^* u ) \f$ you can simply
* call this function with u.adjoint().
+ *
+ * \sa rankUpdate(const MatrixBase<DerivedU>&, const MatrixBase<DerivedV>&, Scalar)
*/
template<typename DerivedU>
- SelfAdjointView& rankKupdate(const MatrixBase<DerivedU>& u, Scalar alpha = Scalar(1));
+ SelfAdjointView& rankUpdate(const MatrixBase<DerivedU>& u, Scalar alpha = Scalar(1));
/////////// Cholesky module ///////////
@@ -231,13 +235,14 @@ struct ei_selfadjoint_product_returntype<Lhs,LhsMode,false,Rhs,0,true>
template<typename Dest> void evalTo(Dest& dst) const
{
- dst.resize(m_lhs.rows(), m_rhs.cols());
dst.setZero();
evalTo(dst,1);
}
template<typename Dest> void evalTo(Dest& dst, Scalar alpha) const
{
+ ei_assert(dst.rows()==m_lhs.rows() && dst.cols()==m_rhs.cols());
+
const ActualLhsType lhs = LhsBlasTraits::extract(m_lhs);
const ActualRhsType rhs = RhsBlasTraits::extract(m_rhs);
diff --git a/Eigen/src/Core/products/SelfadjointMatrixVector.h b/Eigen/src/Core/products/SelfadjointMatrixVector.h
index 8e0dc9526..279836f8a 100644
--- a/Eigen/src/Core/products/SelfadjointMatrixVector.h
+++ b/Eigen/src/Core/products/SelfadjointMatrixVector.h
@@ -63,9 +63,6 @@ static EIGEN_DONT_INLINE void ei_product_selfadjoint_vector(
rhs = r;
}
- for (int i=0;i<size;i++)
- res[i] = 0;
-
int bound = std::max(0,size-8) & 0xfffffffE;
if (FirstTriangular)
bound = size - bound;
diff --git a/Eigen/src/Core/products/SelfadjointProduct.h b/Eigen/src/Core/products/SelfadjointProduct.h
index d971720d6..f6420d10c 100644
--- a/Eigen/src/Core/products/SelfadjointProduct.h
+++ b/Eigen/src/Core/products/SelfadjointProduct.h
@@ -126,7 +126,7 @@ struct ei_selfadjoint_product<Scalar,MatStorageOrder, ColMajor, AAT, UpLo>
template<typename MatrixType, unsigned int UpLo>
template<typename DerivedU>
SelfAdjointView<MatrixType,UpLo>& SelfAdjointView<MatrixType,UpLo>
-::rankKupdate(const MatrixBase<DerivedU>& u, Scalar alpha)
+::rankUpdate(const MatrixBase<DerivedU>& u, Scalar alpha)
{
typedef ei_blas_traits<DerivedU> UBlasTraits;
typedef typename UBlasTraits::DirectLinearAccessType ActualUType;
diff --git a/Eigen/src/Core/products/SelfadjointRank2Update.h b/Eigen/src/Core/products/SelfadjointRank2Update.h
index 65a321fde..f7dcc003f 100644
--- a/Eigen/src/Core/products/SelfadjointRank2Update.h
+++ b/Eigen/src/Core/products/SelfadjointRank2Update.h
@@ -41,7 +41,7 @@ struct ei_selfadjoint_rank2_update_selector<Scalar,UType,VType,LowerTriangular>
// std::cerr << "lower \n" << u.transpose() << "\n" << v.transpose() << "\n\n";
for (int i=0; i<size; ++i)
{
-// std::cerr <<
+// std::cerr <<
Map<Matrix<Scalar,Dynamic,1> >(mat+stride*i+i, size-i) +=
(alpha * ei_conj(u.coeff(i))) * v.end(size-i)
+ (alpha * ei_conj(v.coeff(i))) * u.end(size-i);
@@ -70,13 +70,13 @@ template<bool Cond, typename T> struct ei_conj_expr_if
template<typename MatrixType, unsigned int UpLo>
template<typename DerivedU, typename DerivedV>
SelfAdjointView<MatrixType,UpLo>& SelfAdjointView<MatrixType,UpLo>
-::rank2update(const MatrixBase<DerivedU>& u, const MatrixBase<DerivedV>& v, Scalar alpha)
+::rankUpdate(const MatrixBase<DerivedU>& u, const MatrixBase<DerivedV>& v, Scalar alpha)
{
typedef ei_blas_traits<DerivedU> UBlasTraits;
typedef typename UBlasTraits::DirectLinearAccessType ActualUType;
typedef typename ei_cleantype<ActualUType>::type _ActualUType;
const ActualUType actualU = UBlasTraits::extract(u.derived());
-
+
typedef ei_blas_traits<DerivedV> VBlasTraits;
typedef typename VBlasTraits::DirectLinearAccessType ActualVType;
typedef typename ei_cleantype<ActualVType>::type _ActualVType;
@@ -87,8 +87,8 @@ SelfAdjointView<MatrixType,UpLo>& SelfAdjointView<MatrixType,UpLo>
enum { IsRowMajor = (ei_traits<MatrixType>::Flags&RowMajorBit)?1:0 };
ei_selfadjoint_rank2_update_selector<Scalar,
- typename ei_conj_expr_if<IsRowMajor ^ UBlasTraits::NeedToConjugate,_ActualUType>::ret,
- typename ei_conj_expr_if<IsRowMajor ^ VBlasTraits::NeedToConjugate,_ActualVType>::ret,
+ typename ei_cleantype<typename ei_conj_expr_if<IsRowMajor ^ UBlasTraits::NeedToConjugate,_ActualUType>::ret>::type,
+ typename ei_cleantype<typename ei_conj_expr_if<IsRowMajor ^ VBlasTraits::NeedToConjugate,_ActualVType>::ret>::type,
(IsRowMajor ? (UpLo==UpperTriangular ? LowerTriangular : UpperTriangular) : UpLo)>
::run(const_cast<Scalar*>(_expression().data()),_expression().stride(),actualU,actualV,actualAlpha);
diff --git a/Eigen/src/QR/Tridiagonalization.h b/Eigen/src/QR/Tridiagonalization.h
index cf737bd30..a003fd59a 100644
--- a/Eigen/src/QR/Tridiagonalization.h
+++ b/Eigen/src/QR/Tridiagonalization.h
@@ -224,21 +224,19 @@ void Tridiagonalization<MatrixType>::_compute(MatrixType& matA, CoeffVectorType&
// Apply similarity transformation to remaining columns,
// i.e., A = H' A H where H = I - h v v' and v = matA.col(i).end(n-i-1)
-
matA.col(i).coeffRef(i+1) = 1;
- // hCoeffs.end(n-i-1) = (matA.corner(BottomRight,n-i-1,n-i-1).template part<LowerTriangular|SelfAdjoint>()
- // * matA.col(i).end(n-i-1)).lazy();
- // TODO map the above code to the function call below:
- ei_product_selfadjoint_vector<Scalar,MatrixType::Flags&RowMajorBit,LowerTriangularBit>
- (n-i-1,matA.corner(BottomRight,n-i-1,n-i-1).data(), matA.stride(), matA.col(i).end(n-i-1).data(), const_cast<Scalar*>(hCoeffs.end(n-i-1).data()));
+// hCoeffs.end(n-i-1) = (matA.corner(BottomRight,n-i-1,n-i-1).template selfadjointView<LowerTriangular>()
+// * (h * matA.col(i).end(n-i-1)));
+
+ hCoeffs.end(n-i-1).setZero();
+ ei_product_selfadjoint_vector<Scalar,MatrixType::Flags&RowMajorBit,LowerTriangular,false,false>
+ (n-i-1,matA.corner(BottomRight,n-i-1,n-i-1).data(), matA.stride(), matA.col(i).end(n-i-1).data(), 1, const_cast<Scalar*>(hCoeffs.end(n-i-1).data()), h);
- hCoeffs.end(n-i-1) = hCoeffs.end(n-i-1)*h
- + (h*ei_conj(h)*Scalar(-0.5)*(matA.col(i).end(n-i-1).dot(hCoeffs.end(n-i-1)))) *
- matA.col(i).end(n-i-1);
+ hCoeffs.end(n-i-1) += (h*Scalar(-0.5)*(matA.col(i).end(n-i-1).dot(hCoeffs.end(n-i-1)))) * matA.col(i).end(n-i-1);
matA.corner(BottomRight, n-i-1, n-i-1).template selfadjointView<LowerTriangular>()
- .rank2update(matA.col(i).end(n-i-1), hCoeffs.end(n-i-1), -1);
+ .rankUpdate(matA.col(i).end(n-i-1), hCoeffs.end(n-i-1), -1);
// note: at that point matA(i+1,i+1) is the (i+1)-th element of the final diagonal
// note: the sequence of the beta values leads to the subdiagonal entries
diff --git a/test/eigensolver_selfadjoint.cpp b/test/eigensolver_selfadjoint.cpp
index 6b5092775..1f10b131b 100644
--- a/test/eigensolver_selfadjoint.cpp
+++ b/test/eigensolver_selfadjoint.cpp
@@ -119,8 +119,8 @@ void test_eigensolver_selfadjoint()
// very important to test a 3x3 matrix since we provide a special path for it
CALL_SUBTEST( selfadjointeigensolver(Matrix3f()) );
CALL_SUBTEST( selfadjointeigensolver(Matrix4d()) );
- CALL_SUBTEST( selfadjointeigensolver(MatrixXf(4,4)) );
- CALL_SUBTEST( selfadjointeigensolver(MatrixXcd(7,7)) );
+ CALL_SUBTEST( selfadjointeigensolver(MatrixXf(10,10)) );
+ CALL_SUBTEST( selfadjointeigensolver(MatrixXcd(17,17)) );
CALL_SUBTEST( selfadjointeigensolver(MatrixXd(19,19)) );
// some trivial but implementation-wise tricky cases
diff --git a/test/product_extra.cpp b/test/product_extra.cpp
index e750be65e..213dbced6 100644
--- a/test/product_extra.cpp
+++ b/test/product_extra.cpp
@@ -62,13 +62,14 @@ template<typename MatrixType> void product_extra(const MatrixType& m)
// all the expressions in this test should be compiled as a single matrix product
// TODO: add internal checks to verify that
- VERIFY_IS_APPROX(m1 * m2.adjoint(), m1 * m2.adjoint().eval());
- VERIFY_IS_APPROX(m1.adjoint() * square.adjoint(), m1.adjoint().eval() * square.adjoint().eval());
- VERIFY_IS_APPROX(m1.adjoint() * m2, m1.adjoint().eval() * m2);
- VERIFY_IS_APPROX( (s1 * m1.adjoint()) * m2, (s1 * m1.adjoint()).eval() * m2);
- VERIFY_IS_APPROX( (- m1.adjoint() * s1) * (s3 * m2), (- m1.adjoint() * s1).eval() * (s3 * m2).eval());
- VERIFY_IS_APPROX( (s2 * m1.adjoint() * s1) * m2, (s2 * m1.adjoint() * s1).eval() * m2);
- VERIFY_IS_APPROX( (-m1*s2) * s1*m2.adjoint(), (-m1*s2).eval() * (s1*m2.adjoint()).eval());
+ VERIFY_IS_APPROX(m3 = (m1 * m2.adjoint()).lazy(), m1 * m2.adjoint().eval());
+ VERIFY_IS_APPROX(m3 = (m1.adjoint() * square.adjoint()).lazy(), m1.adjoint().eval() * square.adjoint().eval());
+ VERIFY_IS_APPROX(m3 = (m1.adjoint() * m2).lazy(), m1.adjoint().eval() * m2);
+ VERIFY_IS_APPROX(m3 = ((s1 * m1.adjoint()) * m2).lazy(), (s1 * m1.adjoint()).eval() * m2);
+ VERIFY_IS_APPROX(m3 = ((- m1.adjoint() * s1) * (s3 * m2)).lazy(), (- m1.adjoint() * s1).eval() * (s3 * m2).eval());
+ VERIFY_IS_APPROX(m3 = ((s2 * m1.adjoint() * s1) * m2).lazy(), (s2 * m1.adjoint() * s1).eval() * m2);
+ VERIFY_IS_APPROX(m3 = ((-m1*s2) * s1*m2.adjoint()).lazy(), (-m1*s2).eval() * (s1*m2.adjoint()).eval());
+
// a very tricky case where a scale factor has to be automatically conjugated:
VERIFY_IS_APPROX( m1.adjoint() * (s1*m2).conjugate(), (m1.adjoint()).eval() * ((s1*m2).conjugate()).eval());
diff --git a/test/product_selfadjoint.cpp b/test/product_selfadjoint.cpp
index efa487ab1..e47358197 100644
--- a/test/product_selfadjoint.cpp
+++ b/test/product_selfadjoint.cpp
@@ -52,42 +52,23 @@ template<typename MatrixType> void product_selfadjoint(const MatrixType& m)
m1 = (m1.adjoint() + m1).eval();
- // lower
- m2 = m1.template triangularView<LowerTriangular>();
- VERIFY_IS_APPROX(v3 = (s1*m2).template selfadjointView<LowerTriangular>() * (s2*v1), (s1*m1) * (s2*v1));
- VERIFY_IS_APPROX(v3 = (s1*m2.conjugate()).template selfadjointView<LowerTriangular>() * (s2*v1), (s1*m1.conjugate()) * (s2*v1));
- VERIFY_IS_APPROX(v3 = (s1*m2).template selfadjointView<LowerTriangular>() * (s2*m4.col(1)), (s1*m1) * (s2*m4.col(1)));
-
- VERIFY_IS_APPROX(v3 = (s1*m2).template selfadjointView<LowerTriangular>() * (s2*v1.conjugate()), (s1*m1) * (s2*v1.conjugate()));
- VERIFY_IS_APPROX(v3 = (s1*m2.conjugate()).template selfadjointView<LowerTriangular>() * (s2*v1.conjugate()), (s1*m1.conjugate()) * (s2*v1.conjugate()));
-
- // upper
- m2 = m1.template triangularView<UpperTriangular>();
- VERIFY_IS_APPROX(v3 = (s1*m2).template selfadjointView<UpperTriangular>() * (s2*v1), (s1*m1) * (s2*v1));
- VERIFY_IS_APPROX(v3 = (s1*m2.conjugate()).template selfadjointView<UpperTriangular>() * (s2*v1), (s1*m1.conjugate()) * (s2*v1));
- VERIFY_IS_APPROX(v3 = (s1*m2.adjoint()).template selfadjointView<LowerTriangular>() * (s2*v1), (s1*m1.adjoint()) * (s2*v1));
- VERIFY_IS_APPROX(v3 = (s1*m2.transpose()).template selfadjointView<LowerTriangular>() * (s2*v1), (s1*m1.transpose()) * (s2*v1));
-
- VERIFY_IS_APPROX(v3 = (s1*m2).template selfadjointView<UpperTriangular>() * (s2*v1.conjugate()), (s1*m1) * (s2*v1.conjugate()));
- VERIFY_IS_APPROX(v3 = (s1*m2.conjugate()).template selfadjointView<UpperTriangular>() * (s2*v1.conjugate()), (s1*m1.conjugate()) * (s2*v1.conjugate()));
-
// rank2 update
m2 = m1.template triangularView<LowerTriangular>();
- m2.template selfadjointView<LowerTriangular>().rank2update(v1,v2);
+ m2.template selfadjointView<LowerTriangular>().rankUpdate(v1,v2);
VERIFY_IS_APPROX(m2, (m1 + v1 * v2.adjoint()+ v2 * v1.adjoint()).template triangularView<LowerTriangular>().toDense());
m2 = m1.template triangularView<UpperTriangular>();
- m2.template selfadjointView<UpperTriangular>().rank2update(-v1,s2*v2,s3);
+ m2.template selfadjointView<UpperTriangular>().rankUpdate(-v1,s2*v2,s3);
VERIFY_IS_APPROX(m2, (m1 + (-s2*s3) * (v1 * v2.adjoint()+ v2 * v1.adjoint())).template triangularView<UpperTriangular>().toDense());
m2 = m1.template triangularView<UpperTriangular>();
- m2.template selfadjointView<UpperTriangular>().rank2update(-r1.adjoint(),r2.adjoint()*s3,s1);
+ m2.template selfadjointView<UpperTriangular>().rankUpdate(-r1.adjoint(),r2.adjoint()*s3,s1);
VERIFY_IS_APPROX(m2, (m1 + (-s3*s1) * (r1.adjoint() * r2 + r2.adjoint() * r1)).template triangularView<UpperTriangular>().toDense());
if (rows>1)
{
m2 = m1.template triangularView<LowerTriangular>();
- m2.block(1,1,rows-1,cols-1).template selfadjointView<LowerTriangular>().rank2update(v1.end(rows-1),v2.start(cols-1));
+ m2.block(1,1,rows-1,cols-1).template selfadjointView<LowerTriangular>().rankUpdate(v1.end(rows-1),v2.start(cols-1));
m3 = m1;
m3.block(1,1,rows-1,cols-1) += v1.end(rows-1) * v2.start(cols-1).adjoint()+ v2.start(cols-1) * v1.end(rows-1).adjoint();
VERIFY_IS_APPROX(m2, m3.template triangularView<LowerTriangular>().toDense());
diff --git a/test/product_symm.cpp b/test/product_symm.cpp
index 3a0cd94d0..54bf91fb9 100644
--- a/test/product_symm.cpp
+++ b/test/product_symm.cpp
@@ -24,25 +24,43 @@
#include "main.h"
-template<typename MatrixType> void symm(const MatrixType& m)
+template<int OtherSize> struct symm_extra {
+ template<typename M1, typename M2, typename Scalar>
+ static void run(M1& m1, M1& m2, M2& rhs2, M2& rhs22, M2& rhs23, Scalar s1, Scalar s2)
+ {
+ m2 = m1.template triangularView<LowerTriangular>();
+ VERIFY_IS_APPROX(rhs22 = (rhs2) * (m2).template selfadjointView<LowerTriangular>(),
+ rhs23 = (rhs2) * (m1));
+ VERIFY_IS_APPROX(rhs22 = (s2*rhs2) * (s1*m2).template selfadjointView<LowerTriangular>(),
+ rhs23 = (s2*rhs2) * (s1*m1));
+ }
+};
+
+template<> struct symm_extra<1> {
+ template<typename M1, typename M2, typename Scalar>
+ static void run(M1& m1, M1& m2, M2& rhs2, M2& rhs22, M2& rhs23, Scalar s1, Scalar s2) {}
+};
+
+template<typename Scalar, int Size, int OtherSize> void symm(int size = Size, int othersize = OtherSize)
{
- typedef typename MatrixType::Scalar Scalar;
typedef typename NumTraits<Scalar>::Real RealScalar;
- typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, Dynamic> Rhs1;
- typedef Matrix<Scalar, Dynamic, MatrixType::RowsAtCompileTime> Rhs2;
- typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, Dynamic,RowMajor> Rhs3;
- int rows = m.rows();
- int cols = m.cols();
+ typedef Matrix<Scalar, Size, Size> MatrixType;
+ typedef Matrix<Scalar, Size, OtherSize> Rhs1;
+ typedef Matrix<Scalar, OtherSize, Size> Rhs2;
+ typedef Matrix<Scalar, Size, OtherSize,RowMajor> Rhs3;
+
+ int rows = size;
+ int cols = size;
MatrixType m1 = MatrixType::Random(rows, cols),
m2 = MatrixType::Random(rows, cols);
m1 = (m1+m1.adjoint()).eval();
- Rhs1 rhs1 = Rhs1::Random(cols, ei_random<int>(1,320)), rhs12, rhs13;
- Rhs2 rhs2 = Rhs2::Random(ei_random<int>(1,320), rows), rhs22, rhs23;
- Rhs3 rhs3 = Rhs3::Random(cols, ei_random<int>(1,320)), rhs32, rhs33;
+ Rhs1 rhs1 = Rhs1::Random(cols, othersize), rhs12(cols, othersize), rhs13(cols, othersize);
+ Rhs2 rhs2 = Rhs2::Random(othersize, rows), rhs22(othersize, rows), rhs23(othersize, rows);
+ Rhs3 rhs3 = Rhs3::Random(cols, othersize), rhs32(cols, othersize), rhs33(cols, othersize);
Scalar s1 = ei_random<Scalar>(),
s2 = ei_random<Scalar>();
@@ -51,46 +69,44 @@ template<typename MatrixType> void symm(const MatrixType& m)
VERIFY_IS_APPROX(rhs12 = (s1*m2).template selfadjointView<LowerTriangular>() * (s2*rhs1),
rhs13 = (s1*m1) * (s2*rhs1));
- m2 = m1.template triangularView<UpperTriangular>();
- VERIFY_IS_APPROX(rhs12 = (s1*m2).template selfadjointView<UpperTriangular>() * (s2*rhs1),
- rhs13 = (s1*m1) * (s2*rhs1));
+ m2 = m1.template triangularView<UpperTriangular>(); rhs12.setRandom(); rhs13 = rhs12;
+ VERIFY_IS_APPROX(rhs12 += (s1*m2).template selfadjointView<UpperTriangular>() * (s2*rhs1),
+ rhs13 += (s1*m1) * (s2*rhs1));
m2 = m1.template triangularView<LowerTriangular>();
- VERIFY_IS_APPROX(rhs22 = (s1*m2).template selfadjointView<LowerTriangular>() * (s2*rhs2.adjoint()),
- rhs23 = (s1*m1) * (s2*rhs2.adjoint()));
+ VERIFY_IS_APPROX(rhs12 = (s1*m2).template selfadjointView<LowerTriangular>() * (s2*rhs2.adjoint()),
+ rhs13 = (s1*m1) * (s2*rhs2.adjoint()));
m2 = m1.template triangularView<UpperTriangular>();
- VERIFY_IS_APPROX(rhs22 = (s1*m2).template selfadjointView<UpperTriangular>() * (s2*rhs2.adjoint()),
- rhs23 = (s1*m1) * (s2*rhs2.adjoint()));
+ VERIFY_IS_APPROX(rhs12 = (s1*m2).template selfadjointView<UpperTriangular>() * (s2*rhs2.adjoint()),
+ rhs13 = (s1*m1) * (s2*rhs2.adjoint()));
m2 = m1.template triangularView<UpperTriangular>();
- VERIFY_IS_APPROX(rhs22 = (s1*m2.adjoint()).template selfadjointView<LowerTriangular>() * (s2*rhs2.adjoint()),
- rhs23 = (s1*m1.adjoint()) * (s2*rhs2.adjoint()));
+ VERIFY_IS_APPROX(rhs12 = (s1*m2.adjoint()).template selfadjointView<LowerTriangular>() * (s2*rhs2.adjoint()),
+ rhs13 = (s1*m1.adjoint()) * (s2*rhs2.adjoint()));
// test row major = <...>
- m2 = m1.template triangularView<LowerTriangular>();
- VERIFY_IS_APPROX(rhs32 = (s1*m2).template selfadjointView<LowerTriangular>() * (s2*rhs3),
- rhs33 = (s1*m1) * (s2 * rhs3));
+ m2 = m1.template triangularView<LowerTriangular>(); rhs12.setRandom(); rhs13 = rhs12;
+ VERIFY_IS_APPROX(rhs12 -= (s1*m2).template selfadjointView<LowerTriangular>() * (s2*rhs3),
+ rhs13 -= (s1*m1) * (s2 * rhs3));
m2 = m1.template triangularView<UpperTriangular>();
- VERIFY_IS_APPROX(rhs32 = (s1*m2.adjoint()).template selfadjointView<LowerTriangular>() * (s2*rhs3).conjugate(),
- rhs33 = (s1*m1.adjoint()) * (s2*rhs3).conjugate());
+ VERIFY_IS_APPROX(rhs12 = (s1*m2.adjoint()).template selfadjointView<LowerTriangular>() * (s2*rhs3).conjugate(),
+ rhs13 = (s1*m1.adjoint()) * (s2*rhs3).conjugate());
// test matrix * selfadjoint
- m2 = m1.template triangularView<LowerTriangular>();
- VERIFY_IS_APPROX(rhs22 = (rhs2) * (m2).template selfadjointView<LowerTriangular>(),
- rhs23 = (rhs2) * (m1));
- VERIFY_IS_APPROX(rhs22 = (s2*rhs2) * (s1*m2).template selfadjointView<LowerTriangular>(),
- rhs23 = (s2*rhs2) * (s1*m1));
+ symm_extra<OtherSize>::run(m1,m2,rhs2,rhs22,rhs23,s1,s2);
+
}
+
void test_product_symm()
{
for(int i = 0; i < g_repeat ; i++)
{
- int s;
- s = ei_random<int>(10,320);
- CALL_SUBTEST( symm(MatrixXf(s, s)) );
- s = ei_random<int>(10,320);
- CALL_SUBTEST( symm(MatrixXcd(s, s)) );
+ CALL_SUBTEST(( symm<float,Dynamic,Dynamic>(ei_random<int>(10,320),ei_random<int>(10,320)) ));
+ CALL_SUBTEST(( symm<std::complex<double>,Dynamic,Dynamic>(ei_random<int>(10,320),ei_random<int>(10,320)) ));
+
+ CALL_SUBTEST(( symm<float,Dynamic,1>(ei_random<int>(10,320)) ));
+ CALL_SUBTEST(( symm<std::complex<double>,Dynamic,1>(ei_random<int>(10,320)) ));
}
}
diff --git a/test/product_syrk.cpp b/test/product_syrk.cpp
index ff8bf933d..be0606c96 100644
--- a/test/product_syrk.cpp
+++ b/test/product_syrk.cpp
@@ -46,27 +46,27 @@ template<typename MatrixType> void syrk(const MatrixType& m)
s2 = ei_random<Scalar>();
m2.setZero();
- VERIFY_IS_APPROX((m2.template selfadjointView<LowerTriangular>().rankKupdate(rhs2,s1)._expression()),
+ VERIFY_IS_APPROX((m2.template selfadjointView<LowerTriangular>().rankUpdate(rhs2,s1)._expression()),
((s1 * rhs2 * rhs2.adjoint()).eval().template triangularView<LowerTriangular>().toDense()));
m2.setZero();
- VERIFY_IS_APPROX(m2.template selfadjointView<UpperTriangular>().rankKupdate(rhs2,s1)._expression(),
+ VERIFY_IS_APPROX(m2.template selfadjointView<UpperTriangular>().rankUpdate(rhs2,s1)._expression(),
(s1 * rhs2 * rhs2.adjoint()).eval().template triangularView<UpperTriangular>().toDense());
m2.setZero();
- VERIFY_IS_APPROX(m2.template selfadjointView<LowerTriangular>().rankKupdate(rhs1.adjoint(),s1)._expression(),
+ VERIFY_IS_APPROX(m2.template selfadjointView<LowerTriangular>().rankUpdate(rhs1.adjoint(),s1)._expression(),
(s1 * rhs1.adjoint() * rhs1).eval().template triangularView<LowerTriangular>().toDense());
m2.setZero();
- VERIFY_IS_APPROX(m2.template selfadjointView<UpperTriangular>().rankKupdate(rhs1.adjoint(),s1)._expression(),
+ VERIFY_IS_APPROX(m2.template selfadjointView<UpperTriangular>().rankUpdate(rhs1.adjoint(),s1)._expression(),
(s1 * rhs1.adjoint() * rhs1).eval().template triangularView<UpperTriangular>().toDense());
m2.setZero();
- VERIFY_IS_APPROX(m2.template selfadjointView<LowerTriangular>().rankKupdate(rhs3.adjoint(),s1)._expression(),
+ VERIFY_IS_APPROX(m2.template selfadjointView<LowerTriangular>().rankUpdate(rhs3.adjoint(),s1)._expression(),
(s1 * rhs3.adjoint() * rhs3).eval().template triangularView<LowerTriangular>().toDense());
m2.setZero();
- VERIFY_IS_APPROX(m2.template selfadjointView<UpperTriangular>().rankKupdate(rhs3.adjoint(),s1)._expression(),
+ VERIFY_IS_APPROX(m2.template selfadjointView<UpperTriangular>().rankUpdate(rhs3.adjoint(),s1)._expression(),
(s1 * rhs3.adjoint() * rhs3).eval().template triangularView<UpperTriangular>().toDense());
}