2007-11-28 23:34:40 +08:00
|
|
|
// This file is part of Eigen, a lightweight C++ template library
|
2009-05-23 02:25:33 +08:00
|
|
|
// for linear algebra.
|
2007-11-28 23:34:40 +08:00
|
|
|
//
|
2008-11-24 21:40:43 +08:00
|
|
|
// Copyright (C) 2006-2008 Benoit Jacob <jacob.benoit.1@gmail.com>
|
2007-11-28 23:34:40 +08:00
|
|
|
//
|
2012-07-14 02:42:47 +08:00
|
|
|
// This Source Code Form is subject to the terms of the Mozilla
|
|
|
|
// Public License v. 2.0. If a copy of the MPL was not distributed
|
|
|
|
// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
|
2009-08-16 04:19:29 +08:00
|
|
|
|
2007-11-28 23:34:40 +08:00
|
|
|
#include "main.h"
|
|
|
|
|
2013-02-15 04:33:42 +08:00
|
|
|
template<bool IsInteger> struct adjoint_specific;
|
|
|
|
|
|
|
|
template<> struct adjoint_specific<true> {
|
|
|
|
template<typename Vec, typename Mat, typename Scalar>
|
|
|
|
static void run(const Vec& v1, const Vec& v2, Vec& v3, const Mat& square, Scalar s1, Scalar s2) {
|
2013-06-11 05:40:56 +08:00
|
|
|
VERIFY(test_isApproxWithRef((s1 * v1 + s2 * v2).dot(v3), numext::conj(s1) * v1.dot(v3) + numext::conj(s2) * v2.dot(v3), 0));
|
2013-02-15 04:33:42 +08:00
|
|
|
VERIFY(test_isApproxWithRef(v3.dot(s1 * v1 + s2 * v2), s1*v3.dot(v1)+s2*v3.dot(v2), 0));
|
|
|
|
|
|
|
|
// check compatibility of dot and adjoint
|
|
|
|
VERIFY(test_isApproxWithRef(v1.dot(square * v2), (square.adjoint() * v1).dot(v2), 0));
|
|
|
|
}
|
|
|
|
};
|
|
|
|
|
|
|
|
template<> struct adjoint_specific<false> {
|
|
|
|
template<typename Vec, typename Mat, typename Scalar>
|
|
|
|
static void run(const Vec& v1, const Vec& v2, Vec& v3, const Mat& square, Scalar s1, Scalar s2) {
|
|
|
|
typedef typename NumTraits<Scalar>::Real RealScalar;
|
2013-07-16 03:21:14 +08:00
|
|
|
using std::abs;
|
2013-02-15 04:33:42 +08:00
|
|
|
|
|
|
|
RealScalar ref = NumTraits<Scalar>::IsInteger ? RealScalar(0) : (std::max)((s1 * v1 + s2 * v2).norm(),v3.norm());
|
2013-06-11 05:40:56 +08:00
|
|
|
VERIFY(test_isApproxWithRef((s1 * v1 + s2 * v2).dot(v3), numext::conj(s1) * v1.dot(v3) + numext::conj(s2) * v2.dot(v3), ref));
|
2013-02-15 04:33:42 +08:00
|
|
|
VERIFY(test_isApproxWithRef(v3.dot(s1 * v1 + s2 * v2), s1*v3.dot(v1)+s2*v3.dot(v2), ref));
|
|
|
|
|
|
|
|
VERIFY_IS_APPROX(v1.squaredNorm(), v1.norm() * v1.norm());
|
|
|
|
// check normalized() and normalize()
|
|
|
|
VERIFY_IS_APPROX(v1, v1.norm() * v1.normalized());
|
|
|
|
v3 = v1;
|
|
|
|
v3.normalize();
|
|
|
|
VERIFY_IS_APPROX(v1, v1.norm() * v3);
|
|
|
|
VERIFY_IS_APPROX(v3, v1.normalized());
|
|
|
|
VERIFY_IS_APPROX(v3.norm(), RealScalar(1));
|
2016-01-22 03:43:42 +08:00
|
|
|
|
|
|
|
// check null inputs
|
|
|
|
VERIFY_IS_APPROX((v1*0).normalized(), (v1*0));
|
2016-01-31 05:14:04 +08:00
|
|
|
#if (!EIGEN_ARCH_i386) || defined(EIGEN_VECTORIZE)
|
2016-01-22 03:43:42 +08:00
|
|
|
RealScalar very_small = (std::numeric_limits<RealScalar>::min)();
|
|
|
|
VERIFY( (v1*very_small).norm() == 0 );
|
|
|
|
VERIFY_IS_APPROX((v1*very_small).normalized(), (v1*very_small));
|
|
|
|
v3 = v1*very_small;
|
|
|
|
v3.normalize();
|
|
|
|
VERIFY_IS_APPROX(v3, (v1*very_small));
|
2016-01-31 05:14:04 +08:00
|
|
|
#endif
|
2013-02-15 04:33:42 +08:00
|
|
|
|
|
|
|
// check compatibility of dot and adjoint
|
|
|
|
ref = NumTraits<Scalar>::IsInteger ? 0 : (std::max)((std::max)(v1.norm(),v2.norm()),(std::max)((square * v2).norm(),(square.adjoint() * v1).norm()));
|
2013-07-16 03:21:14 +08:00
|
|
|
VERIFY(internal::isMuchSmallerThan(abs(v1.dot(square * v2) - (square.adjoint() * v1).dot(v2)), ref, test_precision<Scalar>()));
|
2013-02-15 04:33:42 +08:00
|
|
|
|
|
|
|
// check that Random().normalized() works: tricky as the random xpr must be evaluated by
|
|
|
|
// normalized() in order to produce a consistent result.
|
|
|
|
VERIFY_IS_APPROX(Vec::Random(v1.size()).normalized().norm(), RealScalar(1));
|
|
|
|
}
|
|
|
|
};
|
|
|
|
|
2007-11-28 23:34:40 +08:00
|
|
|
template<typename MatrixType> void adjoint(const MatrixType& m)
|
|
|
|
{
|
|
|
|
/* this test covers the following files:
|
|
|
|
Transpose.h Conjugate.h Dot.h
|
|
|
|
*/
|
2012-11-06 22:25:50 +08:00
|
|
|
using std::abs;
|
2007-11-28 23:34:40 +08:00
|
|
|
typedef typename MatrixType::Scalar Scalar;
|
2008-08-23 23:14:20 +08:00
|
|
|
typedef typename NumTraits<Scalar>::Real RealScalar;
|
2008-03-13 01:17:36 +08:00
|
|
|
typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
|
2008-08-23 23:14:20 +08:00
|
|
|
typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> SquareMatrixType;
|
2015-01-27 00:09:01 +08:00
|
|
|
const Index PacketSize = internal::packet_traits<Scalar>::size;
|
2010-06-20 23:37:56 +08:00
|
|
|
|
|
|
|
Index rows = m.rows();
|
|
|
|
Index cols = m.cols();
|
2008-03-05 21:18:19 +08:00
|
|
|
|
2008-09-02 01:31:21 +08:00
|
|
|
MatrixType m1 = MatrixType::Random(rows, cols),
|
|
|
|
m2 = MatrixType::Random(rows, cols),
|
2007-11-28 23:34:40 +08:00
|
|
|
m3(rows, cols),
|
2008-09-02 01:31:21 +08:00
|
|
|
square = SquareMatrixType::Random(rows, rows);
|
|
|
|
VectorType v1 = VectorType::Random(rows),
|
|
|
|
v2 = VectorType::Random(rows),
|
|
|
|
v3 = VectorType::Random(rows),
|
2008-07-21 08:34:46 +08:00
|
|
|
vzero = VectorType::Zero(rows);
|
2007-11-28 23:34:40 +08:00
|
|
|
|
2010-10-25 22:15:22 +08:00
|
|
|
Scalar s1 = internal::random<Scalar>(),
|
|
|
|
s2 = internal::random<Scalar>();
|
2008-03-05 21:18:19 +08:00
|
|
|
|
2007-11-28 23:34:40 +08:00
|
|
|
// check basic compatibility of adjoint, transpose, conjugate
|
2007-12-03 18:23:08 +08:00
|
|
|
VERIFY_IS_APPROX(m1.transpose().conjugate().adjoint(), m1);
|
|
|
|
VERIFY_IS_APPROX(m1.adjoint().conjugate().transpose(), m1);
|
2008-03-05 21:18:19 +08:00
|
|
|
|
2007-11-28 23:34:40 +08:00
|
|
|
// check multiplicative behavior
|
2007-12-03 18:23:08 +08:00
|
|
|
VERIFY_IS_APPROX((m1.adjoint() * m2).adjoint(), m2.adjoint() * m1);
|
2013-06-11 05:40:56 +08:00
|
|
|
VERIFY_IS_APPROX((s1 * m1).adjoint(), numext::conj(s1) * m1.adjoint());
|
2008-03-05 21:18:19 +08:00
|
|
|
|
2013-02-15 04:33:42 +08:00
|
|
|
// check basic properties of dot, squaredNorm
|
2013-06-11 05:40:56 +08:00
|
|
|
VERIFY_IS_APPROX(numext::conj(v1.dot(v2)), v2.dot(v1));
|
|
|
|
VERIFY_IS_APPROX(numext::real(v1.dot(v1)), v1.squaredNorm());
|
2011-06-15 12:30:46 +08:00
|
|
|
|
2013-02-15 04:33:42 +08:00
|
|
|
adjoint_specific<NumTraits<Scalar>::IsInteger>::run(v1, v2, v3, square, s1, s2);
|
|
|
|
|
|
|
|
VERIFY_IS_MUCH_SMALLER_THAN(abs(vzero.dot(v1)), static_cast<RealScalar>(1));
|
2011-02-19 00:39:04 +08:00
|
|
|
|
2007-12-13 01:48:20 +08:00
|
|
|
// like in testBasicStuff, test operator() to check const-qualification
|
2010-10-25 22:15:22 +08:00
|
|
|
Index r = internal::random<Index>(0, rows-1),
|
|
|
|
c = internal::random<Index>(0, cols-1);
|
2013-06-11 05:40:56 +08:00
|
|
|
VERIFY_IS_APPROX(m1.conjugate()(r,c), numext::conj(m1(r,c)));
|
|
|
|
VERIFY_IS_APPROX(m1.adjoint()(c,r), numext::conj(m1(r,c)));
|
2008-03-05 21:18:19 +08:00
|
|
|
|
2008-10-29 23:24:08 +08:00
|
|
|
// check inplace transpose
|
|
|
|
m3 = m1;
|
|
|
|
m3.transposeInPlace();
|
|
|
|
VERIFY_IS_APPROX(m3,m1.transpose());
|
|
|
|
m3.transposeInPlace();
|
|
|
|
VERIFY_IS_APPROX(m3,m1);
|
2015-01-27 00:09:01 +08:00
|
|
|
|
|
|
|
if(PacketSize<m3.rows() && PacketSize<m3.cols())
|
|
|
|
{
|
|
|
|
m3 = m1;
|
|
|
|
Index i = internal::random<Index>(0,m3.rows()-PacketSize);
|
|
|
|
Index j = internal::random<Index>(0,m3.cols()-PacketSize);
|
|
|
|
m3.template block<PacketSize,PacketSize>(i,j).transposeInPlace();
|
|
|
|
VERIFY_IS_APPROX( (m3.template block<PacketSize,PacketSize>(i,j)), (m1.template block<PacketSize,PacketSize>(i,j).transpose()) );
|
|
|
|
m3.template block<PacketSize,PacketSize>(i,j).transposeInPlace();
|
|
|
|
VERIFY_IS_APPROX(m3,m1);
|
|
|
|
}
|
2009-03-31 21:55:40 +08:00
|
|
|
|
|
|
|
// check inplace adjoint
|
|
|
|
m3 = m1;
|
|
|
|
m3.adjointInPlace();
|
|
|
|
VERIFY_IS_APPROX(m3,m1.adjoint());
|
|
|
|
m3.transposeInPlace();
|
|
|
|
VERIFY_IS_APPROX(m3,m1.conjugate());
|
|
|
|
|
2011-01-27 16:59:19 +08:00
|
|
|
// check mixed dot product
|
|
|
|
typedef Matrix<RealScalar, MatrixType::RowsAtCompileTime, 1> RealVectorType;
|
|
|
|
RealVectorType rv1 = RealVectorType::Random(rows);
|
|
|
|
VERIFY_IS_APPROX(v1.dot(rv1.template cast<Scalar>()), v1.dot(rv1));
|
|
|
|
VERIFY_IS_APPROX(rv1.template cast<Scalar>().dot(v1), rv1.dot(v1));
|
2019-01-17 18:33:43 +08:00
|
|
|
|
|
|
|
VERIFY( is_same_type(m1,m1.template conjugateIf<false>()) );
|
|
|
|
VERIFY( is_same_type(m1.conjugate(),m1.template conjugateIf<true>()) );
|
2007-11-28 23:34:40 +08:00
|
|
|
}
|
|
|
|
|
2018-07-17 21:52:58 +08:00
|
|
|
template<int>
|
|
|
|
void adjoint_extra()
|
|
|
|
{
|
|
|
|
MatrixXcf a(10,10), b(10,10);
|
|
|
|
VERIFY_RAISES_ASSERT(a = a.transpose());
|
|
|
|
VERIFY_RAISES_ASSERT(a = a.transpose() + b);
|
|
|
|
VERIFY_RAISES_ASSERT(a = b + a.transpose());
|
|
|
|
VERIFY_RAISES_ASSERT(a = a.conjugate().transpose());
|
|
|
|
VERIFY_RAISES_ASSERT(a = a.adjoint());
|
|
|
|
VERIFY_RAISES_ASSERT(a = a.adjoint() + b);
|
|
|
|
VERIFY_RAISES_ASSERT(a = b + a.adjoint());
|
|
|
|
|
|
|
|
// no assertion should be triggered for these cases:
|
|
|
|
a.transpose() = a.transpose();
|
|
|
|
a.transpose() += a.transpose();
|
|
|
|
a.transpose() += a.transpose() + b;
|
|
|
|
a.transpose() = a.adjoint();
|
|
|
|
a.transpose() += a.adjoint();
|
|
|
|
a.transpose() += a.adjoint() + b;
|
|
|
|
|
|
|
|
// regression tests for check_for_aliasing
|
|
|
|
MatrixXd c(10,10);
|
|
|
|
c = 1.0 * MatrixXd::Ones(10,10) + c;
|
|
|
|
c = MatrixXd::Ones(10,10) * 1.0 + c;
|
|
|
|
c = c + MatrixXd::Ones(10,10) .cwiseProduct( MatrixXd::Zero(10,10) );
|
|
|
|
c = MatrixXd::Ones(10,10) * MatrixXd::Zero(10,10);
|
2019-01-16 21:33:45 +08:00
|
|
|
|
|
|
|
// regression for bug 1646
|
|
|
|
for (int j = 0; j < 10; ++j) {
|
|
|
|
c.col(j).head(j) = c.row(j).head(j);
|
|
|
|
}
|
|
|
|
|
2019-01-16 23:27:00 +08:00
|
|
|
for (int j = 0; j < 10; ++j) {
|
|
|
|
c.col(j) = c.row(j);
|
|
|
|
}
|
|
|
|
|
2019-01-16 21:33:45 +08:00
|
|
|
a.conservativeResize(1,1);
|
|
|
|
a = a.transpose();
|
|
|
|
|
|
|
|
a.conservativeResize(0,0);
|
|
|
|
a = a.transpose();
|
2018-07-17 21:52:58 +08:00
|
|
|
}
|
|
|
|
|
2018-07-17 20:46:15 +08:00
|
|
|
EIGEN_DECLARE_TEST(adjoint)
|
2007-11-28 23:34:40 +08:00
|
|
|
{
|
2009-07-15 05:00:53 +08:00
|
|
|
for(int i = 0; i < g_repeat; i++) {
|
2009-10-29 06:19:29 +08:00
|
|
|
CALL_SUBTEST_1( adjoint(Matrix<float, 1, 1>()) );
|
|
|
|
CALL_SUBTEST_2( adjoint(Matrix3d()) );
|
|
|
|
CALL_SUBTEST_3( adjoint(Matrix4f()) );
|
2015-01-27 00:09:01 +08:00
|
|
|
|
2011-07-12 20:41:00 +08:00
|
|
|
CALL_SUBTEST_4( adjoint(MatrixXcf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE/2), internal::random<int>(1,EIGEN_TEST_MAX_SIZE/2))) );
|
|
|
|
CALL_SUBTEST_5( adjoint(MatrixXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
|
|
|
|
CALL_SUBTEST_6( adjoint(MatrixXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
|
2015-01-27 00:09:01 +08:00
|
|
|
|
|
|
|
// Complement for 128 bits vectorization:
|
|
|
|
CALL_SUBTEST_8( adjoint(Matrix2d()) );
|
|
|
|
CALL_SUBTEST_9( adjoint(Matrix<int,4,4>()) );
|
|
|
|
|
|
|
|
// 256 bits vectorization:
|
|
|
|
CALL_SUBTEST_10( adjoint(Matrix<float,8,8>()) );
|
|
|
|
CALL_SUBTEST_11( adjoint(Matrix<double,4,4>()) );
|
|
|
|
CALL_SUBTEST_12( adjoint(Matrix<int,8,8>()) );
|
2009-07-14 03:14:47 +08:00
|
|
|
}
|
2011-07-12 20:41:00 +08:00
|
|
|
// test a large static matrix only once
|
2009-10-29 06:19:29 +08:00
|
|
|
CALL_SUBTEST_7( adjoint(Matrix<float, 100, 100>()) );
|
2009-09-07 18:46:16 +08:00
|
|
|
|
2018-07-17 21:52:58 +08:00
|
|
|
CALL_SUBTEST_13( adjoint_extra<0>() );
|
2007-11-28 23:34:40 +08:00
|
|
|
}
|
2007-12-03 02:32:59 +08:00
|
|
|
|