• 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 
12 template<typename T>
get_random_size()13 int get_random_size()
14 {
15   const int factor = NumTraits<T>::ReadCost;
16   const int max_test_size = EIGEN_TEST_MAX_SIZE>2*factor ? EIGEN_TEST_MAX_SIZE/factor : EIGEN_TEST_MAX_SIZE;
17   return internal::random<int>(1,max_test_size);
18 }
19 
20 template<typename Scalar, int Mode, int TriOrder, int OtherOrder, int ResOrder, int OtherCols>
trmm(int rows=get_random_size<Scalar> (),int cols=get_random_size<Scalar> (),int otherCols=OtherCols==Dynamic?get_random_size<Scalar> ():OtherCols)21 void trmm(int rows=get_random_size<Scalar>(),
22           int cols=get_random_size<Scalar>(),
23           int otherCols = OtherCols==Dynamic?get_random_size<Scalar>():OtherCols)
24 {
25   typedef Matrix<Scalar,Dynamic,Dynamic,TriOrder> TriMatrix;
26   typedef Matrix<Scalar,Dynamic,OtherCols,OtherCols==1?ColMajor:OtherOrder> OnTheRight;
27   typedef Matrix<Scalar,OtherCols,Dynamic,OtherCols==1?RowMajor:OtherOrder> OnTheLeft;
28 
29   typedef Matrix<Scalar,Dynamic,OtherCols,OtherCols==1?ColMajor:ResOrder> ResXS;
30   typedef Matrix<Scalar,OtherCols,Dynamic,OtherCols==1?RowMajor:ResOrder> ResSX;
31 
32   TriMatrix  mat(rows,cols), tri(rows,cols), triTr(cols,rows);
33 
34   OnTheRight  ge_right(cols,otherCols);
35   OnTheLeft   ge_left(otherCols,rows);
36   ResSX       ge_sx, ge_sx_save;
37   ResXS       ge_xs, ge_xs_save;
38 
39   Scalar s1 = internal::random<Scalar>(),
40          s2 = internal::random<Scalar>();
41 
42   mat.setRandom();
43   tri = mat.template triangularView<Mode>();
44   triTr = mat.transpose().template triangularView<Mode>();
45   ge_right.setRandom();
46   ge_left.setRandom();
47 
48   VERIFY_IS_APPROX( ge_xs = mat.template triangularView<Mode>() * ge_right, tri * ge_right);
49   VERIFY_IS_APPROX( ge_sx = ge_left * mat.template triangularView<Mode>(), ge_left * tri);
50 
51   VERIFY_IS_APPROX( ge_xs.noalias() = mat.template triangularView<Mode>() * ge_right, tri * ge_right);
52   VERIFY_IS_APPROX( ge_sx.noalias() = ge_left * mat.template triangularView<Mode>(), ge_left * tri);
53 
54   VERIFY_IS_APPROX( ge_xs.noalias() = (s1*mat.adjoint()).template triangularView<Mode>() * (s2*ge_left.transpose()), s1*triTr.conjugate() * (s2*ge_left.transpose()));
55   VERIFY_IS_APPROX( ge_sx.noalias() = ge_right.transpose() * mat.adjoint().template triangularView<Mode>(), ge_right.transpose() * triTr.conjugate());
56 
57   VERIFY_IS_APPROX( ge_xs.noalias() = (s1*mat.adjoint()).template triangularView<Mode>() * (s2*ge_left.adjoint()), s1*triTr.conjugate() * (s2*ge_left.adjoint()));
58   VERIFY_IS_APPROX( ge_sx.noalias() = ge_right.adjoint() * mat.adjoint().template triangularView<Mode>(), ge_right.adjoint() * triTr.conjugate());
59 
60   ge_xs_save = ge_xs;
61   VERIFY_IS_APPROX( (ge_xs_save + s1*triTr.conjugate() * (s2*ge_left.adjoint())).eval(), ge_xs.noalias() += (s1*mat.adjoint()).template triangularView<Mode>() * (s2*ge_left.adjoint()) );
62   ge_sx.setRandom();
63   ge_sx_save = ge_sx;
64   VERIFY_IS_APPROX( ge_sx_save - (ge_right.adjoint() * (-s1 * triTr).conjugate()).eval(), ge_sx.noalias() -= (ge_right.adjoint() * (-s1 * mat).adjoint().template triangularView<Mode>()).eval());
65 
66   VERIFY_IS_APPROX( ge_xs = (s1*mat).adjoint().template triangularView<Mode>() * ge_left.adjoint(), numext::conj(s1) * triTr.conjugate() * ge_left.adjoint());
67 
68   // TODO check with sub-matrix expressions ?
69 }
70 
71 template<typename Scalar, int Mode, int TriOrder>
trmv(int rows=get_random_size<Scalar> (),int cols=get_random_size<Scalar> ())72 void trmv(int rows=get_random_size<Scalar>(), int cols=get_random_size<Scalar>())
73 {
74   trmm<Scalar,Mode,TriOrder,ColMajor,ColMajor,1>(rows,cols,1);
75 }
76 
77 template<typename Scalar, int Mode, int TriOrder, int OtherOrder, int ResOrder>
trmm(int rows=get_random_size<Scalar> (),int cols=get_random_size<Scalar> (),int otherCols=get_random_size<Scalar> ())78 void trmm(int rows=get_random_size<Scalar>(), int cols=get_random_size<Scalar>(), int otherCols = get_random_size<Scalar>())
79 {
80   trmm<Scalar,Mode,TriOrder,OtherOrder,ResOrder,Dynamic>(rows,cols,otherCols);
81 }
82 
83 #define CALL_ALL_ORDERS(NB,SCALAR,MODE)                                             \
84   EIGEN_CAT(CALL_SUBTEST_,NB)((trmm<SCALAR, MODE, ColMajor,ColMajor,ColMajor>()));  \
85   EIGEN_CAT(CALL_SUBTEST_,NB)((trmm<SCALAR, MODE, ColMajor,ColMajor,RowMajor>()));  \
86   EIGEN_CAT(CALL_SUBTEST_,NB)((trmm<SCALAR, MODE, ColMajor,RowMajor,ColMajor>()));  \
87   EIGEN_CAT(CALL_SUBTEST_,NB)((trmm<SCALAR, MODE, ColMajor,RowMajor,RowMajor>()));  \
88   EIGEN_CAT(CALL_SUBTEST_,NB)((trmm<SCALAR, MODE, RowMajor,ColMajor,ColMajor>()));  \
89   EIGEN_CAT(CALL_SUBTEST_,NB)((trmm<SCALAR, MODE, RowMajor,ColMajor,RowMajor>()));  \
90   EIGEN_CAT(CALL_SUBTEST_,NB)((trmm<SCALAR, MODE, RowMajor,RowMajor,ColMajor>()));  \
91   EIGEN_CAT(CALL_SUBTEST_,NB)((trmm<SCALAR, MODE, RowMajor,RowMajor,RowMajor>()));  \
92   \
93   EIGEN_CAT(CALL_SUBTEST_1,NB)((trmv<SCALAR, MODE, ColMajor>()));                   \
94   EIGEN_CAT(CALL_SUBTEST_1,NB)((trmv<SCALAR, MODE, RowMajor>()));
95 
96 
97 #define CALL_ALL(NB,SCALAR)                 \
98   CALL_ALL_ORDERS(EIGEN_CAT(1,NB),SCALAR,Upper)          \
99   CALL_ALL_ORDERS(EIGEN_CAT(2,NB),SCALAR,UnitUpper)      \
100   CALL_ALL_ORDERS(EIGEN_CAT(3,NB),SCALAR,StrictlyUpper)  \
101   CALL_ALL_ORDERS(EIGEN_CAT(1,NB),SCALAR,Lower)          \
102   CALL_ALL_ORDERS(EIGEN_CAT(2,NB),SCALAR,UnitLower)      \
103   CALL_ALL_ORDERS(EIGEN_CAT(3,NB),SCALAR,StrictlyLower)
104 
105 
test_product_trmm()106 void test_product_trmm()
107 {
108   for(int i = 0; i < g_repeat ; i++)
109   {
110     CALL_ALL(1,float);                //  EIGEN_SUFFIXES;11;111;21;121;31;131
111     CALL_ALL(2,double);               //  EIGEN_SUFFIXES;12;112;22;122;32;132
112     CALL_ALL(3,std::complex<float>);  //  EIGEN_SUFFIXES;13;113;23;123;33;133
113     CALL_ALL(4,std::complex<double>); //  EIGEN_SUFFIXES;14;114;24;124;34;134
114   }
115 }
116