blob: f37bad78775593bbea7c76c21400e12e4f294990 [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 Jacob3698d8c2008-02-28 15:44:45 +00006// Eigen is free software; you can redistribute it and/or
7// modify it under the terms of the GNU Lesser General Public
Gael Guennebaud8e0d5482008-03-05 13:18:19 +00008// License as published by the Free Software Foundation; either
Benoit Jacob3698d8c2008-02-28 15:44:45 +00009// version 3 of the License, or (at your option) any later version.
10//
11// Alternatively, you can redistribute it and/or
12// modify it under the terms of the GNU General Public License as
Gael Guennebaud8e0d5482008-03-05 13:18:19 +000013// published by the Free Software Foundation; either version 2 of
Benoit Jacob3698d8c2008-02-28 15:44:45 +000014// the License, or (at your option) any later version.
Benoit Jacob2fdd0672007-11-28 15:34:40 +000015//
16// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
17// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
Benoit Jacob3698d8c2008-02-28 15:44:45 +000018// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
19// GNU General Public License for more details.
Benoit Jacob2fdd0672007-11-28 15:34:40 +000020//
Gael Guennebaud8e0d5482008-03-05 13:18:19 +000021// You should have received a copy of the GNU Lesser General Public
Benoit Jacob3698d8c2008-02-28 15:44:45 +000022// License and a copy of the GNU General Public License along with
23// Eigen. If not, see <http://www.gnu.org/licenses/>.
Gael Guennebaud239ada92009-08-15 22:19:29 +020024
Benoit Jacobe2775862010-04-28 18:51:38 -040025#define EIGEN_NO_STATIC_ASSERT
26
Benoit Jacob2fdd0672007-11-28 15:34:40 +000027#include "main.h"
28
29template<typename MatrixType> void adjoint(const MatrixType& m)
30{
31 /* this test covers the following files:
32 Transpose.h Conjugate.h Dot.h
33 */
Hauke Heibelf1679c72010-06-20 17:37:56 +020034 typedef typename MatrixType::Index Index;
Benoit Jacob2fdd0672007-11-28 15:34:40 +000035 typedef typename MatrixType::Scalar Scalar;
Gael Guennebaud2120fed2008-08-23 15:14:20 +000036 typedef typename NumTraits<Scalar>::Real RealScalar;
Benoit Jacob2ee68a02008-03-12 17:17:36 +000037 typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
Gael Guennebaud2120fed2008-08-23 15:14:20 +000038 typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> SquareMatrixType;
Hauke Heibelf1679c72010-06-20 17:37:56 +020039
40 Index rows = m.rows();
41 Index cols = m.cols();
Gael Guennebaud8e0d5482008-03-05 13:18:19 +000042
Gael Guennebaud2120fed2008-08-23 15:14:20 +000043 RealScalar largerEps = test_precision<RealScalar>();
Benoit Jacob47160402010-10-25 10:15:22 -040044 if (internal::is_same_type<RealScalar,float>::ret)
Gael Guennebaud2db58882009-01-12 11:55:52 +000045 largerEps = RealScalar(1e-3f);
Gael Guennebaud2120fed2008-08-23 15:14:20 +000046
Benoit Jacob46fe7a32008-09-01 17:31:21 +000047 MatrixType m1 = MatrixType::Random(rows, cols),
48 m2 = MatrixType::Random(rows, cols),
Benoit Jacob2fdd0672007-11-28 15:34:40 +000049 m3(rows, cols),
Gael Guennebaudc10f0692008-07-21 00:34:46 +000050 mzero = MatrixType::Zero(rows, cols),
Gael Guennebaud2120fed2008-08-23 15:14:20 +000051 identity = SquareMatrixType::Identity(rows, rows),
Benoit Jacob46fe7a32008-09-01 17:31:21 +000052 square = SquareMatrixType::Random(rows, rows);
53 VectorType v1 = VectorType::Random(rows),
54 v2 = VectorType::Random(rows),
55 v3 = VectorType::Random(rows),
Gael Guennebaudc10f0692008-07-21 00:34:46 +000056 vzero = VectorType::Zero(rows);
Benoit Jacob2fdd0672007-11-28 15:34:40 +000057
Benoit Jacob47160402010-10-25 10:15:22 -040058 Scalar s1 = internal::random<Scalar>(),
59 s2 = internal::random<Scalar>();
Gael Guennebaud8e0d5482008-03-05 13:18:19 +000060
Benoit Jacob2fdd0672007-11-28 15:34:40 +000061 // check basic compatibility of adjoint, transpose, conjugate
Benoit Jacob346c00f2007-12-03 10:23:08 +000062 VERIFY_IS_APPROX(m1.transpose().conjugate().adjoint(), m1);
63 VERIFY_IS_APPROX(m1.adjoint().conjugate().transpose(), m1);
Gael Guennebaud8e0d5482008-03-05 13:18:19 +000064
Benoit Jacob2fdd0672007-11-28 15:34:40 +000065 // check multiplicative behavior
Benoit Jacob346c00f2007-12-03 10:23:08 +000066 VERIFY_IS_APPROX((m1.adjoint() * m2).adjoint(), m2.adjoint() * m1);
Benoit Jacob47160402010-10-25 10:15:22 -040067 VERIFY_IS_APPROX((s1 * m1).adjoint(), internal::conj(s1) * m1.adjoint());
Gael Guennebaud8e0d5482008-03-05 13:18:19 +000068
Benoit Jacob2fdd0672007-11-28 15:34:40 +000069 // check basic properties of dot, norm, norm2
70 typedef typename NumTraits<Scalar>::Real RealScalar;
Benoit Jacob47160402010-10-25 10:15:22 -040071 VERIFY(internal::isApprox((s1 * v1 + s2 * v2).dot(v3), internal::conj(s1) * v1.dot(v3) + internal::conj(s2) * v2.dot(v3), largerEps));
72 VERIFY(internal::isApprox(v3.dot(s1 * v1 + s2 * v2), s1*v3.dot(v1)+s2*v3.dot(v2), largerEps));
73 VERIFY_IS_APPROX(internal::conj(v1.dot(v2)), v2.dot(v1));
74 VERIFY_IS_APPROX(internal::abs(v1.dot(v1)), v1.squaredNorm());
Benoit Jacobe2775862010-04-28 18:51:38 -040075 if(!NumTraits<Scalar>::IsInteger)
Gael Guennebaud36c8a642009-01-28 22:11:56 +000076 VERIFY_IS_APPROX(v1.squaredNorm(), v1.norm() * v1.norm());
Benoit Jacob47160402010-10-25 10:15:22 -040077 VERIFY_IS_MUCH_SMALLER_THAN(internal::abs(vzero.dot(v1)), static_cast<RealScalar>(1));
Gael Guennebaud8e0d5482008-03-05 13:18:19 +000078
Benoit Jacob2fdd0672007-11-28 15:34:40 +000079 // check compatibility of dot and adjoint
Benoit Jacob47160402010-10-25 10:15:22 -040080 VERIFY(internal::isApprox(v1.dot(square * v2), (square.adjoint() * v1).dot(v2), largerEps));
Gael Guennebaud8e0d5482008-03-05 13:18:19 +000081
Benoit Jacobfc7b2b52007-12-12 17:48:20 +000082 // like in testBasicStuff, test operator() to check const-qualification
Benoit Jacob47160402010-10-25 10:15:22 -040083 Index r = internal::random<Index>(0, rows-1),
84 c = internal::random<Index>(0, cols-1);
85 VERIFY_IS_APPROX(m1.conjugate()(r,c), internal::conj(m1(r,c)));
86 VERIFY_IS_APPROX(m1.adjoint()(c,r), internal::conj(m1(r,c)));
Gael Guennebaud8e0d5482008-03-05 13:18:19 +000087
Benoit Jacobe2775862010-04-28 18:51:38 -040088 if(!NumTraits<Scalar>::IsInteger)
Benoit Jacob13ad8872008-08-12 02:14:02 +000089 {
90 // check that Random().normalized() works: tricky as the random xpr must be evaluated by
91 // normalized() in order to produce a consistent result.
92 VERIFY_IS_APPROX(VectorType::Random(rows).normalized().norm(), RealScalar(1));
93 }
Gael Guennebaudebe14aa2008-10-29 15:24:08 +000094
95 // check inplace transpose
96 m3 = m1;
97 m3.transposeInPlace();
98 VERIFY_IS_APPROX(m3,m1.transpose());
99 m3.transposeInPlace();
100 VERIFY_IS_APPROX(m3,m1);
Benoit Jacobbf596d02009-03-31 13:55:40 +0000101
102 // check inplace adjoint
103 m3 = m1;
104 m3.adjointInPlace();
105 VERIFY_IS_APPROX(m3,m1.adjoint());
106 m3.transposeInPlace();
107 VERIFY_IS_APPROX(m3,m1.conjugate());
108
Benoit Jacob2fdd0672007-11-28 15:34:40 +0000109}
110
Gael Guennebaud522e24f2008-05-22 12:18:55 +0000111void test_adjoint()
Benoit Jacob2fdd0672007-11-28 15:34:40 +0000112{
Gael Guennebaud7a9519a2009-07-14 23:00:53 +0200113 for(int i = 0; i < g_repeat; i++) {
Benoit Jacob2840ac72009-10-28 18:19:29 -0400114 CALL_SUBTEST_1( adjoint(Matrix<float, 1, 1>()) );
115 CALL_SUBTEST_2( adjoint(Matrix3d()) );
116 CALL_SUBTEST_3( adjoint(Matrix4f()) );
117 CALL_SUBTEST_4( adjoint(MatrixXcf(4, 4)) );
118 CALL_SUBTEST_5( adjoint(MatrixXi(8, 12)) );
119 CALL_SUBTEST_6( adjoint(MatrixXf(21, 21)) );
Gael Guennebaudf5d23172009-07-13 21:14:47 +0200120 }
Gael Guennebaud7a9519a2009-07-14 23:00:53 +0200121 // test a large matrix only once
Benoit Jacob2840ac72009-10-28 18:19:29 -0400122 CALL_SUBTEST_7( adjoint(Matrix<float, 100, 100>()) );
Gael Guennebaudb56bb442009-09-07 12:46:16 +0200123
Benoit Jacob2840ac72009-10-28 18:19:29 -0400124#ifdef EIGEN_TEST_PART_4
Gael Guennebaud239ada92009-08-15 22:19:29 +0200125 {
126 MatrixXcf a(10,10), b(10,10);
127 VERIFY_RAISES_ASSERT(a = a.transpose());
128 VERIFY_RAISES_ASSERT(a = a.transpose() + b);
129 VERIFY_RAISES_ASSERT(a = b + a.transpose());
130 VERIFY_RAISES_ASSERT(a = a.conjugate().transpose());
131 VERIFY_RAISES_ASSERT(a = a.adjoint());
132 VERIFY_RAISES_ASSERT(a = a.adjoint() + b);
133 VERIFY_RAISES_ASSERT(a = b + a.adjoint());
Gael Guennebaud6db67742009-12-16 11:41:16 +0100134
135 // no assertion should be triggered for these cases:
136 a.transpose() = a.transpose();
137 a.transpose() += a.transpose();
138 a.transpose() += a.transpose() + b;
139 a.transpose() = a.adjoint();
140 a.transpose() += a.adjoint();
141 a.transpose() += a.adjoint() + b;
Gael Guennebaud239ada92009-08-15 22:19:29 +0200142 }
Benoit Jacob2840ac72009-10-28 18:19:29 -0400143#endif
Benoit Jacob2fdd0672007-11-28 15:34:40 +0000144}
Benoit Jacobe05f2912007-12-02 18:32:59 +0000145