blob: b63e843c66cf65b7cfd6ed01fe8777039e980a86 [file] [log] [blame]
Benoit Jacob2fdd0672007-11-28 15:34:40 +00001// This file is part of Eigen, a lightweight C++ template library
Benoit Jacob6347b1d2009-05-22 20:25:33 +02002// for linear algebra.
Benoit Jacob2fdd0672007-11-28 15:34:40 +00003//
Benoit Jacob00f89a82008-11-24 13:40:43 +00004// Copyright (C) 2006-2008 Benoit Jacob <jacob.benoit.1@gmail.com>
Benoit Jacob2fdd0672007-11-28 15:34:40 +00005//
Benoit Jacob69124cf2012-07-13 14:42:47 -04006// This Source Code Form is subject to the terms of the Mozilla
7// Public License v. 2.0. If a copy of the MPL was not distributed
8// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
Gael Guennebaud239ada92009-08-15 22:19:29 +02009
Benoit Jacobe2775862010-04-28 18:51:38 -040010#define EIGEN_NO_STATIC_ASSERT
11
Benoit Jacob2fdd0672007-11-28 15:34:40 +000012#include "main.h"
13
Gael Guennebauda0fb8852013-02-14 21:33:42 +010014template<bool IsInteger> struct adjoint_specific;
15
16template<> struct adjoint_specific<true> {
17 template<typename Vec, typename Mat, typename Scalar>
18 static void run(const Vec& v1, const Vec& v2, Vec& v3, const Mat& square, Scalar s1, Scalar s2) {
Gael Guennebaud62670c82013-06-10 23:40:56 +020019 VERIFY(test_isApproxWithRef((s1 * v1 + s2 * v2).dot(v3), numext::conj(s1) * v1.dot(v3) + numext::conj(s2) * v2.dot(v3), 0));
Gael Guennebauda0fb8852013-02-14 21:33:42 +010020 VERIFY(test_isApproxWithRef(v3.dot(s1 * v1 + s2 * v2), s1*v3.dot(v1)+s2*v3.dot(v2), 0));
21
22 // check compatibility of dot and adjoint
23 VERIFY(test_isApproxWithRef(v1.dot(square * v2), (square.adjoint() * v1).dot(v2), 0));
24 }
25};
26
27template<> struct adjoint_specific<false> {
28 template<typename Vec, typename Mat, typename Scalar>
29 static void run(const Vec& v1, const Vec& v2, Vec& v3, const Mat& square, Scalar s1, Scalar s2) {
30 typedef typename NumTraits<Scalar>::Real RealScalar;
31
32 RealScalar ref = NumTraits<Scalar>::IsInteger ? RealScalar(0) : (std::max)((s1 * v1 + s2 * v2).norm(),v3.norm());
Gael Guennebaud62670c82013-06-10 23:40:56 +020033 VERIFY(test_isApproxWithRef((s1 * v1 + s2 * v2).dot(v3), numext::conj(s1) * v1.dot(v3) + numext::conj(s2) * v2.dot(v3), ref));
Gael Guennebauda0fb8852013-02-14 21:33:42 +010034 VERIFY(test_isApproxWithRef(v3.dot(s1 * v1 + s2 * v2), s1*v3.dot(v1)+s2*v3.dot(v2), ref));
35
36 VERIFY_IS_APPROX(v1.squaredNorm(), v1.norm() * v1.norm());
37 // check normalized() and normalize()
38 VERIFY_IS_APPROX(v1, v1.norm() * v1.normalized());
39 v3 = v1;
40 v3.normalize();
41 VERIFY_IS_APPROX(v1, v1.norm() * v3);
42 VERIFY_IS_APPROX(v3, v1.normalized());
43 VERIFY_IS_APPROX(v3.norm(), RealScalar(1));
44
45 // check compatibility of dot and adjoint
46 ref = NumTraits<Scalar>::IsInteger ? 0 : (std::max)((std::max)(v1.norm(),v2.norm()),(std::max)((square * v2).norm(),(square.adjoint() * v1).norm()));
47 VERIFY(test_isApproxWithRef(v1.dot(square * v2), (square.adjoint() * v1).dot(v2), ref));
48
49 // check that Random().normalized() works: tricky as the random xpr must be evaluated by
50 // normalized() in order to produce a consistent result.
51 VERIFY_IS_APPROX(Vec::Random(v1.size()).normalized().norm(), RealScalar(1));
52 }
53};
54
Benoit Jacob2fdd0672007-11-28 15:34:40 +000055template<typename MatrixType> void adjoint(const MatrixType& m)
56{
57 /* this test covers the following files:
58 Transpose.h Conjugate.h Dot.h
59 */
Gael Guennebauda76fbbf2012-11-06 15:25:50 +010060 using std::abs;
Hauke Heibelf1679c72010-06-20 17:37:56 +020061 typedef typename MatrixType::Index Index;
Benoit Jacob2fdd0672007-11-28 15:34:40 +000062 typedef typename MatrixType::Scalar Scalar;
Gael Guennebaud2120fed2008-08-23 15:14:20 +000063 typedef typename NumTraits<Scalar>::Real RealScalar;
Benoit Jacob2ee68a02008-03-12 17:17:36 +000064 typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
Gael Guennebaud2120fed2008-08-23 15:14:20 +000065 typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> SquareMatrixType;
Hauke Heibelf1679c72010-06-20 17:37:56 +020066
67 Index rows = m.rows();
68 Index cols = m.cols();
Gael Guennebaud8e0d5482008-03-05 13:18:19 +000069
Benoit Jacob46fe7a32008-09-01 17:31:21 +000070 MatrixType m1 = MatrixType::Random(rows, cols),
71 m2 = MatrixType::Random(rows, cols),
Benoit Jacob2fdd0672007-11-28 15:34:40 +000072 m3(rows, cols),
Benoit Jacob46fe7a32008-09-01 17:31:21 +000073 square = SquareMatrixType::Random(rows, rows);
74 VectorType v1 = VectorType::Random(rows),
75 v2 = VectorType::Random(rows),
76 v3 = VectorType::Random(rows),
Gael Guennebaudc10f0692008-07-21 00:34:46 +000077 vzero = VectorType::Zero(rows);
Benoit Jacob2fdd0672007-11-28 15:34:40 +000078
Benoit Jacob47160402010-10-25 10:15:22 -040079 Scalar s1 = internal::random<Scalar>(),
80 s2 = internal::random<Scalar>();
Gael Guennebaud8e0d5482008-03-05 13:18:19 +000081
Benoit Jacob2fdd0672007-11-28 15:34:40 +000082 // check basic compatibility of adjoint, transpose, conjugate
Benoit Jacob346c00f2007-12-03 10:23:08 +000083 VERIFY_IS_APPROX(m1.transpose().conjugate().adjoint(), m1);
84 VERIFY_IS_APPROX(m1.adjoint().conjugate().transpose(), m1);
Gael Guennebaud8e0d5482008-03-05 13:18:19 +000085
Benoit Jacob2fdd0672007-11-28 15:34:40 +000086 // check multiplicative behavior
Benoit Jacob346c00f2007-12-03 10:23:08 +000087 VERIFY_IS_APPROX((m1.adjoint() * m2).adjoint(), m2.adjoint() * m1);
Gael Guennebaud62670c82013-06-10 23:40:56 +020088 VERIFY_IS_APPROX((s1 * m1).adjoint(), numext::conj(s1) * m1.adjoint());
Gael Guennebaud8e0d5482008-03-05 13:18:19 +000089
Gael Guennebauda0fb8852013-02-14 21:33:42 +010090 // check basic properties of dot, squaredNorm
Gael Guennebaud62670c82013-06-10 23:40:56 +020091 VERIFY_IS_APPROX(numext::conj(v1.dot(v2)), v2.dot(v1));
92 VERIFY_IS_APPROX(numext::real(v1.dot(v1)), v1.squaredNorm());
Gael Guennebauda0fb8852013-02-14 21:33:42 +010093
94 adjoint_specific<NumTraits<Scalar>::IsInteger>::run(v1, v2, v3, square, s1, s2);
95
Gael Guennebauda76fbbf2012-11-06 15:25:50 +010096 VERIFY_IS_MUCH_SMALLER_THAN(abs(vzero.dot(v1)), static_cast<RealScalar>(1));
Benoit Jacobd2673d82011-06-15 00:30:46 -040097
Benoit Jacobfc7b2b52007-12-12 17:48:20 +000098 // like in testBasicStuff, test operator() to check const-qualification
Benoit Jacob47160402010-10-25 10:15:22 -040099 Index r = internal::random<Index>(0, rows-1),
100 c = internal::random<Index>(0, cols-1);
Gael Guennebaud62670c82013-06-10 23:40:56 +0200101 VERIFY_IS_APPROX(m1.conjugate()(r,c), numext::conj(m1(r,c)));
102 VERIFY_IS_APPROX(m1.adjoint()(c,r), numext::conj(m1(r,c)));
Gael Guennebaud8e0d5482008-03-05 13:18:19 +0000103
Gael Guennebaudebe14aa2008-10-29 15:24:08 +0000104 // check inplace transpose
105 m3 = m1;
106 m3.transposeInPlace();
107 VERIFY_IS_APPROX(m3,m1.transpose());
108 m3.transposeInPlace();
109 VERIFY_IS_APPROX(m3,m1);
Benoit Jacobbf596d02009-03-31 13:55:40 +0000110
111 // check inplace adjoint
112 m3 = m1;
113 m3.adjointInPlace();
114 VERIFY_IS_APPROX(m3,m1.adjoint());
115 m3.transposeInPlace();
116 VERIFY_IS_APPROX(m3,m1.conjugate());
117
Gael Guennebaud0bfb78c2011-01-27 09:59:19 +0100118 // check mixed dot product
119 typedef Matrix<RealScalar, MatrixType::RowsAtCompileTime, 1> RealVectorType;
120 RealVectorType rv1 = RealVectorType::Random(rows);
121 VERIFY_IS_APPROX(v1.dot(rv1.template cast<Scalar>()), v1.dot(rv1));
122 VERIFY_IS_APPROX(rv1.template cast<Scalar>().dot(v1), rv1.dot(v1));
Benoit Jacob2fdd0672007-11-28 15:34:40 +0000123}
124
Gael Guennebaud522e24f2008-05-22 12:18:55 +0000125void test_adjoint()
Benoit Jacob2fdd0672007-11-28 15:34:40 +0000126{
Gael Guennebaud7a9519a2009-07-14 23:00:53 +0200127 for(int i = 0; i < g_repeat; i++) {
Benoit Jacob2840ac72009-10-28 18:19:29 -0400128 CALL_SUBTEST_1( adjoint(Matrix<float, 1, 1>()) );
129 CALL_SUBTEST_2( adjoint(Matrix3d()) );
130 CALL_SUBTEST_3( adjoint(Matrix4f()) );
Gael Guennebauda8f66fe2011-07-12 14:41:00 +0200131 CALL_SUBTEST_4( adjoint(MatrixXcf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE/2), internal::random<int>(1,EIGEN_TEST_MAX_SIZE/2))) );
132 CALL_SUBTEST_5( adjoint(MatrixXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
133 CALL_SUBTEST_6( adjoint(MatrixXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
Gael Guennebaudf5d23172009-07-13 21:14:47 +0200134 }
Gael Guennebauda8f66fe2011-07-12 14:41:00 +0200135 // test a large static matrix only once
Benoit Jacob2840ac72009-10-28 18:19:29 -0400136 CALL_SUBTEST_7( adjoint(Matrix<float, 100, 100>()) );
Gael Guennebaudb56bb442009-09-07 12:46:16 +0200137
Benoit Jacob2840ac72009-10-28 18:19:29 -0400138#ifdef EIGEN_TEST_PART_4
Gael Guennebaud239ada92009-08-15 22:19:29 +0200139 {
140 MatrixXcf a(10,10), b(10,10);
141 VERIFY_RAISES_ASSERT(a = a.transpose());
142 VERIFY_RAISES_ASSERT(a = a.transpose() + b);
143 VERIFY_RAISES_ASSERT(a = b + a.transpose());
144 VERIFY_RAISES_ASSERT(a = a.conjugate().transpose());
145 VERIFY_RAISES_ASSERT(a = a.adjoint());
146 VERIFY_RAISES_ASSERT(a = a.adjoint() + b);
147 VERIFY_RAISES_ASSERT(a = b + a.adjoint());
Gael Guennebaud6db67742009-12-16 11:41:16 +0100148
149 // no assertion should be triggered for these cases:
150 a.transpose() = a.transpose();
151 a.transpose() += a.transpose();
152 a.transpose() += a.transpose() + b;
153 a.transpose() = a.adjoint();
154 a.transpose() += a.adjoint();
155 a.transpose() += a.adjoint() + b;
Gael Guennebaud239ada92009-08-15 22:19:29 +0200156 }
Benoit Jacob2840ac72009-10-28 18:19:29 -0400157#endif
Benoit Jacob2fdd0672007-11-28 15:34:40 +0000158}
Benoit Jacobe05f2912007-12-02 18:32:59 +0000159