• 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 Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> ColVectorType;
17    typedef Matrix<Scalar, 1, MatrixType::ColsAtCompileTime> RowVectorType;
18  
19    Index rows = m.rows();
20    Index cols = m.cols();
21  
22    MatrixType m1 = MatrixType::Random(rows, cols),
23               m2 = MatrixType::Random(rows, cols),
24               m3(rows, cols);
25  
26    ColVectorType cv1 = ColVectorType::Random(rows);
27    RowVectorType rv1 = RowVectorType::Random(cols);
28  
29    Scalar  s1 = internal::random<Scalar>(),
30            s2 = internal::random<Scalar>();
31  
32    // scalar addition
33    VERIFY_IS_APPROX(m1.array() + s1, s1 + m1.array());
34    VERIFY_IS_APPROX((m1.array() + s1).matrix(), MatrixType::Constant(rows,cols,s1) + m1);
35    VERIFY_IS_APPROX(((m1*Scalar(2)).array() - s2).matrix(), (m1+m1) - MatrixType::Constant(rows,cols,s2) );
36    m3 = m1;
37    m3.array() += s2;
38    VERIFY_IS_APPROX(m3, (m1.array() + s2).matrix());
39    m3 = m1;
40    m3.array() -= s1;
41    VERIFY_IS_APPROX(m3, (m1.array() - s1).matrix());
42  
43    // reductions
44    VERIFY_IS_MUCH_SMALLER_THAN(m1.colwise().sum().sum() - m1.sum(), m1.squaredNorm());
45    VERIFY_IS_MUCH_SMALLER_THAN(m1.rowwise().sum().sum() - m1.sum(), m1.squaredNorm());
46    VERIFY_IS_MUCH_SMALLER_THAN(m1.colwise().sum() + m2.colwise().sum() - (m1+m2).colwise().sum(), (m1+m2).squaredNorm());
47    VERIFY_IS_MUCH_SMALLER_THAN(m1.rowwise().sum() - m2.rowwise().sum() - (m1-m2).rowwise().sum(), (m1-m2).squaredNorm());
48    VERIFY_IS_APPROX(m1.colwise().sum(), m1.colwise().redux(internal::scalar_sum_op<Scalar,Scalar>()));
49  
50    // vector-wise ops
51    m3 = m1;
52    VERIFY_IS_APPROX(m3.colwise() += cv1, m1.colwise() + cv1);
53    m3 = m1;
54    VERIFY_IS_APPROX(m3.colwise() -= cv1, m1.colwise() - cv1);
55    m3 = m1;
56    VERIFY_IS_APPROX(m3.rowwise() += rv1, m1.rowwise() + rv1);
57    m3 = m1;
58    VERIFY_IS_APPROX(m3.rowwise() -= rv1, m1.rowwise() - rv1);
59  
60    // empty objects
61    VERIFY_IS_APPROX(m1.block(0,0,0,cols).colwise().sum(),  RowVectorType::Zero(cols));
62    VERIFY_IS_APPROX(m1.block(0,0,rows,0).rowwise().prod(), ColVectorType::Ones(rows));
63  
64    // verify the const accessors exist
65    const Scalar& ref_m1 = m.matrix().array().coeffRef(0);
66    const Scalar& ref_m2 = m.matrix().array().coeffRef(0,0);
67    const Scalar& ref_a1 = m.array().matrix().coeffRef(0);
68    const Scalar& ref_a2 = m.array().matrix().coeffRef(0,0);
69    VERIFY(&ref_a1 == &ref_m1);
70    VERIFY(&ref_a2 == &ref_m2);
71  
72    // Check write accessors:
73    m1.array().coeffRef(0,0) = 1;
74    VERIFY_IS_APPROX(m1(0,0),Scalar(1));
75    m1.array()(0,0) = 2;
76    VERIFY_IS_APPROX(m1(0,0),Scalar(2));
77    m1.array().matrix().coeffRef(0,0) = 3;
78    VERIFY_IS_APPROX(m1(0,0),Scalar(3));
79    m1.array().matrix()(0,0) = 4;
80    VERIFY_IS_APPROX(m1(0,0),Scalar(4));
81  }
82  
comparisons(const MatrixType & m)83  template<typename MatrixType> void comparisons(const MatrixType& m)
84  {
85    using std::abs;
86    typedef typename MatrixType::Index Index;
87    typedef typename MatrixType::Scalar Scalar;
88    typedef typename NumTraits<Scalar>::Real RealScalar;
89  
90    Index rows = m.rows();
91    Index cols = m.cols();
92  
93    Index r = internal::random<Index>(0, rows-1),
94          c = internal::random<Index>(0, cols-1);
95  
96    MatrixType m1 = MatrixType::Random(rows, cols),
97               m2 = MatrixType::Random(rows, cols),
98               m3(rows, cols);
99  
100    VERIFY(((m1.array() + Scalar(1)) > m1.array()).all());
101    VERIFY(((m1.array() - Scalar(1)) < m1.array()).all());
102    if (rows*cols>1)
103    {
104      m3 = m1;
105      m3(r,c) += 1;
106      VERIFY(! (m1.array() < m3.array()).all() );
107      VERIFY(! (m1.array() > m3.array()).all() );
108    }
109  
110    // comparisons to scalar
111    VERIFY( (m1.array() != (m1(r,c)+1) ).any() );
112    VERIFY( (m1.array() > (m1(r,c)-1) ).any() );
113    VERIFY( (m1.array() < (m1(r,c)+1) ).any() );
114    VERIFY( (m1.array() == m1(r,c) ).any() );
115    VERIFY( m1.cwiseEqual(m1(r,c)).any() );
116  
117    // test Select
118    VERIFY_IS_APPROX( (m1.array()<m2.array()).select(m1,m2), m1.cwiseMin(m2) );
119    VERIFY_IS_APPROX( (m1.array()>m2.array()).select(m1,m2), m1.cwiseMax(m2) );
120    Scalar mid = (m1.cwiseAbs().minCoeff() + m1.cwiseAbs().maxCoeff())/Scalar(2);
121    for (int j=0; j<cols; ++j)
122    for (int i=0; i<rows; ++i)
123      m3(i,j) = abs(m1(i,j))<mid ? 0 : m1(i,j);
124    VERIFY_IS_APPROX( (m1.array().abs()<MatrixType::Constant(rows,cols,mid).array())
125                          .select(MatrixType::Zero(rows,cols),m1), m3);
126    // shorter versions:
127    VERIFY_IS_APPROX( (m1.array().abs()<MatrixType::Constant(rows,cols,mid).array())
128                          .select(0,m1), m3);
129    VERIFY_IS_APPROX( (m1.array().abs()>=MatrixType::Constant(rows,cols,mid).array())
130                          .select(m1,0), m3);
131    // even shorter version:
132    VERIFY_IS_APPROX( (m1.array().abs()<mid).select(0,m1), m3);
133  
134    // count
135    VERIFY(((m1.array().abs()+1)>RealScalar(0.1)).count() == rows*cols);
136  
137    // and/or
138    VERIFY( ((m1.array()<RealScalar(0)).matrix() && (m1.array()>RealScalar(0)).matrix()).count() == 0);
139    VERIFY( ((m1.array()<RealScalar(0)).matrix() || (m1.array()>=RealScalar(0)).matrix()).count() == rows*cols);
140    RealScalar a = m1.cwiseAbs().mean();
141    VERIFY( ((m1.array()<-a).matrix() || (m1.array()>a).matrix()).count() == (m1.cwiseAbs().array()>a).count());
142  
143    typedef Matrix<typename MatrixType::Index, Dynamic, 1> VectorOfIndices;
144  
145    // TODO allows colwise/rowwise for array
146    VERIFY_IS_APPROX(((m1.array().abs()+1)>RealScalar(0.1)).matrix().colwise().count(), VectorOfIndices::Constant(cols,rows).transpose());
147    VERIFY_IS_APPROX(((m1.array().abs()+1)>RealScalar(0.1)).matrix().rowwise().count(), VectorOfIndices::Constant(rows, cols));
148  }
149  
lpNorm(const VectorType & v)150  template<typename VectorType> void lpNorm(const VectorType& v)
151  {
152    using std::sqrt;
153    typedef typename VectorType::RealScalar RealScalar;
154    VectorType u = VectorType::Random(v.size());
155  
156    if(v.size()==0)
157    {
158      VERIFY_IS_APPROX(u.template lpNorm<Infinity>(), RealScalar(0));
159      VERIFY_IS_APPROX(u.template lpNorm<1>(), RealScalar(0));
160      VERIFY_IS_APPROX(u.template lpNorm<2>(), RealScalar(0));
161      VERIFY_IS_APPROX(u.template lpNorm<5>(), RealScalar(0));
162    }
163    else
164    {
165      VERIFY_IS_APPROX(u.template lpNorm<Infinity>(), u.cwiseAbs().maxCoeff());
166    }
167  
168    VERIFY_IS_APPROX(u.template lpNorm<1>(), u.cwiseAbs().sum());
169    VERIFY_IS_APPROX(u.template lpNorm<2>(), sqrt(u.array().abs().square().sum()));
170    VERIFY_IS_APPROX(numext::pow(u.template lpNorm<5>(), typename VectorType::RealScalar(5)), u.array().abs().pow(5).sum());
171  }
172  
cwise_min_max(const MatrixType & m)173  template<typename MatrixType> void cwise_min_max(const MatrixType& m)
174  {
175    typedef typename MatrixType::Index Index;
176    typedef typename MatrixType::Scalar Scalar;
177  
178    Index rows = m.rows();
179    Index cols = m.cols();
180  
181    MatrixType m1 = MatrixType::Random(rows, cols);
182  
183    // min/max with array
184    Scalar maxM1 = m1.maxCoeff();
185    Scalar minM1 = m1.minCoeff();
186  
187    VERIFY_IS_APPROX(MatrixType::Constant(rows,cols, minM1), m1.cwiseMin(MatrixType::Constant(rows,cols, minM1)));
188    VERIFY_IS_APPROX(m1, m1.cwiseMin(MatrixType::Constant(rows,cols, maxM1)));
189  
190    VERIFY_IS_APPROX(MatrixType::Constant(rows,cols, maxM1), m1.cwiseMax(MatrixType::Constant(rows,cols, maxM1)));
191    VERIFY_IS_APPROX(m1, m1.cwiseMax(MatrixType::Constant(rows,cols, minM1)));
192  
193    // min/max with scalar input
194    VERIFY_IS_APPROX(MatrixType::Constant(rows,cols, minM1), m1.cwiseMin( minM1));
195    VERIFY_IS_APPROX(m1, m1.cwiseMin(maxM1));
196    VERIFY_IS_APPROX(-m1, (-m1).cwiseMin(-minM1));
197    VERIFY_IS_APPROX(-m1.array(), ((-m1).array().min)( -minM1));
198  
199    VERIFY_IS_APPROX(MatrixType::Constant(rows,cols, maxM1), m1.cwiseMax( maxM1));
200    VERIFY_IS_APPROX(m1, m1.cwiseMax(minM1));
201    VERIFY_IS_APPROX(-m1, (-m1).cwiseMax(-maxM1));
202    VERIFY_IS_APPROX(-m1.array(), ((-m1).array().max)(-maxM1));
203  
204    VERIFY_IS_APPROX(MatrixType::Constant(rows,cols, minM1).array(), (m1.array().min)( minM1));
205    VERIFY_IS_APPROX(m1.array(), (m1.array().min)( maxM1));
206  
207    VERIFY_IS_APPROX(MatrixType::Constant(rows,cols, maxM1).array(), (m1.array().max)( maxM1));
208    VERIFY_IS_APPROX(m1.array(), (m1.array().max)( minM1));
209  
210  }
211  
resize(const MatrixTraits & t)212  template<typename MatrixTraits> void resize(const MatrixTraits& t)
213  {
214    typedef typename MatrixTraits::Index Index;
215    typedef typename MatrixTraits::Scalar Scalar;
216    typedef Matrix<Scalar,Dynamic,Dynamic> MatrixType;
217    typedef Array<Scalar,Dynamic,Dynamic> Array2DType;
218    typedef Matrix<Scalar,Dynamic,1> VectorType;
219    typedef Array<Scalar,Dynamic,1> Array1DType;
220  
221    Index rows = t.rows(), cols = t.cols();
222  
223    MatrixType m(rows,cols);
224    VectorType v(rows);
225    Array2DType a2(rows,cols);
226    Array1DType a1(rows);
227  
228    m.array().resize(rows+1,cols+1);
229    VERIFY(m.rows()==rows+1 && m.cols()==cols+1);
230    a2.matrix().resize(rows+1,cols+1);
231    VERIFY(a2.rows()==rows+1 && a2.cols()==cols+1);
232    v.array().resize(cols);
233    VERIFY(v.size()==cols);
234    a1.matrix().resize(cols);
235    VERIFY(a1.size()==cols);
236  }
237  
238  template<int>
regression_bug_654()239  void regression_bug_654()
240  {
241    ArrayXf a = RowVectorXf(3);
242    VectorXf v = Array<float,1,Dynamic>(3);
243  }
244  
245  // Check propagation of LvalueBit through Array/Matrix-Wrapper
246  template<int>
regrrssion_bug_1410()247  void regrrssion_bug_1410()
248  {
249    const Matrix4i M;
250    const Array4i A;
251    ArrayWrapper<const Matrix4i> MA = M.array();
252    MA.row(0);
253    MatrixWrapper<const Array4i> AM = A.matrix();
254    AM.row(0);
255  
256    VERIFY((internal::traits<ArrayWrapper<const Matrix4i> >::Flags&LvalueBit)==0);
257    VERIFY((internal::traits<MatrixWrapper<const Array4i> >::Flags&LvalueBit)==0);
258  
259    VERIFY((internal::traits<ArrayWrapper<Matrix4i> >::Flags&LvalueBit)==LvalueBit);
260    VERIFY((internal::traits<MatrixWrapper<Array4i> >::Flags&LvalueBit)==LvalueBit);
261  }
262  
test_array_for_matrix()263  void test_array_for_matrix()
264  {
265    for(int i = 0; i < g_repeat; i++) {
266      CALL_SUBTEST_1( array_for_matrix(Matrix<float, 1, 1>()) );
267      CALL_SUBTEST_2( array_for_matrix(Matrix2f()) );
268      CALL_SUBTEST_3( array_for_matrix(Matrix4d()) );
269      CALL_SUBTEST_4( array_for_matrix(MatrixXcf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
270      CALL_SUBTEST_5( array_for_matrix(MatrixXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
271      CALL_SUBTEST_6( array_for_matrix(MatrixXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
272    }
273    for(int i = 0; i < g_repeat; i++) {
274      CALL_SUBTEST_1( comparisons(Matrix<float, 1, 1>()) );
275      CALL_SUBTEST_2( comparisons(Matrix2f()) );
276      CALL_SUBTEST_3( comparisons(Matrix4d()) );
277      CALL_SUBTEST_5( comparisons(MatrixXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
278      CALL_SUBTEST_6( comparisons(MatrixXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
279    }
280    for(int i = 0; i < g_repeat; i++) {
281      CALL_SUBTEST_1( cwise_min_max(Matrix<float, 1, 1>()) );
282      CALL_SUBTEST_2( cwise_min_max(Matrix2f()) );
283      CALL_SUBTEST_3( cwise_min_max(Matrix4d()) );
284      CALL_SUBTEST_5( cwise_min_max(MatrixXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
285      CALL_SUBTEST_6( cwise_min_max(MatrixXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
286    }
287    for(int i = 0; i < g_repeat; i++) {
288      CALL_SUBTEST_1( lpNorm(Matrix<float, 1, 1>()) );
289      CALL_SUBTEST_2( lpNorm(Vector2f()) );
290      CALL_SUBTEST_7( lpNorm(Vector3d()) );
291      CALL_SUBTEST_8( lpNorm(Vector4f()) );
292      CALL_SUBTEST_5( lpNorm(VectorXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
293      CALL_SUBTEST_4( lpNorm(VectorXcf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
294    }
295    CALL_SUBTEST_5( lpNorm(VectorXf(0)) );
296    CALL_SUBTEST_4( lpNorm(VectorXcf(0)) );
297    for(int i = 0; i < g_repeat; i++) {
298      CALL_SUBTEST_4( resize(MatrixXcf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
299      CALL_SUBTEST_5( resize(MatrixXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
300      CALL_SUBTEST_6( resize(MatrixXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
301    }
302    CALL_SUBTEST_6( regression_bug_654<0>() );
303    CALL_SUBTEST_6( regrrssion_bug_1410<0>() );
304  }
305