product_mmtr.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) 2010-2017 Gael Guennebaud <gael.guennebaud@inria.fr>
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 #include "main.h"
11 
12 #define CHECK_MMTR(DEST, TRI, OP) { \
13  ref3 = DEST; \
14  ref2 = ref1 = DEST; \
15  DEST.template triangularView<TRI>() OP; \
16  ref1 OP; \
17  ref2.template triangularView<TRI>() \
18  = ref1.template triangularView<TRI>(); \
19  VERIFY_IS_APPROX(DEST,ref2); \
20  \
21  DEST = ref3; \
22  ref3 = ref2; \
23  ref3.diagonal() = DEST.diagonal(); \
24  DEST.template triangularView<TRI|ZeroDiag>() OP; \
25  VERIFY_IS_APPROX(DEST,ref3); \
26  }
27 
28 template<typename Scalar> void mmtr(int size)
29 {
30  typedef Matrix<Scalar,Dynamic,Dynamic,ColMajor> MatrixColMaj;
31  typedef Matrix<Scalar,Dynamic,Dynamic,RowMajor> MatrixRowMaj;
32 
33  DenseIndex othersize = internal::random<DenseIndex>(1,200);
34 
35  MatrixColMaj matc = MatrixColMaj::Zero(size, size);
36  MatrixRowMaj matr = MatrixRowMaj::Zero(size, size);
37  MatrixColMaj ref1(size, size), ref2(size, size), ref3(size,size);
38 
39  MatrixColMaj soc(size,othersize); soc.setRandom();
40  MatrixColMaj osc(othersize,size); osc.setRandom();
41  MatrixRowMaj sor(size,othersize); sor.setRandom();
42  MatrixRowMaj osr(othersize,size); osr.setRandom();
43  MatrixColMaj sqc(size,size); sqc.setRandom();
44  MatrixRowMaj sqr(size,size); sqr.setRandom();
45 
46  Scalar s = internal::random<Scalar>();
47 
48  CHECK_MMTR(matc, Lower, = s*soc*sor.adjoint());
49  CHECK_MMTR(matc, Upper, = s*(soc*soc.adjoint()));
50  CHECK_MMTR(matr, Lower, = s*soc*soc.adjoint());
51  CHECK_MMTR(matr, Upper, = soc*(s*sor.adjoint()));
52 
53  CHECK_MMTR(matc, Lower, += s*soc*soc.adjoint());
54  CHECK_MMTR(matc, Upper, += s*(soc*sor.transpose()));
55  CHECK_MMTR(matr, Lower, += s*sor*soc.adjoint());
56  CHECK_MMTR(matr, Upper, += soc*(s*soc.adjoint()));
57 
58  CHECK_MMTR(matc, Lower, -= s*soc*soc.adjoint());
59  CHECK_MMTR(matc, Upper, -= s*(osc.transpose()*osc.conjugate()));
60  CHECK_MMTR(matr, Lower, -= s*soc*soc.adjoint());
61  CHECK_MMTR(matr, Upper, -= soc*(s*soc.adjoint()));
62 
63  CHECK_MMTR(matc, Lower, -= s*sqr*sqc.template triangularView<Upper>());
64  CHECK_MMTR(matc, Upper, = s*sqc*sqr.template triangularView<Upper>());
65  CHECK_MMTR(matc, Lower, += s*sqr*sqc.template triangularView<Lower>());
66  CHECK_MMTR(matc, Upper, = s*sqc*sqc.template triangularView<Lower>());
67 
68  CHECK_MMTR(matc, Lower, = (s*sqr).template triangularView<Upper>()*sqc);
69  CHECK_MMTR(matc, Upper, -= (s*sqc).template triangularView<Upper>()*sqc);
70  CHECK_MMTR(matc, Lower, = (s*sqr).template triangularView<Lower>()*sqc);
71  CHECK_MMTR(matc, Upper, += (s*sqc).template triangularView<Lower>()*sqc);
72 
73  // check aliasing
74  ref2 = ref1 = matc;
75  ref1 = sqc.adjoint() * matc * sqc;
76  ref2.template triangularView<Upper>() = ref1.template triangularView<Upper>();
77  matc.template triangularView<Upper>() = sqc.adjoint() * matc * sqc;
78  VERIFY_IS_APPROX(matc, ref2);
79 
80  ref2 = ref1 = matc;
81  ref1 = sqc * matc * sqc.adjoint();
82  ref2.template triangularView<Lower>() = ref1.template triangularView<Lower>();
83  matc.template triangularView<Lower>() = sqc * matc * sqc.adjoint();
84  VERIFY_IS_APPROX(matc, ref2);
85 
86  // destination with a non-default inner-stride
87  // see bug 1741
88  {
89  typedef Matrix<Scalar,Dynamic,Dynamic> MatrixX;
90  MatrixX buffer(2*size,2*size);
92  buffer.setZero();
93  CHECK_MMTR(map1, Lower, = s*soc*sor.adjoint());
94  }
95 }
96 
97 EIGEN_DECLARE_TEST(product_mmtr)
98 {
99  for(int i = 0; i < g_repeat ; i++)
100  {
101  CALL_SUBTEST_1((mmtr<float>(internal::random<int>(1,EIGEN_TEST_MAX_SIZE))));
102  CALL_SUBTEST_2((mmtr<double>(internal::random<int>(1,EIGEN_TEST_MAX_SIZE))));
103  CALL_SUBTEST_3((mmtr<std::complex<float> >(internal::random<int>(1,EIGEN_TEST_MAX_SIZE/2))));
104  CALL_SUBTEST_4((mmtr<std::complex<double> >(internal::random<int>(1,EIGEN_TEST_MAX_SIZE/2))));
105  }
106 }
Eigen::Stride
Holds strides information for Map.
Definition: Stride.h:48
s
RealScalar s
Definition: level1_cplx_impl.h:126
Eigen::Upper
@ Upper
Definition: Constants.h:211
buffer
Definition: pytypes.h:2270
size
Scalar Scalar int size
Definition: benchVecAdd.cpp:17
CALL_SUBTEST_4
#define CALL_SUBTEST_4(FUNC)
Definition: split_test_helper.h:22
CALL_SUBTEST_3
#define CALL_SUBTEST_3(FUNC)
Definition: split_test_helper.h:16
CALL_SUBTEST_1
#define CALL_SUBTEST_1(FUNC)
Definition: split_test_helper.h:4
EIGEN_DECLARE_TEST
EIGEN_DECLARE_TEST(product_mmtr)
Definition: product_mmtr.cpp:97
Eigen::g_repeat
static int g_repeat
Definition: main.h:169
Eigen::Lower
@ Lower
Definition: Constants.h:209
CALL_SUBTEST_2
#define CALL_SUBTEST_2(FUNC)
Definition: split_test_helper.h:10
Eigen::Map
A matrix or vector expression mapping an existing array of data.
Definition: Map.h:94
VERIFY_IS_APPROX
#define VERIFY_IS_APPROX(a, b)
Definition: integer_types.cpp:15
main.h
Eigen::DenseIndex
EIGEN_DEFAULT_DENSE_INDEX_TYPE DenseIndex
Definition: Meta.h:66
EIGEN_TEST_MAX_SIZE
#define EIGEN_TEST_MAX_SIZE
Definition: boostmultiprec.cpp:16
CHECK_MMTR
#define CHECK_MMTR(DEST, TRI, OP)
Definition: product_mmtr.cpp:12
Eigen::Matrix
The matrix class, also used for vectors and row-vectors.
Definition: 3rdparty/Eigen/Eigen/src/Core/Matrix.h:178
triangularView< Lower >
A triangularView< Lower >().adjoint().solveInPlace(B)
i
int i
Definition: BiCGSTAB_step_by_step.cpp:9
mmtr
void mmtr(int size)
Definition: product_mmtr.cpp:28
Scalar
SCALAR Scalar
Definition: bench_gemm.cpp:46


gtsam
Author(s):
autogenerated on Sat Nov 16 2024 04:03:42