blob: a308f73660eb9995dc919a222346cb67ef59649e [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"
26#include <Eigen/Array>
27
Gael Guennebaud59dc1da2008-09-03 17:16:28 +000028template<typename MatrixType> void array(const MatrixType& m)
Gael Guennebaud8cef5412008-06-21 17:28:07 +000029{
30 /* this test covers the following files:
31 Array.cpp
32 */
33
34 typedef typename MatrixType::Scalar Scalar;
Gael Guennebaud2120fed2008-08-23 15:14:20 +000035 typedef typename NumTraits<Scalar>::Real RealScalar;
Gael Guennebaud627595a2009-06-10 11:20:30 +020036 typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> ColVectorType;
37 typedef Matrix<Scalar, 1, MatrixType::ColsAtCompileTime> RowVectorType;
Gael Guennebaud8cef5412008-06-21 17:28:07 +000038
39 int rows = m.rows();
40 int cols = m.cols();
41
Benoit Jacob46fe7a32008-09-01 17:31:21 +000042 MatrixType m1 = MatrixType::Random(rows, cols),
43 m2 = MatrixType::Random(rows, cols),
Gael Guennebaud8cef5412008-06-21 17:28:07 +000044 m3(rows, cols);
45
Gael Guennebaud627595a2009-06-10 11:20:30 +020046 ColVectorType cv1 = ColVectorType::Random(rows);
47 RowVectorType rv1 = RowVectorType::Random(cols);
48
Benoit Jacob46fe7a32008-09-01 17:31:21 +000049 Scalar s1 = ei_random<Scalar>(),
50 s2 = ei_random<Scalar>();
Gael Guennebaud8cef5412008-06-21 17:28:07 +000051
Gael Guennebaud59dc1da2008-09-03 17:16:28 +000052 // scalar addition
Benoit Jacobf5791ee2008-07-08 00:49:10 +000053 VERIFY_IS_APPROX(m1.cwise() + s1, s1 + m1.cwise());
Gael Guennebaudc10f0692008-07-21 00:34:46 +000054 VERIFY_IS_APPROX(m1.cwise() + s1, MatrixType::Constant(rows,cols,s1) + m1);
55 VERIFY_IS_APPROX((m1*Scalar(2)).cwise() - s2, (m1+m1) - MatrixType::Constant(rows,cols,s2) );
Gael Guennebaud8cef5412008-06-21 17:28:07 +000056 m3 = m1;
Benoit Jacobf5791ee2008-07-08 00:49:10 +000057 m3.cwise() += s2;
58 VERIFY_IS_APPROX(m3, m1.cwise() + s2);
Gael Guennebaud8cef5412008-06-21 17:28:07 +000059 m3 = m1;
Benoit Jacobf5791ee2008-07-08 00:49:10 +000060 m3.cwise() -= s1;
61 VERIFY_IS_APPROX(m3, m1.cwise() - s1);
Gael Guennebaud05ad0832008-07-19 13:03:23 +000062
Gael Guennebaud59dc1da2008-09-03 17:16:28 +000063 // reductions
Gael Guennebaud05ad0832008-07-19 13:03:23 +000064 VERIFY_IS_APPROX(m1.colwise().sum().sum(), m1.sum());
65 VERIFY_IS_APPROX(m1.rowwise().sum().sum(), m1.sum());
Gael Guennebaud2120fed2008-08-23 15:14:20 +000066 if (!ei_isApprox(m1.sum(), (m1+m2).sum()))
67 VERIFY_IS_NOT_APPROX(((m1+m2).rowwise().sum()).sum(), m1.sum());
Gael Guennebaud05ad0832008-07-19 13:03:23 +000068 VERIFY_IS_APPROX(m1.colwise().sum(), m1.colwise().redux(ei_scalar_sum_op<Scalar>()));
Gael Guennebaud627595a2009-06-10 11:20:30 +020069
70 // vector-wise ops
71 m3 = m1;
72 VERIFY_IS_APPROX(m3.colwise() += cv1, m1.colwise() + cv1);
73 m3 = m1;
74 VERIFY_IS_APPROX(m3.colwise() -= cv1, m1.colwise() - cv1);
75 m3 = m1;
76 VERIFY_IS_APPROX(m3.rowwise() += rv1, m1.rowwise() + rv1);
77 m3 = m1;
78 VERIFY_IS_APPROX(m3.rowwise() -= rv1, m1.rowwise() - rv1);
Gael Guennebaud8cef5412008-06-21 17:28:07 +000079}
80
81template<typename MatrixType> void comparisons(const MatrixType& m)
82{
83 typedef typename MatrixType::Scalar Scalar;
Benoit Jacob874ff5a2009-01-26 17:56:04 +000084 typedef typename NumTraits<Scalar>::Real RealScalar;
Gael Guennebaud8cef5412008-06-21 17:28:07 +000085 typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
86
87 int rows = m.rows();
88 int cols = m.cols();
89
90 int r = ei_random<int>(0, rows-1),
91 c = ei_random<int>(0, cols-1);
92
Gael Guennebaudc10f0692008-07-21 00:34:46 +000093 MatrixType m1 = MatrixType::Random(rows, cols),
94 m2 = MatrixType::Random(rows, cols),
Gael Guennebaud8cef5412008-06-21 17:28:07 +000095 m3(rows, cols);
96
Benoit Jacobf5791ee2008-07-08 00:49:10 +000097 VERIFY(((m1.cwise() + Scalar(1)).cwise() > m1).all());
98 VERIFY(((m1.cwise() - Scalar(1)).cwise() < m1).all());
Gael Guennebaud8cef5412008-06-21 17:28:07 +000099 if (rows*cols>1)
100 {
101 m3 = m1;
102 m3(r,c) += 1;
Benoit Jacobf5791ee2008-07-08 00:49:10 +0000103 VERIFY(! (m1.cwise() < m3).all() );
104 VERIFY(! (m1.cwise() > 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
108 VERIFY( (m1.cwise() != (m1(r,c)+1) ).any() );
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) ).any() );
Gael Guennebaud627595a2009-06-10 11:20:30 +0200112
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000113 // test Select
114 VERIFY_IS_APPROX( (m1.cwise()<m2).select(m1,m2), m1.cwise().min(m2) );
115 VERIFY_IS_APPROX( (m1.cwise()>m2).select(m1,m2), m1.cwise().max(m2) );
116 Scalar mid = (m1.cwise().abs().minCoeff() + m1.cwise().abs().maxCoeff())/Scalar(2);
117 for (int j=0; j<cols; ++j)
118 for (int i=0; i<rows; ++i)
119 m3(i,j) = ei_abs(m1(i,j))<mid ? 0 : m1(i,j);
120 VERIFY_IS_APPROX( (m1.cwise().abs().cwise()<MatrixType::Constant(rows,cols,mid))
121 .select(MatrixType::Zero(rows,cols),m1), m3);
Gael Guennebaude14aa8c2008-09-03 17:56:06 +0000122 // shorter versions:
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000123 VERIFY_IS_APPROX( (m1.cwise().abs().cwise()<MatrixType::Constant(rows,cols,mid))
124 .select(0,m1), m3);
125 VERIFY_IS_APPROX( (m1.cwise().abs().cwise()>=MatrixType::Constant(rows,cols,mid))
126 .select(m1,0), m3);
Gael Guennebaude14aa8c2008-09-03 17:56:06 +0000127 // even shorter version:
128 VERIFY_IS_APPROX( (m1.cwise().abs().cwise()<mid).select(0,m1), m3);
Gael Guennebaud627595a2009-06-10 11:20:30 +0200129
Gael Guennebaud56c7e162009-01-24 15:22:44 +0000130 // count
Benoit Jacob874ff5a2009-01-26 17:56:04 +0000131 VERIFY(((m1.cwise().abs().cwise()+1).cwise()>RealScalar(0.1)).count() == rows*cols);
132 VERIFY_IS_APPROX(((m1.cwise().abs().cwise()+1).cwise()>RealScalar(0.1)).colwise().count(), RowVectorXi::Constant(cols,rows));
133 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 +0000134}
135
Benoit Jacobe8009992008-11-03 22:47:00 +0000136template<typename VectorType> void lpNorm(const VectorType& v)
137{
138 VectorType u = VectorType::Random(v.size());
139
140 VERIFY_IS_APPROX(u.template lpNorm<Infinity>(), u.cwise().abs().maxCoeff());
141 VERIFY_IS_APPROX(u.template lpNorm<1>(), u.cwise().abs().sum());
142 VERIFY_IS_APPROX(u.template lpNorm<2>(), ei_sqrt(u.cwise().abs().cwise().square().sum()));
143 VERIFY_IS_APPROX(ei_pow(u.template lpNorm<5>(), typename VectorType::RealScalar(5)), u.cwise().abs().cwise().pow(5).sum());
144}
145
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000146void test_array()
147{
148 for(int i = 0; i < g_repeat; i++) {
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000149 CALL_SUBTEST( array(Matrix<float, 1, 1>()) );
150 CALL_SUBTEST( array(Matrix2f()) );
151 CALL_SUBTEST( array(Matrix4d()) );
152 CALL_SUBTEST( array(MatrixXcf(3, 3)) );
153 CALL_SUBTEST( array(MatrixXf(8, 12)) );
154 CALL_SUBTEST( array(MatrixXi(8, 12)) );
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000155 }
156 for(int i = 0; i < g_repeat; i++) {
157 CALL_SUBTEST( comparisons(Matrix<float, 1, 1>()) );
158 CALL_SUBTEST( comparisons(Matrix2f()) );
159 CALL_SUBTEST( comparisons(Matrix4d()) );
160 CALL_SUBTEST( comparisons(MatrixXf(8, 12)) );
161 CALL_SUBTEST( comparisons(MatrixXi(8, 12)) );
162 }
Benoit Jacobe8009992008-11-03 22:47:00 +0000163 for(int i = 0; i < g_repeat; i++) {
164 CALL_SUBTEST( lpNorm(Matrix<float, 1, 1>()) );
165 CALL_SUBTEST( lpNorm(Vector2f()) );
166 CALL_SUBTEST( lpNorm(Vector3d()) );
167 CALL_SUBTEST( lpNorm(Vector4f()) );
168 CALL_SUBTEST( lpNorm(VectorXf(16)) );
169 CALL_SUBTEST( lpNorm(VectorXcd(10)) );
170 }
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000171}