• Home
  • Line#
  • Scopes#
  • Navigate#
  • Raw
  • Download
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra. Eigen itself is part of the KDE project.
3 //
4 // Copyright (C) 2008 Gael Guennebaud <g.gael@free.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 #include <Eigen/Array>
12 
array(const MatrixType & m)13 template<typename MatrixType> void array(const MatrixType& m)
14 {
15   /* this test covers the following files:
16      Array.cpp
17   */
18 
19   typedef typename MatrixType::Scalar Scalar;
20   typedef typename NumTraits<Scalar>::Real RealScalar;
21   typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
22 
23   int rows = m.rows();
24   int cols = m.cols();
25 
26   MatrixType m1 = MatrixType::Random(rows, cols),
27              m2 = MatrixType::Random(rows, cols),
28              m3(rows, cols);
29 
30   Scalar  s1 = ei_random<Scalar>(),
31           s2 = ei_random<Scalar>();
32 
33   // scalar addition
34   VERIFY_IS_APPROX(m1.cwise() + s1, s1 + m1.cwise());
35   VERIFY_IS_APPROX(m1.cwise() + s1, MatrixType::Constant(rows,cols,s1) + m1);
36   VERIFY_IS_APPROX((m1*Scalar(2)).cwise() - s2, (m1+m1) - MatrixType::Constant(rows,cols,s2) );
37   m3 = m1;
38   m3.cwise() += s2;
39   VERIFY_IS_APPROX(m3, m1.cwise() + s2);
40   m3 = m1;
41   m3.cwise() -= s1;
42   VERIFY_IS_APPROX(m3, m1.cwise() - s1);
43 
44   // reductions
45   VERIFY_IS_APPROX(m1.colwise().sum().sum(), m1.sum());
46   VERIFY_IS_APPROX(m1.rowwise().sum().sum(), m1.sum());
47   if (!ei_isApprox(m1.sum(), (m1+m2).sum()))
48     VERIFY_IS_NOT_APPROX(((m1+m2).rowwise().sum()).sum(), m1.sum());
49   VERIFY_IS_APPROX(m1.colwise().sum(), m1.colwise().redux(internal::scalar_sum_op<Scalar>()));
50 }
51 
comparisons(const MatrixType & m)52 template<typename MatrixType> void comparisons(const MatrixType& m)
53 {
54   typedef typename MatrixType::Scalar Scalar;
55   typedef typename NumTraits<Scalar>::Real RealScalar;
56   typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
57 
58   int rows = m.rows();
59   int cols = m.cols();
60 
61   int r = ei_random<int>(0, rows-1),
62       c = ei_random<int>(0, cols-1);
63 
64   MatrixType m1 = MatrixType::Random(rows, cols),
65              m2 = MatrixType::Random(rows, cols),
66              m3(rows, cols);
67 
68   VERIFY(((m1.cwise() + Scalar(1)).cwise() > m1).all());
69   VERIFY(((m1.cwise() - Scalar(1)).cwise() < m1).all());
70   if (rows*cols>1)
71   {
72     m3 = m1;
73     m3(r,c) += 1;
74     VERIFY(! (m1.cwise() < m3).all() );
75     VERIFY(! (m1.cwise() > m3).all() );
76   }
77 
78   // comparisons to scalar
79   VERIFY( (m1.cwise() != (m1(r,c)+1) ).any() );
80   VERIFY( (m1.cwise() > (m1(r,c)-1) ).any() );
81   VERIFY( (m1.cwise() < (m1(r,c)+1) ).any() );
82   VERIFY( (m1.cwise() == m1(r,c) ).any() );
83 
84   // test Select
85   VERIFY_IS_APPROX( (m1.cwise()<m2).select(m1,m2), m1.cwise().min(m2) );
86   VERIFY_IS_APPROX( (m1.cwise()>m2).select(m1,m2), m1.cwise().max(m2) );
87   Scalar mid = (m1.cwise().abs().minCoeff() + m1.cwise().abs().maxCoeff())/Scalar(2);
88   for (int j=0; j<cols; ++j)
89   for (int i=0; i<rows; ++i)
90     m3(i,j) = ei_abs(m1(i,j))<mid ? 0 : m1(i,j);
91   VERIFY_IS_APPROX( (m1.cwise().abs().cwise()<MatrixType::Constant(rows,cols,mid))
92                         .select(MatrixType::Zero(rows,cols),m1), m3);
93   // shorter versions:
94   VERIFY_IS_APPROX( (m1.cwise().abs().cwise()<MatrixType::Constant(rows,cols,mid))
95                         .select(0,m1), m3);
96   VERIFY_IS_APPROX( (m1.cwise().abs().cwise()>=MatrixType::Constant(rows,cols,mid))
97                         .select(m1,0), m3);
98   // even shorter version:
99   VERIFY_IS_APPROX( (m1.cwise().abs().cwise()<mid).select(0,m1), m3);
100 
101   // count
102   VERIFY(((m1.cwise().abs().cwise()+1).cwise()>RealScalar(0.1)).count() == rows*cols);
103   VERIFY_IS_APPROX(((m1.cwise().abs().cwise()+1).cwise()>RealScalar(0.1)).colwise().count().template cast<int>(), RowVectorXi::Constant(cols,rows));
104   VERIFY_IS_APPROX(((m1.cwise().abs().cwise()+1).cwise()>RealScalar(0.1)).rowwise().count().template cast<int>(), VectorXi::Constant(rows, cols));
105 }
106 
lpNorm(const VectorType & v)107 template<typename VectorType> void lpNorm(const VectorType& v)
108 {
109   VectorType u = VectorType::Random(v.size());
110 
111   VERIFY_IS_APPROX(u.template lpNorm<Infinity>(), u.cwise().abs().maxCoeff());
112   VERIFY_IS_APPROX(u.template lpNorm<1>(), u.cwise().abs().sum());
113   VERIFY_IS_APPROX(u.template lpNorm<2>(), ei_sqrt(u.cwise().abs().cwise().square().sum()));
114   VERIFY_IS_APPROX(ei_pow(u.template lpNorm<5>(), typename VectorType::RealScalar(5)), u.cwise().abs().cwise().pow(5).sum());
115 }
116 
test_eigen2_array()117 void test_eigen2_array()
118 {
119   for(int i = 0; i < g_repeat; i++) {
120     CALL_SUBTEST_1( array(Matrix<float, 1, 1>()) );
121     CALL_SUBTEST_2( array(Matrix2f()) );
122     CALL_SUBTEST_3( array(Matrix4d()) );
123     CALL_SUBTEST_4( array(MatrixXcf(3, 3)) );
124     CALL_SUBTEST_5( array(MatrixXf(8, 12)) );
125     CALL_SUBTEST_6( array(MatrixXi(8, 12)) );
126   }
127   for(int i = 0; i < g_repeat; i++) {
128     CALL_SUBTEST_1( comparisons(Matrix<float, 1, 1>()) );
129     CALL_SUBTEST_2( comparisons(Matrix2f()) );
130     CALL_SUBTEST_3( comparisons(Matrix4d()) );
131     CALL_SUBTEST_5( comparisons(MatrixXf(8, 12)) );
132     CALL_SUBTEST_6( comparisons(MatrixXi(8, 12)) );
133   }
134   for(int i = 0; i < g_repeat; i++) {
135     CALL_SUBTEST_1( lpNorm(Matrix<float, 1, 1>()) );
136     CALL_SUBTEST_2( lpNorm(Vector2f()) );
137     CALL_SUBTEST_3( lpNorm(Vector3d()) );
138     CALL_SUBTEST_4( lpNorm(Vector4f()) );
139     CALL_SUBTEST_5( lpNorm(VectorXf(16)) );
140     CALL_SUBTEST_7( lpNorm(VectorXcd(10)) );
141   }
142 }
143