• Home
  • Line#
  • Scopes#
  • Navigate#
  • Raw
  • Download
1 // This file is triangularView of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2008-2009 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 
trmv(const MatrixType & m)12 template<typename MatrixType> void trmv(const MatrixType& m)
13 {
14   typedef typename MatrixType::Index Index;
15   typedef typename MatrixType::Scalar Scalar;
16   typedef typename NumTraits<Scalar>::Real RealScalar;
17   typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
18 
19   RealScalar largerEps = 10*test_precision<RealScalar>();
20 
21   Index rows = m.rows();
22   Index cols = m.cols();
23 
24   MatrixType m1 = MatrixType::Random(rows, cols),
25              m3(rows, cols);
26   VectorType v1 = VectorType::Random(rows);
27 
28   Scalar s1 = internal::random<Scalar>();
29 
30   m1 = MatrixType::Random(rows, cols);
31 
32   // check with a column-major matrix
33   m3 = m1.template triangularView<Eigen::Lower>();
34   VERIFY((m3 * v1).isApprox(m1.template triangularView<Eigen::Lower>() * v1, largerEps));
35   m3 = m1.template triangularView<Eigen::Upper>();
36   VERIFY((m3 * v1).isApprox(m1.template triangularView<Eigen::Upper>() * v1, largerEps));
37   m3 = m1.template triangularView<Eigen::UnitLower>();
38   VERIFY((m3 * v1).isApprox(m1.template triangularView<Eigen::UnitLower>() * v1, largerEps));
39   m3 = m1.template triangularView<Eigen::UnitUpper>();
40   VERIFY((m3 * v1).isApprox(m1.template triangularView<Eigen::UnitUpper>() * v1, largerEps));
41 
42   // check conjugated and scalar multiple expressions (col-major)
43   m3 = m1.template triangularView<Eigen::Lower>();
44   VERIFY(((s1*m3).conjugate() * v1).isApprox((s1*m1).conjugate().template triangularView<Eigen::Lower>() * v1, largerEps));
45   m3 = m1.template triangularView<Eigen::Upper>();
46   VERIFY((m3.conjugate() * v1.conjugate()).isApprox(m1.conjugate().template triangularView<Eigen::Upper>() * v1.conjugate(), largerEps));
47 
48   // check with a row-major matrix
49   m3 = m1.template triangularView<Eigen::Upper>();
50   VERIFY((m3.transpose() * v1).isApprox(m1.transpose().template triangularView<Eigen::Lower>() * v1, largerEps));
51   m3 = m1.template triangularView<Eigen::Lower>();
52   VERIFY((m3.transpose() * v1).isApprox(m1.transpose().template triangularView<Eigen::Upper>() * v1, largerEps));
53   m3 = m1.template triangularView<Eigen::UnitUpper>();
54   VERIFY((m3.transpose() * v1).isApprox(m1.transpose().template triangularView<Eigen::UnitLower>() * v1, largerEps));
55   m3 = m1.template triangularView<Eigen::UnitLower>();
56   VERIFY((m3.transpose() * v1).isApprox(m1.transpose().template triangularView<Eigen::UnitUpper>() * v1, largerEps));
57 
58   // check conjugated and scalar multiple expressions (row-major)
59   m3 = m1.template triangularView<Eigen::Upper>();
60   VERIFY((m3.adjoint() * v1).isApprox(m1.adjoint().template triangularView<Eigen::Lower>() * v1, largerEps));
61   m3 = m1.template triangularView<Eigen::Lower>();
62   VERIFY((m3.adjoint() * (s1*v1.conjugate())).isApprox(m1.adjoint().template triangularView<Eigen::Upper>() * (s1*v1.conjugate()), largerEps));
63   m3 = m1.template triangularView<Eigen::UnitUpper>();
64 
65   // check transposed cases:
66   m3 = m1.template triangularView<Eigen::Lower>();
67   VERIFY((v1.transpose() * m3).isApprox(v1.transpose() * m1.template triangularView<Eigen::Lower>(), largerEps));
68   VERIFY((v1.adjoint() * m3).isApprox(v1.adjoint() * m1.template triangularView<Eigen::Lower>(), largerEps));
69   VERIFY((v1.adjoint() * m3.adjoint()).isApprox(v1.adjoint() * m1.template triangularView<Eigen::Lower>().adjoint(), largerEps));
70 
71   // TODO check with sub-matrices
72 }
73 
test_product_trmv()74 void test_product_trmv()
75 {
76   int s = 0;
77   for(int i = 0; i < g_repeat ; i++) {
78     CALL_SUBTEST_1( trmv(Matrix<float, 1, 1>()) );
79     CALL_SUBTEST_2( trmv(Matrix<float, 2, 2>()) );
80     CALL_SUBTEST_3( trmv(Matrix3d()) );
81     s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/2);
82     CALL_SUBTEST_4( trmv(MatrixXcf(s,s)) );
83     s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/2);
84     CALL_SUBTEST_5( trmv(MatrixXcd(s,s)) );
85     s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE);
86     CALL_SUBTEST_6( trmv(Matrix<float,Dynamic,Dynamic,RowMajor>(s, s)) );
87   }
88   TEST_SET_BUT_UNUSED_VARIABLE(s);
89 }
90