blob: 3006503389ad52cbff73e8913a5fa79784732d0c [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/>.
Benoit Jacobbf596d02009-03-31 13:55:40 +000024#define EIGEN_NO_ASSERTION_CHECKING
Benoit Jacob2fdd0672007-11-28 15:34:40 +000025#include "main.h"
26
27template<typename MatrixType> void adjoint(const MatrixType& m)
28{
29 /* this test covers the following files:
30 Transpose.h Conjugate.h Dot.h
31 */
32
33 typedef typename MatrixType::Scalar Scalar;
Gael Guennebaud2120fed2008-08-23 15:14:20 +000034 typedef typename NumTraits<Scalar>::Real RealScalar;
Benoit Jacob2ee68a02008-03-12 17:17:36 +000035 typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
Gael Guennebaud2120fed2008-08-23 15:14:20 +000036 typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> SquareMatrixType;
Benoit Jacob2fdd0672007-11-28 15:34:40 +000037 int rows = m.rows();
38 int cols = m.cols();
Gael Guennebaud8e0d5482008-03-05 13:18:19 +000039
Gael Guennebaud2120fed2008-08-23 15:14:20 +000040 RealScalar largerEps = test_precision<RealScalar>();
41 if (ei_is_same_type<RealScalar,float>::ret)
Gael Guennebaud2db58882009-01-12 11:55:52 +000042 largerEps = RealScalar(1e-3f);
Gael Guennebaud2120fed2008-08-23 15:14:20 +000043
Benoit Jacob46fe7a32008-09-01 17:31:21 +000044 MatrixType m1 = MatrixType::Random(rows, cols),
45 m2 = MatrixType::Random(rows, cols),
Benoit Jacob2fdd0672007-11-28 15:34:40 +000046 m3(rows, cols),
Gael Guennebaudc10f0692008-07-21 00:34:46 +000047 mzero = MatrixType::Zero(rows, cols),
Gael Guennebaud2120fed2008-08-23 15:14:20 +000048 identity = SquareMatrixType::Identity(rows, rows),
Benoit Jacob46fe7a32008-09-01 17:31:21 +000049 square = SquareMatrixType::Random(rows, rows);
50 VectorType v1 = VectorType::Random(rows),
51 v2 = VectorType::Random(rows),
52 v3 = VectorType::Random(rows),
Gael Guennebaudc10f0692008-07-21 00:34:46 +000053 vzero = VectorType::Zero(rows);
Benoit Jacob2fdd0672007-11-28 15:34:40 +000054
Benoit Jacob46fe7a32008-09-01 17:31:21 +000055 Scalar s1 = ei_random<Scalar>(),
56 s2 = ei_random<Scalar>();
Gael Guennebaud8e0d5482008-03-05 13:18:19 +000057
Benoit Jacob2fdd0672007-11-28 15:34:40 +000058 // check basic compatibility of adjoint, transpose, conjugate
Benoit Jacob346c00f2007-12-03 10:23:08 +000059 VERIFY_IS_APPROX(m1.transpose().conjugate().adjoint(), m1);
60 VERIFY_IS_APPROX(m1.adjoint().conjugate().transpose(), m1);
Gael Guennebaud8e0d5482008-03-05 13:18:19 +000061
Benoit Jacob2fdd0672007-11-28 15:34:40 +000062 // check multiplicative behavior
Benoit Jacob346c00f2007-12-03 10:23:08 +000063 VERIFY_IS_APPROX((m1.adjoint() * m2).adjoint(), m2.adjoint() * m1);
Benoit Jacob69078862008-02-28 12:38:12 +000064 VERIFY_IS_APPROX((s1 * m1).adjoint(), ei_conj(s1) * m1.adjoint());
Gael Guennebaud8e0d5482008-03-05 13:18:19 +000065
Benoit Jacob2fdd0672007-11-28 15:34:40 +000066 // check basic properties of dot, norm, norm2
67 typedef typename NumTraits<Scalar>::Real RealScalar;
Benoit Jacob523cded2009-08-03 17:20:45 +020068 VERIFY(ei_isApprox((s1 * v1 + s2 * v2).dot(v3), ei_conj(s1) * v1.dot(v3) + ei_conj(s2) * v2.dot(v3), largerEps));
69 VERIFY(ei_isApprox(v3.dot(s1 * v1 + s2 * v2), s1*v3.dot(v1)+s2*v3.dot(v2), largerEps));
Gael Guennebaud2120fed2008-08-23 15:14:20 +000070 VERIFY_IS_APPROX(ei_conj(v1.dot(v2)), v2.dot(v1));
Benoit Jacob3d90c132008-11-03 19:14:17 +000071 VERIFY_IS_APPROX(ei_abs(v1.dot(v1)), v1.squaredNorm());
Benoit Jacobe05f2912007-12-02 18:32:59 +000072 if(NumTraits<Scalar>::HasFloatingPoint)
Gael Guennebaud36c8a642009-01-28 22:11:56 +000073 VERIFY_IS_APPROX(v1.squaredNorm(), v1.norm() * v1.norm());
Gael Guennebaud2120fed2008-08-23 15:14:20 +000074 VERIFY_IS_MUCH_SMALLER_THAN(ei_abs(vzero.dot(v1)), static_cast<RealScalar>(1));
Benoit Jacobe05f2912007-12-02 18:32:59 +000075 if(NumTraits<Scalar>::HasFloatingPoint)
Gael Guennebaud36c8a642009-01-28 22:11:56 +000076 {
Gael Guennebaud2120fed2008-08-23 15:14:20 +000077 VERIFY_IS_MUCH_SMALLER_THAN(vzero.norm(), static_cast<RealScalar>(1));
Gael Guennebaud36c8a642009-01-28 22:11:56 +000078 VERIFY_IS_APPROX(v1.norm(), v1.stableNorm());
Gael Guennebaud32b08ac2009-07-17 16:22:39 +020079 VERIFY_IS_APPROX(v1.blueNorm(), v1.stableNorm());
80 VERIFY_IS_APPROX(v1.hypotNorm(), v1.stableNorm());
Gael Guennebaud36c8a642009-01-28 22:11:56 +000081 }
Gael Guennebaud8e0d5482008-03-05 13:18:19 +000082
Benoit Jacob2fdd0672007-11-28 15:34:40 +000083 // check compatibility of dot and adjoint
Gael Guennebaud2120fed2008-08-23 15:14:20 +000084 VERIFY(ei_isApprox(v1.dot(square * v2), (square.adjoint() * v1).dot(v2), largerEps));
Gael Guennebaud8e0d5482008-03-05 13:18:19 +000085
Benoit Jacobfc7b2b52007-12-12 17:48:20 +000086 // like in testBasicStuff, test operator() to check const-qualification
Benoit Jacob69078862008-02-28 12:38:12 +000087 int r = ei_random<int>(0, rows-1),
88 c = ei_random<int>(0, cols-1);
89 VERIFY_IS_APPROX(m1.conjugate()(r,c), ei_conj(m1(r,c)));
90 VERIFY_IS_APPROX(m1.adjoint()(c,r), ei_conj(m1(r,c)));
Gael Guennebaud8e0d5482008-03-05 13:18:19 +000091
Benoit Jacob13ad8872008-08-12 02:14:02 +000092 if(NumTraits<Scalar>::HasFloatingPoint)
93 {
94 // check that Random().normalized() works: tricky as the random xpr must be evaluated by
95 // normalized() in order to produce a consistent result.
96 VERIFY_IS_APPROX(VectorType::Random(rows).normalized().norm(), RealScalar(1));
97 }
Gael Guennebaudebe14aa2008-10-29 15:24:08 +000098
99 // check inplace transpose
100 m3 = m1;
101 m3.transposeInPlace();
102 VERIFY_IS_APPROX(m3,m1.transpose());
103 m3.transposeInPlace();
104 VERIFY_IS_APPROX(m3,m1);
Benoit Jacobbf596d02009-03-31 13:55:40 +0000105
106 // check inplace adjoint
107 m3 = m1;
108 m3.adjointInPlace();
109 VERIFY_IS_APPROX(m3,m1.adjoint());
110 m3.transposeInPlace();
111 VERIFY_IS_APPROX(m3,m1.conjugate());
112
Gael Guennebaud15784212009-07-15 14:20:45 +0200113
Benoit Jacob2fdd0672007-11-28 15:34:40 +0000114}
115
Gael Guennebaud522e24f2008-05-22 12:18:55 +0000116void test_adjoint()
Benoit Jacob2fdd0672007-11-28 15:34:40 +0000117{
Gael Guennebaud7a9519a2009-07-14 23:00:53 +0200118 for(int i = 0; i < g_repeat; i++) {
119 CALL_SUBTEST( adjoint(Matrix<float, 1, 1>()) );
120 CALL_SUBTEST( adjoint(Matrix3d()) );
121 CALL_SUBTEST( adjoint(Matrix4f()) );
122 CALL_SUBTEST( adjoint(MatrixXcf(4, 4)) );
123 CALL_SUBTEST( adjoint(MatrixXi(8, 12)) );
124 CALL_SUBTEST( adjoint(MatrixXf(21, 21)) );
Gael Guennebaudf5d23172009-07-13 21:14:47 +0200125 }
Gael Guennebaud7a9519a2009-07-14 23:00:53 +0200126 // test a large matrix only once
127 CALL_SUBTEST( adjoint(Matrix<float, 100, 100>()) );
Benoit Jacob2fdd0672007-11-28 15:34:40 +0000128}
Benoit Jacobe05f2912007-12-02 18:32:59 +0000129