blob: 9b349db499c07d6249ea5453d455814b6b415486 [file] [log] [blame]
Gael Guennebaud8cef5412008-06-21 17:28:07 +00001// This file is part of Eigen, a lightweight C++ template library
Benoit Jacob6347b1d2009-05-22 20:25:33 +02002// for linear algebra.
Gael Guennebaud8cef5412008-06-21 17:28:07 +00003//
Gael Guennebaudec5c6082009-07-10 16:10:03 +02004// Copyright (C) 2008-2009 Gael Guennebaud <g.gael@free.fr>
Gael Guennebaud8cef5412008-06-21 17:28:07 +00005//
6// Eigen is free software; you can redistribute it and/or
7// modify it under the terms of the GNU Lesser General Public
8// License as published by the Free Software Foundation; either
9// version 3 of the License, or (at your option) any later version.
10//
11// Alternatively, you can redistribute it and/or
12// modify it under the terms of the GNU General Public License as
13// published by the Free Software Foundation; either version 2 of
14// the License, or (at your option) any later version.
15//
16// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
17// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
18// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
19// GNU General Public License for more details.
20//
21// You should have received a copy of the GNU Lesser General Public
22// License and a copy of the GNU General Public License along with
23// Eigen. If not, see <http://www.gnu.org/licenses/>.
24
Gael Guennebaude3d890b2009-11-18 18:15:19 +010025#define EIGEN2_SUPPORT
Gael Guennebaud8cef5412008-06-21 17:28:07 +000026#include "main.h"
27#include <Eigen/Array>
28
Gael Guennebaud59dc1da2008-09-03 17:16:28 +000029template<typename MatrixType> void array(const MatrixType& m)
Gael Guennebaud8cef5412008-06-21 17:28:07 +000030{
31 /* this test covers the following files:
32 Array.cpp
33 */
34
35 typedef typename MatrixType::Scalar Scalar;
Gael Guennebaud2120fed2008-08-23 15:14:20 +000036 typedef typename NumTraits<Scalar>::Real RealScalar;
Gael Guennebaud627595a2009-06-10 11:20:30 +020037 typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> ColVectorType;
38 typedef Matrix<Scalar, 1, MatrixType::ColsAtCompileTime> RowVectorType;
Gael Guennebaud8cef5412008-06-21 17:28:07 +000039
40 int rows = m.rows();
41 int cols = m.cols();
42
Benoit Jacob46fe7a32008-09-01 17:31:21 +000043 MatrixType m1 = MatrixType::Random(rows, cols),
44 m2 = MatrixType::Random(rows, cols),
Gael Guennebaud8cef5412008-06-21 17:28:07 +000045 m3(rows, cols);
46
Gael Guennebaud627595a2009-06-10 11:20:30 +020047 ColVectorType cv1 = ColVectorType::Random(rows);
48 RowVectorType rv1 = RowVectorType::Random(cols);
49
Benoit Jacob46fe7a32008-09-01 17:31:21 +000050 Scalar s1 = ei_random<Scalar>(),
51 s2 = ei_random<Scalar>();
Gael Guennebaud8cef5412008-06-21 17:28:07 +000052
Gael Guennebaud59dc1da2008-09-03 17:16:28 +000053 // scalar addition
Benoit Jacobf5791ee2008-07-08 00:49:10 +000054 VERIFY_IS_APPROX(m1.cwise() + s1, s1 + m1.cwise());
Gael Guennebaudc10f0692008-07-21 00:34:46 +000055 VERIFY_IS_APPROX(m1.cwise() + s1, MatrixType::Constant(rows,cols,s1) + m1);
56 VERIFY_IS_APPROX((m1*Scalar(2)).cwise() - s2, (m1+m1) - MatrixType::Constant(rows,cols,s2) );
Gael Guennebaud8cef5412008-06-21 17:28:07 +000057 m3 = m1;
Benoit Jacobf5791ee2008-07-08 00:49:10 +000058 m3.cwise() += s2;
59 VERIFY_IS_APPROX(m3, m1.cwise() + s2);
Gael Guennebaud8cef5412008-06-21 17:28:07 +000060 m3 = m1;
Benoit Jacobf5791ee2008-07-08 00:49:10 +000061 m3.cwise() -= s1;
62 VERIFY_IS_APPROX(m3, m1.cwise() - s1);
Gael Guennebaud05ad0832008-07-19 13:03:23 +000063
Gael Guennebaud59dc1da2008-09-03 17:16:28 +000064 // reductions
Gael Guennebaud05ad0832008-07-19 13:03:23 +000065 VERIFY_IS_APPROX(m1.colwise().sum().sum(), m1.sum());
66 VERIFY_IS_APPROX(m1.rowwise().sum().sum(), m1.sum());
Gael Guennebaud2120fed2008-08-23 15:14:20 +000067 if (!ei_isApprox(m1.sum(), (m1+m2).sum()))
68 VERIFY_IS_NOT_APPROX(((m1+m2).rowwise().sum()).sum(), m1.sum());
Gael Guennebaud05ad0832008-07-19 13:03:23 +000069 VERIFY_IS_APPROX(m1.colwise().sum(), m1.colwise().redux(ei_scalar_sum_op<Scalar>()));
Gael Guennebaud627595a2009-06-10 11:20:30 +020070
71 // vector-wise ops
72 m3 = m1;
73 VERIFY_IS_APPROX(m3.colwise() += cv1, m1.colwise() + cv1);
74 m3 = m1;
75 VERIFY_IS_APPROX(m3.colwise() -= cv1, m1.colwise() - cv1);
76 m3 = m1;
77 VERIFY_IS_APPROX(m3.rowwise() += rv1, m1.rowwise() + rv1);
78 m3 = m1;
79 VERIFY_IS_APPROX(m3.rowwise() -= rv1, m1.rowwise() - rv1);
Gael Guennebaud8cef5412008-06-21 17:28:07 +000080}
81
82template<typename MatrixType> void comparisons(const MatrixType& m)
83{
84 typedef typename MatrixType::Scalar Scalar;
Benoit Jacob874ff5a2009-01-26 17:56:04 +000085 typedef typename NumTraits<Scalar>::Real RealScalar;
Gael Guennebaud8cef5412008-06-21 17:28:07 +000086 typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
87
88 int rows = m.rows();
89 int cols = m.cols();
90
91 int r = ei_random<int>(0, rows-1),
92 c = ei_random<int>(0, cols-1);
93
Gael Guennebaudc10f0692008-07-21 00:34:46 +000094 MatrixType m1 = MatrixType::Random(rows, cols),
95 m2 = MatrixType::Random(rows, cols),
Gael Guennebaud8cef5412008-06-21 17:28:07 +000096 m3(rows, cols);
97
Benoit Jacobf5791ee2008-07-08 00:49:10 +000098 VERIFY(((m1.cwise() + Scalar(1)).cwise() > m1).all());
99 VERIFY(((m1.cwise() - Scalar(1)).cwise() < m1).all());
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000100 if (rows*cols>1)
101 {
102 m3 = m1;
103 m3(r,c) += 1;
Benoit Jacobf5791ee2008-07-08 00:49:10 +0000104 VERIFY(! (m1.cwise() < m3).all() );
105 VERIFY(! (m1.cwise() > m3).all() );
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000106 }
Gael Guennebaud627595a2009-06-10 11:20:30 +0200107
Gael Guennebaude14aa8c2008-09-03 17:56:06 +0000108 // comparisons to scalar
109 VERIFY( (m1.cwise() != (m1(r,c)+1) ).any() );
110 VERIFY( (m1.cwise() > (m1(r,c)-1) ).any() );
111 VERIFY( (m1.cwise() < (m1(r,c)+1) ).any() );
112 VERIFY( (m1.cwise() == m1(r,c) ).any() );
Gael Guennebaud627595a2009-06-10 11:20:30 +0200113
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000114 // test Select
115 VERIFY_IS_APPROX( (m1.cwise()<m2).select(m1,m2), m1.cwise().min(m2) );
116 VERIFY_IS_APPROX( (m1.cwise()>m2).select(m1,m2), m1.cwise().max(m2) );
117 Scalar mid = (m1.cwise().abs().minCoeff() + m1.cwise().abs().maxCoeff())/Scalar(2);
118 for (int j=0; j<cols; ++j)
119 for (int i=0; i<rows; ++i)
120 m3(i,j) = ei_abs(m1(i,j))<mid ? 0 : m1(i,j);
121 VERIFY_IS_APPROX( (m1.cwise().abs().cwise()<MatrixType::Constant(rows,cols,mid))
122 .select(MatrixType::Zero(rows,cols),m1), m3);
Gael Guennebaude14aa8c2008-09-03 17:56:06 +0000123 // shorter versions:
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000124 VERIFY_IS_APPROX( (m1.cwise().abs().cwise()<MatrixType::Constant(rows,cols,mid))
125 .select(0,m1), m3);
126 VERIFY_IS_APPROX( (m1.cwise().abs().cwise()>=MatrixType::Constant(rows,cols,mid))
127 .select(m1,0), m3);
Gael Guennebaude14aa8c2008-09-03 17:56:06 +0000128 // even shorter version:
129 VERIFY_IS_APPROX( (m1.cwise().abs().cwise()<mid).select(0,m1), m3);
Gael Guennebaud627595a2009-06-10 11:20:30 +0200130
Gael Guennebaud56c7e162009-01-24 15:22:44 +0000131 // count
Benoit Jacob874ff5a2009-01-26 17:56:04 +0000132 VERIFY(((m1.cwise().abs().cwise()+1).cwise()>RealScalar(0.1)).count() == rows*cols);
133 VERIFY_IS_APPROX(((m1.cwise().abs().cwise()+1).cwise()>RealScalar(0.1)).colwise().count(), RowVectorXi::Constant(cols,rows));
134 VERIFY_IS_APPROX(((m1.cwise().abs().cwise()+1).cwise()>RealScalar(0.1)).rowwise().count(), VectorXi::Constant(rows, cols));
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000135}
136
Benoit Jacobe8009992008-11-03 22:47:00 +0000137template<typename VectorType> void lpNorm(const VectorType& v)
138{
139 VectorType u = VectorType::Random(v.size());
140
141 VERIFY_IS_APPROX(u.template lpNorm<Infinity>(), u.cwise().abs().maxCoeff());
142 VERIFY_IS_APPROX(u.template lpNorm<1>(), u.cwise().abs().sum());
143 VERIFY_IS_APPROX(u.template lpNorm<2>(), ei_sqrt(u.cwise().abs().cwise().square().sum()));
144 VERIFY_IS_APPROX(ei_pow(u.template lpNorm<5>(), typename VectorType::RealScalar(5)), u.cwise().abs().cwise().pow(5).sum());
145}
146
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000147void test_array()
148{
149 for(int i = 0; i < g_repeat; i++) {
Benoit Jacob2840ac72009-10-28 18:19:29 -0400150 CALL_SUBTEST_1( array(Matrix<float, 1, 1>()) );
151 CALL_SUBTEST_2( array(Matrix2f()) );
152 CALL_SUBTEST_3( array(Matrix4d()) );
153 CALL_SUBTEST_4( array(MatrixXcf(3, 3)) );
154 CALL_SUBTEST_5( array(MatrixXf(8, 12)) );
155 CALL_SUBTEST_6( array(MatrixXi(8, 12)) );
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000156 }
157 for(int i = 0; i < g_repeat; i++) {
Benoit Jacob2840ac72009-10-28 18:19:29 -0400158 CALL_SUBTEST_1( comparisons(Matrix<float, 1, 1>()) );
159 CALL_SUBTEST_2( comparisons(Matrix2f()) );
160 CALL_SUBTEST_3( comparisons(Matrix4d()) );
161 CALL_SUBTEST_5( comparisons(MatrixXf(8, 12)) );
162 CALL_SUBTEST_6( comparisons(MatrixXi(8, 12)) );
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000163 }
Benoit Jacobe8009992008-11-03 22:47:00 +0000164 for(int i = 0; i < g_repeat; i++) {
Benoit Jacob2840ac72009-10-28 18:19:29 -0400165 CALL_SUBTEST_1( lpNorm(Matrix<float, 1, 1>()) );
166 CALL_SUBTEST_2( lpNorm(Vector2f()) );
167 CALL_SUBTEST_7( lpNorm(Vector3d()) );
168 CALL_SUBTEST_8( lpNorm(Vector4f()) );
169 CALL_SUBTEST_5( lpNorm(VectorXf(16)) );
170 CALL_SUBTEST_4( lpNorm(VectorXcf(10)) );
Benoit Jacobe8009992008-11-03 22:47:00 +0000171 }
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000172}