• Home
  • Line#
  • Scopes#
  • Navigate#
  • Raw
  • Download
1 // This file is part 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 
product_selfadjoint(const MatrixType & m)12 template<typename MatrixType> void product_selfadjoint(const MatrixType& m)
13 {
14   typedef typename MatrixType::Index Index;
15   typedef typename MatrixType::Scalar Scalar;
16   typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
17   typedef Matrix<Scalar, 1, MatrixType::RowsAtCompileTime> RowVectorType;
18 
19   typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, Dynamic, RowMajor> RhsMatrixType;
20 
21   Index rows = m.rows();
22   Index cols = m.cols();
23 
24   MatrixType m1 = MatrixType::Random(rows, cols),
25              m2 = MatrixType::Random(rows, cols),
26              m3;
27   VectorType v1 = VectorType::Random(rows),
28              v2 = VectorType::Random(rows),
29              v3(rows);
30   RowVectorType r1 = RowVectorType::Random(rows),
31                 r2 = RowVectorType::Random(rows);
32   RhsMatrixType m4 = RhsMatrixType::Random(rows,10);
33 
34   Scalar s1 = internal::random<Scalar>(),
35          s2 = internal::random<Scalar>(),
36          s3 = internal::random<Scalar>();
37 
38   m1 = (m1.adjoint() + m1).eval();
39 
40   // rank2 update
41   m2 = m1.template triangularView<Lower>();
42   m2.template selfadjointView<Lower>().rankUpdate(v1,v2);
43   VERIFY_IS_APPROX(m2, (m1 + v1 * v2.adjoint()+ v2 * v1.adjoint()).template triangularView<Lower>().toDenseMatrix());
44 
45   m2 = m1.template triangularView<Upper>();
46   m2.template selfadjointView<Upper>().rankUpdate(-v1,s2*v2,s3);
47   VERIFY_IS_APPROX(m2, (m1 + (s3*(-v1)*(s2*v2).adjoint()+numext::conj(s3)*(s2*v2)*(-v1).adjoint())).template triangularView<Upper>().toDenseMatrix());
48 
49   m2 = m1.template triangularView<Upper>();
50   m2.template selfadjointView<Upper>().rankUpdate(-s2*r1.adjoint(),r2.adjoint()*s3,s1);
51   VERIFY_IS_APPROX(m2, (m1 + s1*(-s2*r1.adjoint())*(r2.adjoint()*s3).adjoint() + numext::conj(s1)*(r2.adjoint()*s3) * (-s2*r1.adjoint()).adjoint()).template triangularView<Upper>().toDenseMatrix());
52 
53   if (rows>1)
54   {
55     m2 = m1.template triangularView<Lower>();
56     m2.block(1,1,rows-1,cols-1).template selfadjointView<Lower>().rankUpdate(v1.tail(rows-1),v2.head(cols-1));
57     m3 = m1;
58     m3.block(1,1,rows-1,cols-1) += v1.tail(rows-1) * v2.head(cols-1).adjoint()+ v2.head(cols-1) * v1.tail(rows-1).adjoint();
59     VERIFY_IS_APPROX(m2, m3.template triangularView<Lower>().toDenseMatrix());
60   }
61 }
62 
test_product_selfadjoint()63 void test_product_selfadjoint()
64 {
65   int s = 0;
66   for(int i = 0; i < g_repeat ; i++) {
67     CALL_SUBTEST_1( product_selfadjoint(Matrix<float, 1, 1>()) );
68     CALL_SUBTEST_2( product_selfadjoint(Matrix<float, 2, 2>()) );
69     CALL_SUBTEST_3( product_selfadjoint(Matrix3d()) );
70 
71     s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/2);
72     CALL_SUBTEST_4( product_selfadjoint(MatrixXcf(s, s)) );
73     TEST_SET_BUT_UNUSED_VARIABLE(s)
74 
75     s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/2);
76     CALL_SUBTEST_5( product_selfadjoint(MatrixXcd(s,s)) );
77     TEST_SET_BUT_UNUSED_VARIABLE(s)
78 
79     s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE);
80     CALL_SUBTEST_6( product_selfadjoint(MatrixXd(s,s)) );
81     TEST_SET_BUT_UNUSED_VARIABLE(s)
82 
83     s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE);
84     CALL_SUBTEST_7( product_selfadjoint(Matrix<float,Dynamic,Dynamic,RowMajor>(s,s)) );
85     TEST_SET_BUT_UNUSED_VARIABLE(s)
86   }
87 }
88