Benoit Jacob | 2fdd067 | 2007-11-28 15:34:40 +0000 | [diff] [blame] | 1 | // This file is part of Eigen, a lightweight C++ template library |
Benoit Jacob | 6347b1d | 2009-05-22 20:25:33 +0200 | [diff] [blame] | 2 | // for linear algebra. |
Benoit Jacob | 2fdd067 | 2007-11-28 15:34:40 +0000 | [diff] [blame] | 3 | // |
Benoit Jacob | 00f89a8 | 2008-11-24 13:40:43 +0000 | [diff] [blame] | 4 | // Copyright (C) 2006-2008 Benoit Jacob <jacob.benoit.1@gmail.com> |
Benoit Jacob | 2fdd067 | 2007-11-28 15:34:40 +0000 | [diff] [blame] | 5 | // |
Benoit Jacob | 69124cf | 2012-07-13 14:42:47 -0400 | [diff] [blame] | 6 | // 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 Guennebaud | 239ada9 | 2009-08-15 22:19:29 +0200 | [diff] [blame] | 9 | |
Benoit Jacob | 2fdd067 | 2007-11-28 15:34:40 +0000 | [diff] [blame] | 10 | #include "main.h" |
| 11 | |
Gael Guennebaud | a0fb885 | 2013-02-14 21:33:42 +0100 | [diff] [blame] | 12 | template<bool IsInteger> struct adjoint_specific; |
| 13 | |
| 14 | template<> struct adjoint_specific<true> { |
| 15 | template<typename Vec, typename Mat, typename Scalar> |
| 16 | static void run(const Vec& v1, const Vec& v2, Vec& v3, const Mat& square, Scalar s1, Scalar s2) { |
Gael Guennebaud | 62670c8 | 2013-06-10 23:40:56 +0200 | [diff] [blame] | 17 | VERIFY(test_isApproxWithRef((s1 * v1 + s2 * v2).dot(v3), numext::conj(s1) * v1.dot(v3) + numext::conj(s2) * v2.dot(v3), 0)); |
Gael Guennebaud | a0fb885 | 2013-02-14 21:33:42 +0100 | [diff] [blame] | 18 | VERIFY(test_isApproxWithRef(v3.dot(s1 * v1 + s2 * v2), s1*v3.dot(v1)+s2*v3.dot(v2), 0)); |
| 19 | |
| 20 | // check compatibility of dot and adjoint |
| 21 | VERIFY(test_isApproxWithRef(v1.dot(square * v2), (square.adjoint() * v1).dot(v2), 0)); |
| 22 | } |
| 23 | }; |
| 24 | |
| 25 | template<> struct adjoint_specific<false> { |
| 26 | template<typename Vec, typename Mat, typename Scalar> |
| 27 | static void run(const Vec& v1, const Vec& v2, Vec& v3, const Mat& square, Scalar s1, Scalar s2) { |
| 28 | typedef typename NumTraits<Scalar>::Real RealScalar; |
Gael Guennebaud | d02e329 | 2013-07-15 21:21:14 +0200 | [diff] [blame] | 29 | using std::abs; |
Gael Guennebaud | a0fb885 | 2013-02-14 21:33:42 +0100 | [diff] [blame] | 30 | |
| 31 | RealScalar ref = NumTraits<Scalar>::IsInteger ? RealScalar(0) : (std::max)((s1 * v1 + s2 * v2).norm(),v3.norm()); |
Gael Guennebaud | 62670c8 | 2013-06-10 23:40:56 +0200 | [diff] [blame] | 32 | VERIFY(test_isApproxWithRef((s1 * v1 + s2 * v2).dot(v3), numext::conj(s1) * v1.dot(v3) + numext::conj(s2) * v2.dot(v3), ref)); |
Gael Guennebaud | a0fb885 | 2013-02-14 21:33:42 +0100 | [diff] [blame] | 33 | VERIFY(test_isApproxWithRef(v3.dot(s1 * v1 + s2 * v2), s1*v3.dot(v1)+s2*v3.dot(v2), ref)); |
| 34 | |
| 35 | VERIFY_IS_APPROX(v1.squaredNorm(), v1.norm() * v1.norm()); |
| 36 | // check normalized() and normalize() |
| 37 | VERIFY_IS_APPROX(v1, v1.norm() * v1.normalized()); |
| 38 | v3 = v1; |
| 39 | v3.normalize(); |
| 40 | VERIFY_IS_APPROX(v1, v1.norm() * v3); |
| 41 | VERIFY_IS_APPROX(v3, v1.normalized()); |
| 42 | VERIFY_IS_APPROX(v3.norm(), RealScalar(1)); |
Gael Guennebaud | ee37eb4 | 2016-01-21 20:43:42 +0100 | [diff] [blame] | 43 | |
| 44 | // check null inputs |
| 45 | VERIFY_IS_APPROX((v1*0).normalized(), (v1*0)); |
Gael Guennebaud | 3ba8a3a | 2016-01-30 22:14:04 +0100 | [diff] [blame] | 46 | #if (!EIGEN_ARCH_i386) || defined(EIGEN_VECTORIZE) |
Gael Guennebaud | ee37eb4 | 2016-01-21 20:43:42 +0100 | [diff] [blame] | 47 | RealScalar very_small = (std::numeric_limits<RealScalar>::min)(); |
Erik Schultheis | d271a7d | 2022-01-26 18:16:19 +0000 | [diff] [blame] | 48 | VERIFY( numext::is_exactly_zero((v1*very_small).norm()) ); |
Gael Guennebaud | ee37eb4 | 2016-01-21 20:43:42 +0100 | [diff] [blame] | 49 | VERIFY_IS_APPROX((v1*very_small).normalized(), (v1*very_small)); |
| 50 | v3 = v1*very_small; |
| 51 | v3.normalize(); |
| 52 | VERIFY_IS_APPROX(v3, (v1*very_small)); |
Gael Guennebaud | 3ba8a3a | 2016-01-30 22:14:04 +0100 | [diff] [blame] | 53 | #endif |
Gael Guennebaud | a0fb885 | 2013-02-14 21:33:42 +0100 | [diff] [blame] | 54 | |
| 55 | // check compatibility of dot and adjoint |
| 56 | ref = NumTraits<Scalar>::IsInteger ? 0 : (std::max)((std::max)(v1.norm(),v2.norm()),(std::max)((square * v2).norm(),(square.adjoint() * v1).norm())); |
Gael Guennebaud | d02e329 | 2013-07-15 21:21:14 +0200 | [diff] [blame] | 57 | VERIFY(internal::isMuchSmallerThan(abs(v1.dot(square * v2) - (square.adjoint() * v1).dot(v2)), ref, test_precision<Scalar>())); |
Gael Guennebaud | a0fb885 | 2013-02-14 21:33:42 +0100 | [diff] [blame] | 58 | |
| 59 | // check that Random().normalized() works: tricky as the random xpr must be evaluated by |
| 60 | // normalized() in order to produce a consistent result. |
| 61 | VERIFY_IS_APPROX(Vec::Random(v1.size()).normalized().norm(), RealScalar(1)); |
| 62 | } |
| 63 | }; |
| 64 | |
Antonio Sánchez | 3234809 | 2022-05-23 14:46:16 +0000 | [diff] [blame] | 65 | template<typename MatrixType, typename Scalar = typename MatrixType::Scalar> |
| 66 | MatrixType RandomMatrix(int rows, int cols, Scalar min, Scalar max) { |
| 67 | MatrixType M = MatrixType(rows, cols); |
| 68 | for (int i=0; i<rows; ++i) { |
| 69 | for (int j=0; j<cols; ++j) { |
| 70 | M(i, j) = Eigen::internal::random<Scalar>(min, max); |
| 71 | } |
| 72 | } |
| 73 | return M; |
| 74 | } |
| 75 | |
Benoit Jacob | 2fdd067 | 2007-11-28 15:34:40 +0000 | [diff] [blame] | 76 | template<typename MatrixType> void adjoint(const MatrixType& m) |
| 77 | { |
| 78 | /* this test covers the following files: |
| 79 | Transpose.h Conjugate.h Dot.h |
| 80 | */ |
Gael Guennebaud | a76fbbf | 2012-11-06 15:25:50 +0100 | [diff] [blame] | 81 | using std::abs; |
Benoit Jacob | 2fdd067 | 2007-11-28 15:34:40 +0000 | [diff] [blame] | 82 | typedef typename MatrixType::Scalar Scalar; |
Gael Guennebaud | 2120fed | 2008-08-23 15:14:20 +0000 | [diff] [blame] | 83 | typedef typename NumTraits<Scalar>::Real RealScalar; |
Benoit Jacob | 2ee68a0 | 2008-03-12 17:17:36 +0000 | [diff] [blame] | 84 | typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType; |
Gael Guennebaud | 2120fed | 2008-08-23 15:14:20 +0000 | [diff] [blame] | 85 | typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> SquareMatrixType; |
Gael Guennebaud | c6eb84a | 2015-01-26 17:09:01 +0100 | [diff] [blame] | 86 | const Index PacketSize = internal::packet_traits<Scalar>::size; |
Hauke Heibel | f1679c7 | 2010-06-20 17:37:56 +0200 | [diff] [blame] | 87 | |
| 88 | Index rows = m.rows(); |
| 89 | Index cols = m.cols(); |
Gael Guennebaud | 8e0d548 | 2008-03-05 13:18:19 +0000 | [diff] [blame] | 90 | |
Antonio Sánchez | 3234809 | 2022-05-23 14:46:16 +0000 | [diff] [blame] | 91 | // Avoid integer overflow by limiting input values. |
| 92 | RealScalar rmin = static_cast<RealScalar>(NumTraits<Scalar>::IsInteger ? NumTraits<Scalar>::IsSigned ? -100 : 0 : -1); |
| 93 | RealScalar rmax = static_cast<RealScalar>(NumTraits<Scalar>::IsInteger ? 100 : 1); |
| 94 | |
| 95 | MatrixType m1 = RandomMatrix<MatrixType>(rows, cols, rmin, rmax), |
| 96 | m2 = RandomMatrix<MatrixType>(rows, cols, rmin, rmax), |
Benoit Jacob | 2fdd067 | 2007-11-28 15:34:40 +0000 | [diff] [blame] | 97 | m3(rows, cols), |
Antonio Sánchez | 3234809 | 2022-05-23 14:46:16 +0000 | [diff] [blame] | 98 | square = RandomMatrix<SquareMatrixType>(rows, rows, rmin, rmax); |
| 99 | VectorType v1 = RandomMatrix<VectorType>(rows, 1, rmin, rmax), |
| 100 | v2 = RandomMatrix<VectorType>(rows, 1, rmin, rmax), |
| 101 | v3 = RandomMatrix<VectorType>(rows, 1, rmin, rmax), |
Gael Guennebaud | c10f069 | 2008-07-21 00:34:46 +0000 | [diff] [blame] | 102 | vzero = VectorType::Zero(rows); |
Benoit Jacob | 2fdd067 | 2007-11-28 15:34:40 +0000 | [diff] [blame] | 103 | |
Antonio Sánchez | 3234809 | 2022-05-23 14:46:16 +0000 | [diff] [blame] | 104 | Scalar s1 = internal::random<Scalar>(rmin, rmax), |
| 105 | s2 = internal::random<Scalar>(rmin, rmax); |
Gael Guennebaud | 8e0d548 | 2008-03-05 13:18:19 +0000 | [diff] [blame] | 106 | |
Benoit Jacob | 2fdd067 | 2007-11-28 15:34:40 +0000 | [diff] [blame] | 107 | // check basic compatibility of adjoint, transpose, conjugate |
Benoit Jacob | 346c00f | 2007-12-03 10:23:08 +0000 | [diff] [blame] | 108 | VERIFY_IS_APPROX(m1.transpose().conjugate().adjoint(), m1); |
| 109 | VERIFY_IS_APPROX(m1.adjoint().conjugate().transpose(), m1); |
Gael Guennebaud | 8e0d548 | 2008-03-05 13:18:19 +0000 | [diff] [blame] | 110 | |
Benoit Jacob | 2fdd067 | 2007-11-28 15:34:40 +0000 | [diff] [blame] | 111 | // check multiplicative behavior |
Benoit Jacob | 346c00f | 2007-12-03 10:23:08 +0000 | [diff] [blame] | 112 | VERIFY_IS_APPROX((m1.adjoint() * m2).adjoint(), m2.adjoint() * m1); |
Gael Guennebaud | 62670c8 | 2013-06-10 23:40:56 +0200 | [diff] [blame] | 113 | VERIFY_IS_APPROX((s1 * m1).adjoint(), numext::conj(s1) * m1.adjoint()); |
Gael Guennebaud | 8e0d548 | 2008-03-05 13:18:19 +0000 | [diff] [blame] | 114 | |
Gael Guennebaud | a0fb885 | 2013-02-14 21:33:42 +0100 | [diff] [blame] | 115 | // check basic properties of dot, squaredNorm |
Gael Guennebaud | 62670c8 | 2013-06-10 23:40:56 +0200 | [diff] [blame] | 116 | VERIFY_IS_APPROX(numext::conj(v1.dot(v2)), v2.dot(v1)); |
| 117 | VERIFY_IS_APPROX(numext::real(v1.dot(v1)), v1.squaredNorm()); |
Gael Guennebaud | a0fb885 | 2013-02-14 21:33:42 +0100 | [diff] [blame] | 118 | |
| 119 | adjoint_specific<NumTraits<Scalar>::IsInteger>::run(v1, v2, v3, square, s1, s2); |
| 120 | |
Gael Guennebaud | a76fbbf | 2012-11-06 15:25:50 +0100 | [diff] [blame] | 121 | VERIFY_IS_MUCH_SMALLER_THAN(abs(vzero.dot(v1)), static_cast<RealScalar>(1)); |
Benoit Jacob | d2673d8 | 2011-06-15 00:30:46 -0400 | [diff] [blame] | 122 | |
Benoit Jacob | fc7b2b5 | 2007-12-12 17:48:20 +0000 | [diff] [blame] | 123 | // like in testBasicStuff, test operator() to check const-qualification |
Benoit Jacob | 4716040 | 2010-10-25 10:15:22 -0400 | [diff] [blame] | 124 | Index r = internal::random<Index>(0, rows-1), |
| 125 | c = internal::random<Index>(0, cols-1); |
Gael Guennebaud | 62670c8 | 2013-06-10 23:40:56 +0200 | [diff] [blame] | 126 | VERIFY_IS_APPROX(m1.conjugate()(r,c), numext::conj(m1(r,c))); |
| 127 | VERIFY_IS_APPROX(m1.adjoint()(c,r), numext::conj(m1(r,c))); |
Gael Guennebaud | 8e0d548 | 2008-03-05 13:18:19 +0000 | [diff] [blame] | 128 | |
Gael Guennebaud | ebe14aa | 2008-10-29 15:24:08 +0000 | [diff] [blame] | 129 | // check inplace transpose |
| 130 | m3 = m1; |
| 131 | m3.transposeInPlace(); |
| 132 | VERIFY_IS_APPROX(m3,m1.transpose()); |
| 133 | m3.transposeInPlace(); |
| 134 | VERIFY_IS_APPROX(m3,m1); |
Gael Guennebaud | c6eb84a | 2015-01-26 17:09:01 +0100 | [diff] [blame] | 135 | |
| 136 | if(PacketSize<m3.rows() && PacketSize<m3.cols()) |
| 137 | { |
| 138 | m3 = m1; |
| 139 | Index i = internal::random<Index>(0,m3.rows()-PacketSize); |
| 140 | Index j = internal::random<Index>(0,m3.cols()-PacketSize); |
| 141 | m3.template block<PacketSize,PacketSize>(i,j).transposeInPlace(); |
| 142 | VERIFY_IS_APPROX( (m3.template block<PacketSize,PacketSize>(i,j)), (m1.template block<PacketSize,PacketSize>(i,j).transpose()) ); |
| 143 | m3.template block<PacketSize,PacketSize>(i,j).transposeInPlace(); |
| 144 | VERIFY_IS_APPROX(m3,m1); |
| 145 | } |
Benoit Jacob | bf596d0 | 2009-03-31 13:55:40 +0000 | [diff] [blame] | 146 | |
| 147 | // check inplace adjoint |
| 148 | m3 = m1; |
| 149 | m3.adjointInPlace(); |
| 150 | VERIFY_IS_APPROX(m3,m1.adjoint()); |
| 151 | m3.transposeInPlace(); |
| 152 | VERIFY_IS_APPROX(m3,m1.conjugate()); |
| 153 | |
Gael Guennebaud | 0bfb78c | 2011-01-27 09:59:19 +0100 | [diff] [blame] | 154 | // check mixed dot product |
| 155 | typedef Matrix<RealScalar, MatrixType::RowsAtCompileTime, 1> RealVectorType; |
Antonio Sánchez | 3234809 | 2022-05-23 14:46:16 +0000 | [diff] [blame] | 156 | RealVectorType rv1 = RandomMatrix<RealVectorType>(rows, 1, rmin, rmax); |
| 157 | |
Gael Guennebaud | 0bfb78c | 2011-01-27 09:59:19 +0100 | [diff] [blame] | 158 | VERIFY_IS_APPROX(v1.dot(rv1.template cast<Scalar>()), v1.dot(rv1)); |
| 159 | VERIFY_IS_APPROX(rv1.template cast<Scalar>().dot(v1), rv1.dot(v1)); |
Gael Guennebaud | 7f32109 | 2019-01-17 11:33:43 +0100 | [diff] [blame] | 160 | |
| 161 | VERIFY( is_same_type(m1,m1.template conjugateIf<false>()) ); |
| 162 | VERIFY( is_same_type(m1.conjugate(),m1.template conjugateIf<true>()) ); |
Benoit Jacob | 2fdd067 | 2007-11-28 15:34:40 +0000 | [diff] [blame] | 163 | } |
| 164 | |
Gael Guennebaud | dff3a92 | 2018-07-17 15:52:58 +0200 | [diff] [blame] | 165 | template<int> |
| 166 | void adjoint_extra() |
| 167 | { |
| 168 | MatrixXcf a(10,10), b(10,10); |
| 169 | VERIFY_RAISES_ASSERT(a = a.transpose()); |
| 170 | VERIFY_RAISES_ASSERT(a = a.transpose() + b); |
| 171 | VERIFY_RAISES_ASSERT(a = b + a.transpose()); |
| 172 | VERIFY_RAISES_ASSERT(a = a.conjugate().transpose()); |
| 173 | VERIFY_RAISES_ASSERT(a = a.adjoint()); |
| 174 | VERIFY_RAISES_ASSERT(a = a.adjoint() + b); |
| 175 | VERIFY_RAISES_ASSERT(a = b + a.adjoint()); |
| 176 | |
| 177 | // no assertion should be triggered for these cases: |
| 178 | a.transpose() = a.transpose(); |
| 179 | a.transpose() += a.transpose(); |
| 180 | a.transpose() += a.transpose() + b; |
| 181 | a.transpose() = a.adjoint(); |
| 182 | a.transpose() += a.adjoint(); |
| 183 | a.transpose() += a.adjoint() + b; |
| 184 | |
| 185 | // regression tests for check_for_aliasing |
| 186 | MatrixXd c(10,10); |
| 187 | c = 1.0 * MatrixXd::Ones(10,10) + c; |
| 188 | c = MatrixXd::Ones(10,10) * 1.0 + c; |
| 189 | c = c + MatrixXd::Ones(10,10) .cwiseProduct( MatrixXd::Zero(10,10) ); |
| 190 | c = MatrixXd::Ones(10,10) * MatrixXd::Zero(10,10); |
Gael Guennebaud | 502f717 | 2019-01-16 14:33:45 +0100 | [diff] [blame] | 191 | |
| 192 | // regression for bug 1646 |
| 193 | for (int j = 0; j < 10; ++j) { |
| 194 | c.col(j).head(j) = c.row(j).head(j); |
| 195 | } |
| 196 | |
Gael Guennebaud | c8e40ed | 2019-01-16 16:27:00 +0100 | [diff] [blame] | 197 | for (int j = 0; j < 10; ++j) { |
| 198 | c.col(j) = c.row(j); |
| 199 | } |
| 200 | |
Gael Guennebaud | 502f717 | 2019-01-16 14:33:45 +0100 | [diff] [blame] | 201 | a.conservativeResize(1,1); |
| 202 | a = a.transpose(); |
| 203 | |
| 204 | a.conservativeResize(0,0); |
| 205 | a = a.transpose(); |
Gael Guennebaud | dff3a92 | 2018-07-17 15:52:58 +0200 | [diff] [blame] | 206 | } |
| 207 | |
Gael Guennebaud | 82f0ce2 | 2018-07-17 14:46:15 +0200 | [diff] [blame] | 208 | EIGEN_DECLARE_TEST(adjoint) |
Benoit Jacob | 2fdd067 | 2007-11-28 15:34:40 +0000 | [diff] [blame] | 209 | { |
Gael Guennebaud | 7a9519a | 2009-07-14 23:00:53 +0200 | [diff] [blame] | 210 | for(int i = 0; i < g_repeat; i++) { |
Benoit Jacob | 2840ac7 | 2009-10-28 18:19:29 -0400 | [diff] [blame] | 211 | CALL_SUBTEST_1( adjoint(Matrix<float, 1, 1>()) ); |
| 212 | CALL_SUBTEST_2( adjoint(Matrix3d()) ); |
| 213 | CALL_SUBTEST_3( adjoint(Matrix4f()) ); |
Gael Guennebaud | c6eb84a | 2015-01-26 17:09:01 +0100 | [diff] [blame] | 214 | |
Gael Guennebaud | a8f66fe | 2011-07-12 14:41:00 +0200 | [diff] [blame] | 215 | CALL_SUBTEST_4( adjoint(MatrixXcf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE/2), internal::random<int>(1,EIGEN_TEST_MAX_SIZE/2))) ); |
| 216 | CALL_SUBTEST_5( adjoint(MatrixXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) ); |
| 217 | CALL_SUBTEST_6( adjoint(MatrixXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) ); |
Gael Guennebaud | c6eb84a | 2015-01-26 17:09:01 +0100 | [diff] [blame] | 218 | |
| 219 | // Complement for 128 bits vectorization: |
| 220 | CALL_SUBTEST_8( adjoint(Matrix2d()) ); |
| 221 | CALL_SUBTEST_9( adjoint(Matrix<int,4,4>()) ); |
| 222 | |
| 223 | // 256 bits vectorization: |
| 224 | CALL_SUBTEST_10( adjoint(Matrix<float,8,8>()) ); |
| 225 | CALL_SUBTEST_11( adjoint(Matrix<double,4,4>()) ); |
| 226 | CALL_SUBTEST_12( adjoint(Matrix<int,8,8>()) ); |
Gael Guennebaud | f5d2317 | 2009-07-13 21:14:47 +0200 | [diff] [blame] | 227 | } |
Gael Guennebaud | a8f66fe | 2011-07-12 14:41:00 +0200 | [diff] [blame] | 228 | // test a large static matrix only once |
Benoit Jacob | 2840ac7 | 2009-10-28 18:19:29 -0400 | [diff] [blame] | 229 | CALL_SUBTEST_7( adjoint(Matrix<float, 100, 100>()) ); |
Gael Guennebaud | b56bb44 | 2009-09-07 12:46:16 +0200 | [diff] [blame] | 230 | |
Gael Guennebaud | dff3a92 | 2018-07-17 15:52:58 +0200 | [diff] [blame] | 231 | CALL_SUBTEST_13( adjoint_extra<0>() ); |
Benoit Jacob | 2fdd067 | 2007-11-28 15:34:40 +0000 | [diff] [blame] | 232 | } |
Benoit Jacob | e05f291 | 2007-12-02 18:32:59 +0000 | [diff] [blame] | 233 | |