• 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 
array_for_matrix(const MatrixType & m)12 template<typename MatrixType> void array_for_matrix(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> ColVectorType;
18   typedef Matrix<Scalar, 1, MatrixType::ColsAtCompileTime> RowVectorType;
19 
20   Index rows = m.rows();
21   Index cols = m.cols();
22 
23   MatrixType m1 = MatrixType::Random(rows, cols),
24              m2 = MatrixType::Random(rows, cols),
25              m3(rows, cols);
26 
27   ColVectorType cv1 = ColVectorType::Random(rows);
28   RowVectorType rv1 = RowVectorType::Random(cols);
29 
30   Scalar  s1 = internal::random<Scalar>(),
31           s2 = internal::random<Scalar>();
32 
33   // scalar addition
34   VERIFY_IS_APPROX(m1.array() + s1, s1 + m1.array());
35   VERIFY_IS_APPROX((m1.array() + s1).matrix(), MatrixType::Constant(rows,cols,s1) + m1);
36   VERIFY_IS_APPROX(((m1*Scalar(2)).array() - s2).matrix(), (m1+m1) - MatrixType::Constant(rows,cols,s2) );
37   m3 = m1;
38   m3.array() += s2;
39   VERIFY_IS_APPROX(m3, (m1.array() + s2).matrix());
40   m3 = m1;
41   m3.array() -= s1;
42   VERIFY_IS_APPROX(m3, (m1.array() - s1).matrix());
43 
44   // reductions
45   VERIFY_IS_MUCH_SMALLER_THAN(m1.colwise().sum().sum() - m1.sum(), m1.cwiseAbs().maxCoeff());
46   VERIFY_IS_MUCH_SMALLER_THAN(m1.rowwise().sum().sum() - m1.sum(), m1.cwiseAbs().maxCoeff());
47   VERIFY_IS_MUCH_SMALLER_THAN(m1.colwise().sum() + m2.colwise().sum() - (m1+m2).colwise().sum(), (m1+m2).cwiseAbs().maxCoeff());
48   VERIFY_IS_MUCH_SMALLER_THAN(m1.rowwise().sum() - m2.rowwise().sum() - (m1-m2).rowwise().sum(), (m1-m2).cwiseAbs().maxCoeff());
49   VERIFY_IS_APPROX(m1.colwise().sum(), m1.colwise().redux(internal::scalar_sum_op<Scalar>()));
50 
51   // vector-wise ops
52   m3 = m1;
53   VERIFY_IS_APPROX(m3.colwise() += cv1, m1.colwise() + cv1);
54   m3 = m1;
55   VERIFY_IS_APPROX(m3.colwise() -= cv1, m1.colwise() - cv1);
56   m3 = m1;
57   VERIFY_IS_APPROX(m3.rowwise() += rv1, m1.rowwise() + rv1);
58   m3 = m1;
59   VERIFY_IS_APPROX(m3.rowwise() -= rv1, m1.rowwise() - rv1);
60 
61   // empty objects
62   VERIFY_IS_APPROX(m1.block(0,0,0,cols).colwise().sum(),  RowVectorType::Zero(cols));
63   VERIFY_IS_APPROX(m1.block(0,0,rows,0).rowwise().prod(), ColVectorType::Ones(rows));
64 
65   // verify the const accessors exist
66   const Scalar& ref_m1 = m.matrix().array().coeffRef(0);
67   const Scalar& ref_m2 = m.matrix().array().coeffRef(0,0);
68   const Scalar& ref_a1 = m.array().matrix().coeffRef(0);
69   const Scalar& ref_a2 = m.array().matrix().coeffRef(0,0);
70   VERIFY(&ref_a1 == &ref_m1);
71   VERIFY(&ref_a2 == &ref_m2);
72 }
73 
comparisons(const MatrixType & m)74 template<typename MatrixType> void comparisons(const MatrixType& m)
75 {
76   typedef typename MatrixType::Index Index;
77   typedef typename MatrixType::Scalar Scalar;
78   typedef typename NumTraits<Scalar>::Real RealScalar;
79   typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
80 
81   Index rows = m.rows();
82   Index cols = m.cols();
83 
84   Index r = internal::random<Index>(0, rows-1),
85         c = internal::random<Index>(0, cols-1);
86 
87   MatrixType m1 = MatrixType::Random(rows, cols),
88              m2 = MatrixType::Random(rows, cols),
89              m3(rows, cols);
90 
91   VERIFY(((m1.array() + Scalar(1)) > m1.array()).all());
92   VERIFY(((m1.array() - Scalar(1)) < m1.array()).all());
93   if (rows*cols>1)
94   {
95     m3 = m1;
96     m3(r,c) += 1;
97     VERIFY(! (m1.array() < m3.array()).all() );
98     VERIFY(! (m1.array() > m3.array()).all() );
99   }
100 
101   // comparisons to scalar
102   VERIFY( (m1.array() != (m1(r,c)+1) ).any() );
103   VERIFY( (m1.array() > (m1(r,c)-1) ).any() );
104   VERIFY( (m1.array() < (m1(r,c)+1) ).any() );
105   VERIFY( (m1.array() == m1(r,c) ).any() );
106 
107   // test Select
108   VERIFY_IS_APPROX( (m1.array()<m2.array()).select(m1,m2), m1.cwiseMin(m2) );
109   VERIFY_IS_APPROX( (m1.array()>m2.array()).select(m1,m2), m1.cwiseMax(m2) );
110   Scalar mid = (m1.cwiseAbs().minCoeff() + m1.cwiseAbs().maxCoeff())/Scalar(2);
111   for (int j=0; j<cols; ++j)
112   for (int i=0; i<rows; ++i)
113     m3(i,j) = internal::abs(m1(i,j))<mid ? 0 : m1(i,j);
114   VERIFY_IS_APPROX( (m1.array().abs()<MatrixType::Constant(rows,cols,mid).array())
115                         .select(MatrixType::Zero(rows,cols),m1), m3);
116   // shorter versions:
117   VERIFY_IS_APPROX( (m1.array().abs()<MatrixType::Constant(rows,cols,mid).array())
118                         .select(0,m1), m3);
119   VERIFY_IS_APPROX( (m1.array().abs()>=MatrixType::Constant(rows,cols,mid).array())
120                         .select(m1,0), m3);
121   // even shorter version:
122   VERIFY_IS_APPROX( (m1.array().abs()<mid).select(0,m1), m3);
123 
124   // count
125   VERIFY(((m1.array().abs()+1)>RealScalar(0.1)).count() == rows*cols);
126 
127   typedef Matrix<typename MatrixType::Index, Dynamic, 1> VectorOfIndices;
128 
129   // TODO allows colwise/rowwise for array
130   VERIFY_IS_APPROX(((m1.array().abs()+1)>RealScalar(0.1)).matrix().colwise().count(), VectorOfIndices::Constant(cols,rows).transpose());
131   VERIFY_IS_APPROX(((m1.array().abs()+1)>RealScalar(0.1)).matrix().rowwise().count(), VectorOfIndices::Constant(rows, cols));
132 }
133 
lpNorm(const VectorType & v)134 template<typename VectorType> void lpNorm(const VectorType& v)
135 {
136   VectorType u = VectorType::Random(v.size());
137 
138   VERIFY_IS_APPROX(u.template lpNorm<Infinity>(), u.cwiseAbs().maxCoeff());
139   VERIFY_IS_APPROX(u.template lpNorm<1>(), u.cwiseAbs().sum());
140   VERIFY_IS_APPROX(u.template lpNorm<2>(), internal::sqrt(u.array().abs().square().sum()));
141   VERIFY_IS_APPROX(internal::pow(u.template lpNorm<5>(), typename VectorType::RealScalar(5)), u.array().abs().pow(5).sum());
142 }
143 
cwise_min_max(const MatrixType & m)144 template<typename MatrixType> void cwise_min_max(const MatrixType& m)
145 {
146   typedef typename MatrixType::Index Index;
147   typedef typename MatrixType::Scalar Scalar;
148 
149   Index rows = m.rows();
150   Index cols = m.cols();
151 
152   MatrixType m1 = MatrixType::Random(rows, cols);
153 
154   // min/max with array
155   Scalar maxM1 = m1.maxCoeff();
156   Scalar minM1 = m1.minCoeff();
157 
158   VERIFY_IS_APPROX(MatrixType::Constant(rows,cols, minM1), m1.cwiseMin(MatrixType::Constant(rows,cols, minM1)));
159   VERIFY_IS_APPROX(m1, m1.cwiseMin(MatrixType::Constant(rows,cols, maxM1)));
160 
161   VERIFY_IS_APPROX(MatrixType::Constant(rows,cols, maxM1), m1.cwiseMax(MatrixType::Constant(rows,cols, maxM1)));
162   VERIFY_IS_APPROX(m1, m1.cwiseMax(MatrixType::Constant(rows,cols, minM1)));
163 
164   // min/max with scalar input
165   VERIFY_IS_APPROX(MatrixType::Constant(rows,cols, minM1), m1.cwiseMin( minM1));
166   VERIFY_IS_APPROX(m1, m1.cwiseMin( maxM1));
167 
168   VERIFY_IS_APPROX(MatrixType::Constant(rows,cols, maxM1), m1.cwiseMax( maxM1));
169   VERIFY_IS_APPROX(m1, m1.cwiseMax( minM1));
170 
171 }
172 
test_array_for_matrix()173 void test_array_for_matrix()
174 {
175   for(int i = 0; i < g_repeat; i++) {
176     CALL_SUBTEST_1( array_for_matrix(Matrix<float, 1, 1>()) );
177     CALL_SUBTEST_2( array_for_matrix(Matrix2f()) );
178     CALL_SUBTEST_3( array_for_matrix(Matrix4d()) );
179     CALL_SUBTEST_4( array_for_matrix(MatrixXcf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
180     CALL_SUBTEST_5( array_for_matrix(MatrixXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
181     CALL_SUBTEST_6( array_for_matrix(MatrixXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
182   }
183   for(int i = 0; i < g_repeat; i++) {
184     CALL_SUBTEST_1( comparisons(Matrix<float, 1, 1>()) );
185     CALL_SUBTEST_2( comparisons(Matrix2f()) );
186     CALL_SUBTEST_3( comparisons(Matrix4d()) );
187     CALL_SUBTEST_5( comparisons(MatrixXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
188     CALL_SUBTEST_6( comparisons(MatrixXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
189   }
190   for(int i = 0; i < g_repeat; i++) {
191     CALL_SUBTEST_1( cwise_min_max(Matrix<float, 1, 1>()) );
192     CALL_SUBTEST_2( cwise_min_max(Matrix2f()) );
193     CALL_SUBTEST_3( cwise_min_max(Matrix4d()) );
194     CALL_SUBTEST_5( cwise_min_max(MatrixXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
195     CALL_SUBTEST_6( cwise_min_max(MatrixXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
196   }
197   for(int i = 0; i < g_repeat; i++) {
198     CALL_SUBTEST_1( lpNorm(Matrix<float, 1, 1>()) );
199     CALL_SUBTEST_2( lpNorm(Vector2f()) );
200     CALL_SUBTEST_7( lpNorm(Vector3d()) );
201     CALL_SUBTEST_8( lpNorm(Vector4f()) );
202     CALL_SUBTEST_5( lpNorm(VectorXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
203     CALL_SUBTEST_4( lpNorm(VectorXcf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
204   }
205 }
206