fix a couple of compilations issues

This commit is contained in:
Gael Guennebaud 2009-08-06 14:10:02 +02:00
parent c2861dd41a
commit 1d4fea48b5
8 changed files with 21 additions and 22 deletions

View File

@ -205,8 +205,8 @@ void LDLT<MatrixType>::compute(const MatrixType& a)
continue;
}
RealScalar Djj = ei_real(m_matrix.coeff(j,j) - (m_matrix.row(j).start(j)
* m_matrix.col(j).start(j).conjugate()).coeff(0,0));
RealScalar Djj = ei_real(m_matrix.coeff(j,j) - m_matrix.row(j).start(j)
.dot(m_matrix.col(j).start(j)));
m_matrix.coeffRef(j,j) = Djj;
// Finish early if the matrix is not full rank.

View File

@ -149,7 +149,7 @@ class GeneralProduct<Lhs, Rhs, InnerProduct>
template<typename Dest> void addTo(Dest& dst, Scalar alpha) const
{
ei_assert(dst.rows()==1 && dst.cols()==1);
dst.coeffRef(0,0) += (m_lhs.cwise()*m_rhs).sum();
dst.coeffRef(0,0) += alpha * (m_lhs.cwise()*m_rhs).sum();
}
};

View File

@ -911,11 +911,7 @@ Transform<Scalar,Dim,Mode>::inverse(TransformTraits hint) const
}
// translation and remaining parts
res.template corner<Dim,1>(TopRight) = - res.template corner<Dim,Dim>(TopLeft) * translation();
if (int(Mode)!=AffineCompact)
{
res.template corner<1,Dim>(BottomLeft).setZero();
res.coeffRef(Dim,Dim) = Scalar(1);
}
makeAffine();
return res;
}
}

View File

@ -174,8 +174,10 @@ umeyama(const MatrixBase<Derived>& src, const MatrixBase<OtherDerived>& dst, boo
const Scalar c = 1/src_var * svd.singularValues().dot(S);
// Eq. (41)
// TODO: lazyness does not make much sense over here, right?
Rt.col(m).segment(0,m) = dst_mean - c*Rt.block(0,0,m,m)*src_mean;
// Note that we first assign dst_mean to the destination so that there no need
// for a temporary.
Rt.col(m).start(m) = dst_mean;
Rt.col(m).start(m) -= (c*Rt.corner(TopLeft,m,m)*src_mean).lazy();
if (with_scaling) Rt.block(0,0,m,m) *= c;

View File

@ -170,7 +170,7 @@ void HessenbergDecomposition<MatrixType>::_compute(MatrixType& matA, CoeffVector
// 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;
// first let's do A = H A
// first let's do A = H' A
matA.corner(BottomRight,n-i-1,n-i-1) -= ((ei_conj(h) * matA.col(i).end(n-i-1)) *
(matA.col(i).end(n-i-1).adjoint() * matA.corner(BottomRight,n-i-1,n-i-1))).lazy();
@ -213,14 +213,15 @@ HessenbergDecomposition<MatrixType>::matrixQ(void) const
{
int n = m_matrix.rows();
MatrixType matQ = MatrixType::Identity(n,n);
Matrix<Scalar,1,MatrixType::ColsAtCompileTime> temp(n);
for (int i = n-2; i>=0; i--)
{
Scalar tmp = m_matrix.coeff(i+1,i);
m_matrix.const_cast_derived().coeffRef(i+1,i) = 1;
matQ.corner(BottomRight,n-i-1,n-i-1) -=
((m_hCoeffs.coeff(i) * m_matrix.col(i).end(n-i-1)) *
(m_matrix.col(i).end(n-i-1).adjoint() * matQ.corner(BottomRight,n-i-1,n-i-1)).lazy()).lazy();
int rs = n-i-1;
temp.end(rs) = (m_matrix.col(i).end(rs).adjoint() * matQ.corner(BottomRight,rs,rs)).lazy();
matQ.corner(BottomRight,rs,rs) -= (m_hCoeffs.coeff(i) * m_matrix.col(i).end(rs) * temp.end(rs)).lazy();
m_matrix.const_cast_derived().coeffRef(i+1,i) = tmp;
}

View File

@ -103,15 +103,15 @@ template<typename MatrixType> void product_notemporary(const MatrixType& m)
VERIFY_EVALUATION_COUNT( m3.col(c0) = ((s1 * m1).adjoint().template selfadjointView<LowerTriangular>() * (-m2.row(c0)*s3).adjoint()).lazy(), 0);
VERIFY_EVALUATION_COUNT( m3.col(c0) -= ((s1 * m1).adjoint().template selfadjointView<UpperTriangular>() * (-m2.row(c0)*s3).adjoint()).lazy(), 0);
VERIFY_EVALUATION_COUNT( m3.block(r0,r0,r1,r1) += ((m1.block(r0,r0,r1,r1).template selfadjointView<UpperTriangular>() * (s1*m2.block(c0,r0,c1,r1)) )).lazy(), 0);
VERIFY_EVALUATION_COUNT( m3.block(r0,r0,r1,r1) = ((m1.block(r0,r0,r1,r1).template selfadjointView<UpperTriangular>() * m2.block(c0,r0,c1,r1) )).lazy(), 0);
VERIFY_EVALUATION_COUNT( m3.block(r0,c0,r1,c1) += ((m1.block(r0,r0,r1,r1).template selfadjointView<UpperTriangular>() * (s1*m2.block(r0,c0,r1,c1)) )).lazy(), 0);
VERIFY_EVALUATION_COUNT( m3.block(r0,c0,r1,c1) = ((m1.block(r0,r0,r1,r1).template selfadjointView<UpperTriangular>() * m2.block(r0,c0,r1,c1) )).lazy(), 0);
VERIFY_EVALUATION_COUNT( m3.template selfadjointView<LowerTriangular>().rankUpdate(m2.adjoint()), 0);
m3.resize(1,1);
VERIFY_EVALUATION_COUNT( m3 = ((m1.block(r0,r0,r1,r1).template selfadjointView<LowerTriangular>() * m2.block(c0,r0,c1,r1) )).lazy(), 0);
VERIFY_EVALUATION_COUNT( m3 = ((m1.block(r0,r0,r1,r1).template selfadjointView<LowerTriangular>() * m2.block(r0,c0,r1,c1) )).lazy(), 0);
m3.resize(1,1);
VERIFY_EVALUATION_COUNT( m3 = ((m1.block(r0,r0,r1,r1).template triangularView<UnitUpperTriangular>() * m2.block(c0,r0,c1,r1) )).lazy(), 0);
VERIFY_EVALUATION_COUNT( m3 = ((m1.block(r0,r0,r1,r1).template triangularView<UnitUpperTriangular>() * m2.block(r0,c0,r1,c1) )).lazy(), 0);
}
void test_product_notemporary()

View File

@ -86,7 +86,7 @@ template<typename MatrixType> void submatrices(const MatrixType& m)
//check row() and col()
VERIFY_IS_APPROX(m1.col(c1).transpose(), m1.transpose().row(c1));
VERIFY_IS_APPROX(m1.col(c1).dot(square.row(r1)), (square.lazy() * m1.conjugate())(r1,c1));
VERIFY_IS_APPROX(m1.col(c1).dot(square.row(r1)), (square.lazy() * m1.conjugate()).eval()(r1,c1));
//check operator(), both constant and non-constant, on row() and col()
m1.row(r1) += s1 * m1.row(r2);
m1.col(c1) += s1 * m1.col(c2);

View File

@ -72,13 +72,13 @@ template<typename MatrixType> void triangular(const MatrixType& m)
// test overloaded operator=
m1.setZero();
m1.template triangularView<Eigen::UpperTriangular>() = (m2.transpose() * m2).lazy();
m3 = m2.transpose() * m2;
m1.template triangularView<Eigen::UpperTriangular>() = m2.transpose() + m2;
m3 = m2.transpose() + m2;
VERIFY_IS_APPROX(m3.template triangularView<Eigen::LowerTriangular>().transpose().toDense(), m1);
// test overloaded operator=
m1.setZero();
m1.template triangularView<Eigen::LowerTriangular>() = (m2.transpose() * m2).lazy();
m1.template triangularView<Eigen::LowerTriangular>() = m2.transpose() + m2;
VERIFY_IS_APPROX(m3.template triangularView<Eigen::LowerTriangular>().toDense(), m1);
m1 = MatrixType::Random(rows, cols);