adjoint.cpp
Go to the documentation of this file.
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2006-2008 Benoit Jacob <jacob.benoit.1@gmail.com>
5 //
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/.
9 
10 #define EIGEN_NO_STATIC_ASSERT
11 
12 #include "main.h"
13 
14 template<bool IsInteger> struct adjoint_specific;
15 
16 template<> 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) {
19  VERIFY(test_isApproxWithRef((s1 * v1 + s2 * v2).dot(v3), numext::conj(s1) * v1.dot(v3) + numext::conj(s2) * v2.dot(v3), 0));
20  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 
27 template<> 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  using std::abs;
32 
33  RealScalar ref = NumTraits<Scalar>::IsInteger ? RealScalar(0) : (std::max)((s1 * v1 + s2 * v2).norm(),v3.norm());
34  VERIFY(test_isApproxWithRef((s1 * v1 + s2 * v2).dot(v3), numext::conj(s1) * v1.dot(v3) + numext::conj(s2) * v2.dot(v3), ref));
35  VERIFY(test_isApproxWithRef(v3.dot(s1 * v1 + s2 * v2), s1*v3.dot(v1)+s2*v3.dot(v2), ref));
36 
37  VERIFY_IS_APPROX(v1.squaredNorm(), v1.norm() * v1.norm());
38  // check normalized() and normalize()
39  VERIFY_IS_APPROX(v1, v1.norm() * v1.normalized());
40  v3 = v1;
41  v3.normalize();
42  VERIFY_IS_APPROX(v1, v1.norm() * v3);
43  VERIFY_IS_APPROX(v3, v1.normalized());
44  VERIFY_IS_APPROX(v3.norm(), RealScalar(1));
45 
46  // check null inputs
47  VERIFY_IS_APPROX((v1*0).normalized(), (v1*0));
48 #if (!EIGEN_ARCH_i386) || defined(EIGEN_VECTORIZE)
49  RealScalar very_small = (std::numeric_limits<RealScalar>::min)();
50  VERIFY( (v1*very_small).norm() == 0 );
51  VERIFY_IS_APPROX((v1*very_small).normalized(), (v1*very_small));
52  v3 = v1*very_small;
53  v3.normalize();
54  VERIFY_IS_APPROX(v3, (v1*very_small));
55 #endif
56 
57  // check compatibility of dot and adjoint
58  ref = NumTraits<Scalar>::IsInteger ? 0 : (std::max)((std::max)(v1.norm(),v2.norm()),(std::max)((square * v2).norm(),(square.adjoint() * v1).norm()));
59  VERIFY(internal::isMuchSmallerThan(abs(v1.dot(square * v2) - (square.adjoint() * v1).dot(v2)), ref, test_precision<Scalar>()));
60 
61  // check that Random().normalized() works: tricky as the random xpr must be evaluated by
62  // normalized() in order to produce a consistent result.
63  VERIFY_IS_APPROX(Vec::Random(v1.size()).normalized().norm(), RealScalar(1));
64  }
65 };
66 
67 template<typename MatrixType> void adjoint(const MatrixType& m)
68 {
69  /* this test covers the following files:
70  Transpose.h Conjugate.h Dot.h
71  */
72  using std::abs;
73  typedef typename MatrixType::Scalar Scalar;
74  typedef typename NumTraits<Scalar>::Real RealScalar;
78 
79  Index rows = m.rows();
80  Index cols = m.cols();
81 
82  MatrixType m1 = MatrixType::Random(rows, cols),
83  m2 = MatrixType::Random(rows, cols),
84  m3(rows, cols),
85  square = SquareMatrixType::Random(rows, rows);
86  VectorType v1 = VectorType::Random(rows),
87  v2 = VectorType::Random(rows),
88  v3 = VectorType::Random(rows),
89  vzero = VectorType::Zero(rows);
90 
91  Scalar s1 = internal::random<Scalar>(),
92  s2 = internal::random<Scalar>();
93 
94  // check basic compatibility of adjoint, transpose, conjugate
95  VERIFY_IS_APPROX(m1.transpose().conjugate().adjoint(), m1);
96  VERIFY_IS_APPROX(m1.adjoint().conjugate().transpose(), m1);
97 
98  // check multiplicative behavior
99  VERIFY_IS_APPROX((m1.adjoint() * m2).adjoint(), m2.adjoint() * m1);
100  VERIFY_IS_APPROX((s1 * m1).adjoint(), numext::conj(s1) * m1.adjoint());
101 
102  // check basic properties of dot, squaredNorm
103  VERIFY_IS_APPROX(numext::conj(v1.dot(v2)), v2.dot(v1));
104  VERIFY_IS_APPROX(numext::real(v1.dot(v1)), v1.squaredNorm());
105 
107 
108  VERIFY_IS_MUCH_SMALLER_THAN(abs(vzero.dot(v1)), static_cast<RealScalar>(1));
109 
110  // like in testBasicStuff, test operator() to check const-qualification
111  Index r = internal::random<Index>(0, rows-1),
112  c = internal::random<Index>(0, cols-1);
113  VERIFY_IS_APPROX(m1.conjugate()(r,c), numext::conj(m1(r,c)));
114  VERIFY_IS_APPROX(m1.adjoint()(c,r), numext::conj(m1(r,c)));
115 
116  // check inplace transpose
117  m3 = m1;
118  m3.transposeInPlace();
119  VERIFY_IS_APPROX(m3,m1.transpose());
120  m3.transposeInPlace();
121  VERIFY_IS_APPROX(m3,m1);
122 
123  if(PacketSize<m3.rows() && PacketSize<m3.cols())
124  {
125  m3 = m1;
126  Index i = internal::random<Index>(0,m3.rows()-PacketSize);
127  Index j = internal::random<Index>(0,m3.cols()-PacketSize);
128  m3.template block<PacketSize,PacketSize>(i,j).transposeInPlace();
129  VERIFY_IS_APPROX( (m3.template block<PacketSize,PacketSize>(i,j)), (m1.template block<PacketSize,PacketSize>(i,j).transpose()) );
130  m3.template block<PacketSize,PacketSize>(i,j).transposeInPlace();
131  VERIFY_IS_APPROX(m3,m1);
132  }
133 
134  // check inplace adjoint
135  m3 = m1;
136  m3.adjointInPlace();
137  VERIFY_IS_APPROX(m3,m1.adjoint());
138  m3.transposeInPlace();
139  VERIFY_IS_APPROX(m3,m1.conjugate());
140 
141  // check mixed dot product
143  RealVectorType rv1 = RealVectorType::Random(rows);
144  VERIFY_IS_APPROX(v1.dot(rv1.template cast<Scalar>()), v1.dot(rv1));
145  VERIFY_IS_APPROX(rv1.template cast<Scalar>().dot(v1), rv1.dot(v1));
146 }
147 
149 {
150  for(int i = 0; i < g_repeat; i++) {
151  CALL_SUBTEST_1( adjoint(Matrix<float, 1, 1>()) );
152  CALL_SUBTEST_2( adjoint(Matrix3d()) );
153  CALL_SUBTEST_3( adjoint(Matrix4f()) );
154 
155  CALL_SUBTEST_4( adjoint(MatrixXcf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE/2), internal::random<int>(1,EIGEN_TEST_MAX_SIZE/2))) );
156  CALL_SUBTEST_5( adjoint(MatrixXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
157  CALL_SUBTEST_6( adjoint(MatrixXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
158 
159  // Complement for 128 bits vectorization:
160  CALL_SUBTEST_8( adjoint(Matrix2d()) );
161  CALL_SUBTEST_9( adjoint(Matrix<int,4,4>()) );
162 
163  // 256 bits vectorization:
164  CALL_SUBTEST_10( adjoint(Matrix<float,8,8>()) );
165  CALL_SUBTEST_11( adjoint(Matrix<double,4,4>()) );
166  CALL_SUBTEST_12( adjoint(Matrix<int,8,8>()) );
167  }
168  // test a large static matrix only once
169  CALL_SUBTEST_7( adjoint(Matrix<float, 100, 100>()) );
170 
171 #ifdef EIGEN_TEST_PART_13
172  {
173  MatrixXcf a(10,10), b(10,10);
174  VERIFY_RAISES_ASSERT(a = a.transpose());
175  VERIFY_RAISES_ASSERT(a = a.transpose() + b);
176  VERIFY_RAISES_ASSERT(a = b + a.transpose());
177  VERIFY_RAISES_ASSERT(a = a.conjugate().transpose());
178  VERIFY_RAISES_ASSERT(a = a.adjoint());
179  VERIFY_RAISES_ASSERT(a = a.adjoint() + b);
180  VERIFY_RAISES_ASSERT(a = b + a.adjoint());
181 
182  // no assertion should be triggered for these cases:
183  a.transpose() = a.transpose();
184  a.transpose() += a.transpose();
185  a.transpose() += a.transpose() + b;
186  a.transpose() = a.adjoint();
187  a.transpose() += a.adjoint();
188  a.transpose() += a.adjoint() + b;
189 
190  // regression tests for check_for_aliasing
191  MatrixXd c(10,10);
192  c = 1.0 * MatrixXd::Ones(10,10) + c;
193  c = MatrixXd::Ones(10,10) * 1.0 + c;
194  c = c + MatrixXd::Ones(10,10) .cwiseProduct( MatrixXd::Zero(10,10) );
195  c = MatrixXd::Ones(10,10) * MatrixXd::Zero(10,10);
196  }
197 #endif
198 }
199 
Matrix3f m
SCALAR Scalar
Definition: bench_gemm.cpp:33
#define VERIFY_RAISES_ASSERT(a)
Definition: main.h:285
#define max(a, b)
Definition: datatypes.h:20
EIGEN_DEVICE_FUNC bool isMuchSmallerThan(const Scalar &x, const OtherScalar &y, const typename NumTraits< Scalar >::Real &precision=NumTraits< Scalar >::dummy_precision())
float real
Definition: datatypes.h:10
Vector v2
Scalar * b
Definition: benchVecAdd.cpp:17
void adjoint(const MatrixType &m)
Definition: adjoint.cpp:67
Vector v1
#define min(a, b)
Definition: datatypes.h:19
MatrixType m2(n_dims)
Scalar Scalar * c
Definition: benchVecAdd.cpp:17
MatrixXf MatrixType
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:150
Array33i a
Scalar Scalar int size
Definition: benchVecAdd.cpp:17
#define VERIFY_IS_APPROX(a, b)
Vector v3
Scalar EIGEN_BLAS_FUNC() dot(int *n, RealScalar *px, int *incx, RealScalar *py, int *incy)
static void run(const Vec &v1, const Vec &v2, Vec &v3, const Mat &square, Scalar s1, Scalar s2)
Definition: adjoint.cpp:29
bool test_isApproxWithRef(const Scalar &a, const Scalar &b, const ScalarRef &ref)
Definition: main.h:553
Matrix3d m1
Definition: IOFormat.cpp:2
static int g_repeat
Definition: main.h:144
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:33
static void run(const Vec &v1, const Vec &v2, Vec &v3, const Mat &square, Scalar s1, Scalar s2)
Definition: adjoint.cpp:18
NumTraits< Scalar >::Real RealScalar
Definition: bench_gemm.cpp:34
#define VERIFY_IS_MUCH_SMALLER_THAN(a, b)
Definition: main.h:335
#define VERIFY(a)
Definition: main.h:325
#define EIGEN_TEST_MAX_SIZE
The matrix class, also used for vectors and row-vectors.
void test_adjoint()
Definition: adjoint.cpp:148
void run(Expr &expr, Dev &dev)
Definition: TensorSyclRun.h:33
#define abs(x)
Definition: datatypes.h:17
std::ptrdiff_t j
EIGEN_DEVICE_FUNC const SquareReturnType square() const
ScalarWithExceptions conj(const ScalarWithExceptions &x)
Definition: exceptions.cpp:74


gtsam
Author(s):
autogenerated on Sat May 8 2021 02:41:36