blob: 72d3584e650cb5900354c3c5496b623c7c532853 [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 Guennebaud28e64b02010-06-24 23:21:58 +02004// Copyright (C) 2008-2009 Gael Guennebaud <gael.guennebaud@inria.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
25#include "main.h"
Gael Guennebaud8cef5412008-06-21 17:28:07 +000026
Benoit Jacobe2775862010-04-28 18:51:38 -040027template<typename ArrayType> void array(const ArrayType& m)
Gael Guennebaud8cef5412008-06-21 17:28:07 +000028{
Hauke Heibelf1679c72010-06-20 17:37:56 +020029 typedef typename ArrayType::Index Index;
Benoit Jacobe2775862010-04-28 18:51:38 -040030 typedef typename ArrayType::Scalar Scalar;
Gael Guennebaud2120fed2008-08-23 15:14:20 +000031 typedef typename NumTraits<Scalar>::Real RealScalar;
Benoit Jacobe2775862010-04-28 18:51:38 -040032 typedef Array<Scalar, ArrayType::RowsAtCompileTime, 1> ColVectorType;
33 typedef Array<Scalar, 1, ArrayType::ColsAtCompileTime> RowVectorType;
Gael Guennebaud8cef5412008-06-21 17:28:07 +000034
Hauke Heibelf1679c72010-06-20 17:37:56 +020035 Index rows = m.rows();
36 Index cols = m.cols();
Gael Guennebaud8cef5412008-06-21 17:28:07 +000037
Benoit Jacobe2775862010-04-28 18:51:38 -040038 ArrayType m1 = ArrayType::Random(rows, cols),
39 m2 = ArrayType::Random(rows, cols),
Gael Guennebaud8cef5412008-06-21 17:28:07 +000040 m3(rows, cols);
41
Gael Guennebaud627595a2009-06-10 11:20:30 +020042 ColVectorType cv1 = ColVectorType::Random(rows);
43 RowVectorType rv1 = RowVectorType::Random(cols);
44
Benoit Jacob47160402010-10-25 10:15:22 -040045 Scalar s1 = internal::random<Scalar>(),
46 s2 = internal::random<Scalar>();
Gael Guennebaud8cef5412008-06-21 17:28:07 +000047
Gael Guennebaud59dc1da2008-09-03 17:16:28 +000048 // scalar addition
Gael Guennebaudc70d5422010-01-18 22:54:20 +010049 VERIFY_IS_APPROX(m1 + s1, s1 + m1);
Benoit Jacobe2775862010-04-28 18:51:38 -040050 VERIFY_IS_APPROX(m1 + s1, ArrayType::Constant(rows,cols,s1) + m1);
Gael Guennebaudc70d5422010-01-18 22:54:20 +010051 VERIFY_IS_APPROX(s1 - m1, (-m1)+s1 );
Benoit Jacobe2775862010-04-28 18:51:38 -040052 VERIFY_IS_APPROX(m1 - s1, m1 - ArrayType::Constant(rows,cols,s1));
53 VERIFY_IS_APPROX(s1 - m1, ArrayType::Constant(rows,cols,s1) - m1);
54 VERIFY_IS_APPROX((m1*Scalar(2)) - s2, (m1+m1) - ArrayType::Constant(rows,cols,s2) );
Gael Guennebaud8cef5412008-06-21 17:28:07 +000055 m3 = m1;
Gael Guennebaudc70d5422010-01-18 22:54:20 +010056 m3 += s2;
57 VERIFY_IS_APPROX(m3, m1 + s2);
Gael Guennebaud8cef5412008-06-21 17:28:07 +000058 m3 = m1;
Gael Guennebaudc70d5422010-01-18 22:54:20 +010059 m3 -= s1;
60 VERIFY_IS_APPROX(m3, m1 - s1);
Gael Guennebaud05ad0832008-07-19 13:03:23 +000061
Gael Guennebaud59dc1da2008-09-03 17:16:28 +000062 // reductions
Gael Guennebaud05ad0832008-07-19 13:03:23 +000063 VERIFY_IS_APPROX(m1.colwise().sum().sum(), m1.sum());
64 VERIFY_IS_APPROX(m1.rowwise().sum().sum(), m1.sum());
Benoit Jacob47160402010-10-25 10:15:22 -040065 if (!internal::isApprox(m1.sum(), (m1+m2).sum()))
Gael Guennebaud2120fed2008-08-23 15:14:20 +000066 VERIFY_IS_NOT_APPROX(((m1+m2).rowwise().sum()).sum(), m1.sum());
Benoit Jacob47160402010-10-25 10:15:22 -040067 VERIFY_IS_APPROX(m1.colwise().sum(), m1.colwise().redux(internal::scalar_sum_op<Scalar>()));
Gael Guennebaud627595a2009-06-10 11:20:30 +020068
69 // vector-wise ops
70 m3 = m1;
71 VERIFY_IS_APPROX(m3.colwise() += cv1, m1.colwise() + cv1);
72 m3 = m1;
73 VERIFY_IS_APPROX(m3.colwise() -= cv1, m1.colwise() - cv1);
74 m3 = m1;
75 VERIFY_IS_APPROX(m3.rowwise() += rv1, m1.rowwise() + rv1);
76 m3 = m1;
77 VERIFY_IS_APPROX(m3.rowwise() -= rv1, m1.rowwise() - rv1);
Gael Guennebaud8cef5412008-06-21 17:28:07 +000078}
79
Benoit Jacobe2775862010-04-28 18:51:38 -040080template<typename ArrayType> void comparisons(const ArrayType& m)
Gael Guennebaud8cef5412008-06-21 17:28:07 +000081{
Hauke Heibelf1679c72010-06-20 17:37:56 +020082 typedef typename ArrayType::Index Index;
Benoit Jacobe2775862010-04-28 18:51:38 -040083 typedef typename ArrayType::Scalar Scalar;
Benoit Jacob874ff5a2009-01-26 17:56:04 +000084 typedef typename NumTraits<Scalar>::Real RealScalar;
Benoit Jacobe2775862010-04-28 18:51:38 -040085 typedef Array<Scalar, ArrayType::RowsAtCompileTime, 1> VectorType;
Gael Guennebaud8cef5412008-06-21 17:28:07 +000086
Hauke Heibelf1679c72010-06-20 17:37:56 +020087 Index rows = m.rows();
88 Index cols = m.cols();
Gael Guennebaud8cef5412008-06-21 17:28:07 +000089
Benoit Jacob47160402010-10-25 10:15:22 -040090 Index r = internal::random<Index>(0, rows-1),
91 c = internal::random<Index>(0, cols-1);
Gael Guennebaud8cef5412008-06-21 17:28:07 +000092
Benoit Jacobe2775862010-04-28 18:51:38 -040093 ArrayType m1 = ArrayType::Random(rows, cols),
94 m2 = ArrayType::Random(rows, cols),
Gael Guennebaud8cef5412008-06-21 17:28:07 +000095 m3(rows, cols);
96
Gael Guennebaudc70d5422010-01-18 22:54:20 +010097 VERIFY(((m1 + Scalar(1)) > m1).all());
98 VERIFY(((m1 - Scalar(1)) < m1).all());
Gael Guennebaud8cef5412008-06-21 17:28:07 +000099 if (rows*cols>1)
100 {
101 m3 = m1;
102 m3(r,c) += 1;
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100103 VERIFY(! (m1 < m3).all() );
104 VERIFY(! (m1 > m3).all() );
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000105 }
Gael Guennebaud627595a2009-06-10 11:20:30 +0200106
Gael Guennebaude14aa8c2008-09-03 17:56:06 +0000107 // comparisons to scalar
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100108 VERIFY( (m1 != (m1(r,c)+1) ).any() );
109 VERIFY( (m1 > (m1(r,c)-1) ).any() );
110 VERIFY( (m1 < (m1(r,c)+1) ).any() );
111 VERIFY( (m1 == m1(r,c) ).any() );
Gael Guennebaud627595a2009-06-10 11:20:30 +0200112
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000113 // test Select
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100114 VERIFY_IS_APPROX( (m1<m2).select(m1,m2), m1.cwiseMin(m2) );
115 VERIFY_IS_APPROX( (m1>m2).select(m1,m2), m1.cwiseMax(m2) );
Gael Guennebaud9f795582009-12-16 19:18:40 +0100116 Scalar mid = (m1.cwiseAbs().minCoeff() + m1.cwiseAbs().maxCoeff())/Scalar(2);
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000117 for (int j=0; j<cols; ++j)
118 for (int i=0; i<rows; ++i)
Benoit Jacob47160402010-10-25 10:15:22 -0400119 m3(i,j) = internal::abs(m1(i,j))<mid ? 0 : m1(i,j);
Benoit Jacobe2775862010-04-28 18:51:38 -0400120 VERIFY_IS_APPROX( (m1.abs()<ArrayType::Constant(rows,cols,mid))
121 .select(ArrayType::Zero(rows,cols),m1), m3);
Gael Guennebaude14aa8c2008-09-03 17:56:06 +0000122 // shorter versions:
Benoit Jacobe2775862010-04-28 18:51:38 -0400123 VERIFY_IS_APPROX( (m1.abs()<ArrayType::Constant(rows,cols,mid))
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000124 .select(0,m1), m3);
Benoit Jacobe2775862010-04-28 18:51:38 -0400125 VERIFY_IS_APPROX( (m1.abs()>=ArrayType::Constant(rows,cols,mid))
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000126 .select(m1,0), m3);
Gael Guennebaude14aa8c2008-09-03 17:56:06 +0000127 // even shorter version:
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100128 VERIFY_IS_APPROX( (m1.abs()<mid).select(0,m1), m3);
Gael Guennebaud627595a2009-06-10 11:20:30 +0200129
Gael Guennebaud56c7e162009-01-24 15:22:44 +0000130 // count
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100131 VERIFY(((m1.abs()+1)>RealScalar(0.1)).count() == rows*cols);
Benoit Jacobaaaade42010-05-30 16:00:58 -0400132
133 typedef Array<typename ArrayType::Index, Dynamic, 1> ArrayOfIndices;
Gael Guennebaud4ebb8042010-06-02 09:45:57 +0200134
Gael Guennebaud9f795582009-12-16 19:18:40 +0100135 // TODO allows colwise/rowwise for array
Benoit Jacobaaaade42010-05-30 16:00:58 -0400136 VERIFY_IS_APPROX(((m1.abs()+1)>RealScalar(0.1)).colwise().count(), ArrayOfIndices::Constant(cols,rows).transpose());
137 VERIFY_IS_APPROX(((m1.abs()+1)>RealScalar(0.1)).rowwise().count(), ArrayOfIndices::Constant(rows, cols));
Benoit Jacobe8009992008-11-03 22:47:00 +0000138}
139
Benoit Jacobe2775862010-04-28 18:51:38 -0400140template<typename ArrayType> void array_real(const ArrayType& m)
Gael Guennebaud0ce5bc02010-01-27 23:23:59 +0100141{
Hauke Heibelf1679c72010-06-20 17:37:56 +0200142 typedef typename ArrayType::Index Index;
Benoit Jacobe2775862010-04-28 18:51:38 -0400143 typedef typename ArrayType::Scalar Scalar;
Gael Guennebaud0ce5bc02010-01-27 23:23:59 +0100144 typedef typename NumTraits<Scalar>::Real RealScalar;
145
Hauke Heibelf1679c72010-06-20 17:37:56 +0200146 Index rows = m.rows();
147 Index cols = m.cols();
Gael Guennebaud0ce5bc02010-01-27 23:23:59 +0100148
Benoit Jacobe2775862010-04-28 18:51:38 -0400149 ArrayType m1 = ArrayType::Random(rows, cols),
150 m2 = ArrayType::Random(rows, cols),
Gael Guennebaud0ce5bc02010-01-27 23:23:59 +0100151 m3(rows, cols);
152
153 VERIFY_IS_APPROX(m1.sin(), std::sin(m1));
Benoit Jacob47160402010-10-25 10:15:22 -0400154 VERIFY_IS_APPROX(m1.sin(), internal::sin(m1));
Benoit Jacobe2775862010-04-28 18:51:38 -0400155 VERIFY_IS_APPROX(m1.cos(), std::cos(m1));
Benoit Jacob47160402010-10-25 10:15:22 -0400156 VERIFY_IS_APPROX(m1.cos(), internal::cos(m1));
Gael Guennebaud4ebb8042010-06-02 09:45:57 +0200157
Benoit Jacob47160402010-10-25 10:15:22 -0400158 VERIFY_IS_APPROX(internal::cos(m1+RealScalar(3)*m2), internal::cos((m1+RealScalar(3)*m2).eval()));
Benoit Jacobe2775862010-04-28 18:51:38 -0400159 VERIFY_IS_APPROX(std::cos(m1+RealScalar(3)*m2), std::cos((m1+RealScalar(3)*m2).eval()));
Gael Guennebaud4ebb8042010-06-02 09:45:57 +0200160
Gael Guennebaud0ce5bc02010-01-27 23:23:59 +0100161 VERIFY_IS_APPROX(m1.abs().sqrt(), std::sqrt(std::abs(m1)));
Benoit Jacob47160402010-10-25 10:15:22 -0400162 VERIFY_IS_APPROX(m1.abs().sqrt(), internal::sqrt(internal::abs(m1)));
163 VERIFY_IS_APPROX(m1.abs(), internal::sqrt(internal::abs2(m1)));
Benoit Jacob5d63d2c2010-04-28 22:42:34 -0400164
Benoit Jacob47160402010-10-25 10:15:22 -0400165 VERIFY_IS_APPROX(internal::abs2(internal::real(m1)) + internal::abs2(internal::imag(m1)), internal::abs2(m1));
166 VERIFY_IS_APPROX(internal::abs2(std::real(m1)) + internal::abs2(std::imag(m1)), internal::abs2(m1));
Benoit Jacob5d63d2c2010-04-28 22:42:34 -0400167 if(!NumTraits<Scalar>::IsComplex)
Benoit Jacob47160402010-10-25 10:15:22 -0400168 VERIFY_IS_APPROX(internal::real(m1), m1);
Gael Guennebaud4ebb8042010-06-02 09:45:57 +0200169
Gael Guennebaud0ce5bc02010-01-27 23:23:59 +0100170 VERIFY_IS_APPROX(m1.abs().log(), std::log(std::abs(m1)));
Benoit Jacob47160402010-10-25 10:15:22 -0400171 VERIFY_IS_APPROX(m1.abs().log(), internal::log(internal::abs(m1)));
Gael Guennebaud4ebb8042010-06-02 09:45:57 +0200172
Gael Guennebaud0ce5bc02010-01-27 23:23:59 +0100173 VERIFY_IS_APPROX(m1.exp(), std::exp(m1));
Benoit Jacobe2775862010-04-28 18:51:38 -0400174 VERIFY_IS_APPROX(m1.exp() * m2.exp(), std::exp(m1+m2));
Benoit Jacob47160402010-10-25 10:15:22 -0400175 VERIFY_IS_APPROX(m1.exp(), internal::exp(m1));
Benoit Jacobe2775862010-04-28 18:51:38 -0400176 VERIFY_IS_APPROX(m1.exp() / m2.exp(), std::exp(m1-m2));
Gael Guennebaud575ac542010-06-19 23:17:07 +0200177
178 VERIFY_IS_APPROX(m1.pow(2), m1.square());
179 VERIFY_IS_APPROX(std::pow(m1,2), m1.square());
180 m3 = m1.abs();
181 VERIFY_IS_APPROX(m3.pow(RealScalar(0.5)), m3.sqrt());
182 VERIFY_IS_APPROX(std::pow(m3,RealScalar(0.5)), m3.sqrt());
Gael Guennebaud0ce5bc02010-01-27 23:23:59 +0100183}
184
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000185void test_array()
186{
187 for(int i = 0; i < g_repeat; i++) {
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100188 CALL_SUBTEST_1( array(Array<float, 1, 1>()) );
189 CALL_SUBTEST_2( array(Array22f()) );
190 CALL_SUBTEST_3( array(Array44d()) );
191 CALL_SUBTEST_4( array(ArrayXXcf(3, 3)) );
192 CALL_SUBTEST_5( array(ArrayXXf(8, 12)) );
193 CALL_SUBTEST_6( array(ArrayXXi(8, 12)) );
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000194 }
195 for(int i = 0; i < g_repeat; i++) {
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100196 CALL_SUBTEST_1( comparisons(Array<float, 1, 1>()) );
197 CALL_SUBTEST_2( comparisons(Array22f()) );
198 CALL_SUBTEST_3( comparisons(Array44d()) );
199 CALL_SUBTEST_5( comparisons(ArrayXXf(8, 12)) );
200 CALL_SUBTEST_6( comparisons(ArrayXXi(8, 12)) );
Benoit Jacobe8009992008-11-03 22:47:00 +0000201 }
Gael Guennebaud0ce5bc02010-01-27 23:23:59 +0100202 for(int i = 0; i < g_repeat; i++) {
203 CALL_SUBTEST_1( array_real(Array<float, 1, 1>()) );
204 CALL_SUBTEST_2( array_real(Array22f()) );
205 CALL_SUBTEST_3( array_real(Array44d()) );
206 CALL_SUBTEST_5( array_real(ArrayXXf(8, 12)) );
207 }
Benoit Jacob5d63d2c2010-04-28 22:42:34 -0400208
Hauke Heibel7bc8e3a2010-10-25 22:13:49 +0200209 VERIFY((internal::is_same< internal::global_math_functions_filtering_base<int>::type, int >::value));
210 VERIFY((internal::is_same< internal::global_math_functions_filtering_base<float>::type, float >::value));
211 VERIFY((internal::is_same< internal::global_math_functions_filtering_base<Array2i>::type, ArrayBase<Array2i> >::value));
Benoit Jacob47160402010-10-25 10:15:22 -0400212 typedef CwiseUnaryOp<internal::scalar_sum_op<double>, ArrayXd > Xpr;
Hauke Heibel7bc8e3a2010-10-25 22:13:49 +0200213 VERIFY((internal::is_same< internal::global_math_functions_filtering_base<Xpr>::type,
Benoit Jacob5d63d2c2010-04-28 22:42:34 -0400214 ArrayBase<Xpr>
Hauke Heibel7bc8e3a2010-10-25 22:13:49 +0200215 >::value));
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000216}