blob: e51dbac2aaf5579060cf7c912e32aee59e548703 [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
25#include "main.h"
Gael Guennebaud8cef5412008-06-21 17:28:07 +000026
Gael Guennebaud59dc1da2008-09-03 17:16:28 +000027template<typename MatrixType> void array(const MatrixType& m)
Gael Guennebaud8cef5412008-06-21 17:28:07 +000028{
Gael Guennebaud8cef5412008-06-21 17:28:07 +000029 typedef typename MatrixType::Scalar Scalar;
Gael Guennebaud2120fed2008-08-23 15:14:20 +000030 typedef typename NumTraits<Scalar>::Real RealScalar;
Gael Guennebaudc70d5422010-01-18 22:54:20 +010031 typedef Array<Scalar, MatrixType::RowsAtCompileTime, 1> ColVectorType;
32 typedef Array<Scalar, 1, MatrixType::ColsAtCompileTime> RowVectorType;
Gael Guennebaud8cef5412008-06-21 17:28:07 +000033
34 int rows = m.rows();
35 int cols = m.cols();
36
Benoit Jacob46fe7a32008-09-01 17:31:21 +000037 MatrixType m1 = MatrixType::Random(rows, cols),
38 m2 = MatrixType::Random(rows, cols),
Gael Guennebaud8cef5412008-06-21 17:28:07 +000039 m3(rows, cols);
40
Gael Guennebaud627595a2009-06-10 11:20:30 +020041 ColVectorType cv1 = ColVectorType::Random(rows);
42 RowVectorType rv1 = RowVectorType::Random(cols);
43
Benoit Jacob46fe7a32008-09-01 17:31:21 +000044 Scalar s1 = ei_random<Scalar>(),
45 s2 = ei_random<Scalar>();
Gael Guennebaud8cef5412008-06-21 17:28:07 +000046
Gael Guennebaud59dc1da2008-09-03 17:16:28 +000047 // scalar addition
Gael Guennebaudc70d5422010-01-18 22:54:20 +010048 VERIFY_IS_APPROX(m1 + s1, s1 + m1);
49 VERIFY_IS_APPROX(m1 + s1, MatrixType::Constant(rows,cols,s1) + m1);
50 VERIFY_IS_APPROX(s1 - m1, (-m1)+s1 );
51 VERIFY_IS_APPROX(m1 - s1, m1 - MatrixType::Constant(rows,cols,s1));
52 VERIFY_IS_APPROX(s1 - m1, MatrixType::Constant(rows,cols,s1) - m1);
53 VERIFY_IS_APPROX((m1*Scalar(2)) - s2, (m1+m1) - MatrixType::Constant(rows,cols,s2) );
Gael Guennebaud8cef5412008-06-21 17:28:07 +000054 m3 = m1;
Gael Guennebaudc70d5422010-01-18 22:54:20 +010055 m3 += s2;
56 VERIFY_IS_APPROX(m3, m1 + s2);
Gael Guennebaud8cef5412008-06-21 17:28:07 +000057 m3 = m1;
Gael Guennebaudc70d5422010-01-18 22:54:20 +010058 m3 -= s1;
59 VERIFY_IS_APPROX(m3, m1 - s1);
Gael Guennebaud05ad0832008-07-19 13:03:23 +000060
Gael Guennebaud59dc1da2008-09-03 17:16:28 +000061 // reductions
Gael Guennebaud05ad0832008-07-19 13:03:23 +000062 VERIFY_IS_APPROX(m1.colwise().sum().sum(), m1.sum());
63 VERIFY_IS_APPROX(m1.rowwise().sum().sum(), m1.sum());
Gael Guennebaud2120fed2008-08-23 15:14:20 +000064 if (!ei_isApprox(m1.sum(), (m1+m2).sum()))
65 VERIFY_IS_NOT_APPROX(((m1+m2).rowwise().sum()).sum(), m1.sum());
Gael Guennebaud05ad0832008-07-19 13:03:23 +000066 VERIFY_IS_APPROX(m1.colwise().sum(), m1.colwise().redux(ei_scalar_sum_op<Scalar>()));
Gael Guennebaud627595a2009-06-10 11:20:30 +020067
68 // vector-wise ops
69 m3 = m1;
70 VERIFY_IS_APPROX(m3.colwise() += cv1, m1.colwise() + cv1);
71 m3 = m1;
72 VERIFY_IS_APPROX(m3.colwise() -= cv1, m1.colwise() - cv1);
73 m3 = m1;
74 VERIFY_IS_APPROX(m3.rowwise() += rv1, m1.rowwise() + rv1);
75 m3 = m1;
76 VERIFY_IS_APPROX(m3.rowwise() -= rv1, m1.rowwise() - rv1);
Gael Guennebaud8cef5412008-06-21 17:28:07 +000077}
78
79template<typename MatrixType> void comparisons(const MatrixType& m)
80{
81 typedef typename MatrixType::Scalar Scalar;
Benoit Jacob874ff5a2009-01-26 17:56:04 +000082 typedef typename NumTraits<Scalar>::Real RealScalar;
Gael Guennebaudc70d5422010-01-18 22:54:20 +010083 typedef Array<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
Gael Guennebaud8cef5412008-06-21 17:28:07 +000084
85 int rows = m.rows();
86 int cols = m.cols();
87
88 int r = ei_random<int>(0, rows-1),
89 c = ei_random<int>(0, cols-1);
90
Gael Guennebaudc10f0692008-07-21 00:34:46 +000091 MatrixType m1 = MatrixType::Random(rows, cols),
92 m2 = MatrixType::Random(rows, cols),
Gael Guennebaud8cef5412008-06-21 17:28:07 +000093 m3(rows, cols);
94
Gael Guennebaudc70d5422010-01-18 22:54:20 +010095 VERIFY(((m1 + Scalar(1)) > m1).all());
96 VERIFY(((m1 - Scalar(1)) < m1).all());
Gael Guennebaud8cef5412008-06-21 17:28:07 +000097 if (rows*cols>1)
98 {
99 m3 = m1;
100 m3(r,c) += 1;
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100101 VERIFY(! (m1 < m3).all() );
102 VERIFY(! (m1 > m3).all() );
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000103 }
Gael Guennebaud627595a2009-06-10 11:20:30 +0200104
Gael Guennebaude14aa8c2008-09-03 17:56:06 +0000105 // comparisons to scalar
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100106 VERIFY( (m1 != (m1(r,c)+1) ).any() );
107 VERIFY( (m1 > (m1(r,c)-1) ).any() );
108 VERIFY( (m1 < (m1(r,c)+1) ).any() );
109 VERIFY( (m1 == m1(r,c) ).any() );
Gael Guennebaud627595a2009-06-10 11:20:30 +0200110
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000111 // test Select
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100112 VERIFY_IS_APPROX( (m1<m2).select(m1,m2), m1.cwiseMin(m2) );
113 VERIFY_IS_APPROX( (m1>m2).select(m1,m2), m1.cwiseMax(m2) );
Gael Guennebaud9f795582009-12-16 19:18:40 +0100114 Scalar mid = (m1.cwiseAbs().minCoeff() + m1.cwiseAbs().maxCoeff())/Scalar(2);
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000115 for (int j=0; j<cols; ++j)
116 for (int i=0; i<rows; ++i)
117 m3(i,j) = ei_abs(m1(i,j))<mid ? 0 : m1(i,j);
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100118 VERIFY_IS_APPROX( (m1.abs()<MatrixType::Constant(rows,cols,mid))
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000119 .select(MatrixType::Zero(rows,cols),m1), m3);
Gael Guennebaude14aa8c2008-09-03 17:56:06 +0000120 // shorter versions:
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100121 VERIFY_IS_APPROX( (m1.abs()<MatrixType::Constant(rows,cols,mid))
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000122 .select(0,m1), m3);
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100123 VERIFY_IS_APPROX( (m1.abs()>=MatrixType::Constant(rows,cols,mid))
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000124 .select(m1,0), m3);
Gael Guennebaude14aa8c2008-09-03 17:56:06 +0000125 // even shorter version:
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100126 VERIFY_IS_APPROX( (m1.abs()<mid).select(0,m1), m3);
Gael Guennebaud627595a2009-06-10 11:20:30 +0200127
Gael Guennebaud56c7e162009-01-24 15:22:44 +0000128 // count
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100129 VERIFY(((m1.abs()+1)>RealScalar(0.1)).count() == rows*cols);
Gael Guennebaud9f795582009-12-16 19:18:40 +0100130 // TODO allows colwise/rowwise for array
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100131 VERIFY_IS_APPROX(((m1.abs()+1)>RealScalar(0.1)).colwise().count(), ArrayXi::Constant(cols,rows).transpose());
132 VERIFY_IS_APPROX(((m1.abs()+1)>RealScalar(0.1)).rowwise().count(), ArrayXi::Constant(rows, cols));
Benoit Jacobe8009992008-11-03 22:47:00 +0000133}
134
Gael Guennebaud0ce5bc02010-01-27 23:23:59 +0100135template<typename MatrixType> void array_real(const MatrixType& m)
136{
137 typedef typename MatrixType::Scalar Scalar;
138 typedef typename NumTraits<Scalar>::Real RealScalar;
139
140 int rows = m.rows();
141 int cols = m.cols();
142
143 MatrixType m1 = MatrixType::Random(rows, cols),
144 m2 = MatrixType::Random(rows, cols),
145 m3(rows, cols);
146
147 VERIFY_IS_APPROX(m1.sin(), std::sin(m1));
148 VERIFY_IS_APPROX(m1.sin(), ei_sin(m1));
149 VERIFY_IS_APPROX(m1.cos(), ei_cos(m1));
150 VERIFY_IS_APPROX(m1.cos(), ei_cos(m1));
151 VERIFY_IS_APPROX(m1.abs().sqrt(), std::sqrt(std::abs(m1)));
152 VERIFY_IS_APPROX(m1.abs().sqrt(), ei_sqrt(ei_abs(m1)));
153 VERIFY_IS_APPROX(m1.abs().log(), std::log(std::abs(m1)));
154 VERIFY_IS_APPROX(m1.abs().log(), ei_log(ei_abs(m1)));
155 VERIFY_IS_APPROX(m1.exp(), std::exp(m1));
156 VERIFY_IS_APPROX(m1.exp(), ei_exp(m1));
157}
158
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000159void test_array()
160{
161 for(int i = 0; i < g_repeat; i++) {
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100162 CALL_SUBTEST_1( array(Array<float, 1, 1>()) );
163 CALL_SUBTEST_2( array(Array22f()) );
164 CALL_SUBTEST_3( array(Array44d()) );
165 CALL_SUBTEST_4( array(ArrayXXcf(3, 3)) );
166 CALL_SUBTEST_5( array(ArrayXXf(8, 12)) );
167 CALL_SUBTEST_6( array(ArrayXXi(8, 12)) );
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000168 }
169 for(int i = 0; i < g_repeat; i++) {
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100170 CALL_SUBTEST_1( comparisons(Array<float, 1, 1>()) );
171 CALL_SUBTEST_2( comparisons(Array22f()) );
172 CALL_SUBTEST_3( comparisons(Array44d()) );
173 CALL_SUBTEST_5( comparisons(ArrayXXf(8, 12)) );
174 CALL_SUBTEST_6( comparisons(ArrayXXi(8, 12)) );
Benoit Jacobe8009992008-11-03 22:47:00 +0000175 }
Gael Guennebaud0ce5bc02010-01-27 23:23:59 +0100176 for(int i = 0; i < g_repeat; i++) {
177 CALL_SUBTEST_1( array_real(Array<float, 1, 1>()) );
178 CALL_SUBTEST_2( array_real(Array22f()) );
179 CALL_SUBTEST_3( array_real(Array44d()) );
180 CALL_SUBTEST_5( array_real(ArrayXXf(8, 12)) );
181 }
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000182}