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