blob: 0fa13b65b93b2b6f31e1a4d552db92500c37f2a3 [file] [log] [blame]
Gael Guennebaud8cef5412008-06-21 17:28:07 +00001// 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// 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 Guennebaud8cef5412008-06-21 17:28:07 +000036 typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
37
38 int rows = m.rows();
39 int cols = m.cols();
40
Benoit Jacob46fe7a32008-09-01 17:31:21 +000041 MatrixType m1 = MatrixType::Random(rows, cols),
42 m2 = MatrixType::Random(rows, cols),
Gael Guennebaud8cef5412008-06-21 17:28:07 +000043 m3(rows, cols);
44
Benoit Jacob46fe7a32008-09-01 17:31:21 +000045 Scalar s1 = ei_random<Scalar>(),
46 s2 = ei_random<Scalar>();
Gael Guennebaud8cef5412008-06-21 17:28:07 +000047
Gael Guennebaud59dc1da2008-09-03 17:16:28 +000048 // scalar addition
Benoit Jacobf5791ee2008-07-08 00:49:10 +000049 VERIFY_IS_APPROX(m1.cwise() + s1, s1 + m1.cwise());
Gael Guennebaudc10f0692008-07-21 00:34:46 +000050 VERIFY_IS_APPROX(m1.cwise() + s1, MatrixType::Constant(rows,cols,s1) + m1);
51 VERIFY_IS_APPROX((m1*Scalar(2)).cwise() - s2, (m1+m1) - MatrixType::Constant(rows,cols,s2) );
Gael Guennebaud8cef5412008-06-21 17:28:07 +000052 m3 = m1;
Benoit Jacobf5791ee2008-07-08 00:49:10 +000053 m3.cwise() += s2;
54 VERIFY_IS_APPROX(m3, m1.cwise() + s2);
Gael Guennebaud8cef5412008-06-21 17:28:07 +000055 m3 = m1;
Benoit Jacobf5791ee2008-07-08 00:49:10 +000056 m3.cwise() -= s1;
57 VERIFY_IS_APPROX(m3, m1.cwise() - s1);
Gael Guennebaud05ad0832008-07-19 13:03:23 +000058
Gael Guennebaud59dc1da2008-09-03 17:16:28 +000059 // reductions
Gael Guennebaud05ad0832008-07-19 13:03:23 +000060 VERIFY_IS_APPROX(m1.colwise().sum().sum(), m1.sum());
61 VERIFY_IS_APPROX(m1.rowwise().sum().sum(), m1.sum());
Gael Guennebaud2120fed2008-08-23 15:14:20 +000062 if (!ei_isApprox(m1.sum(), (m1+m2).sum()))
63 VERIFY_IS_NOT_APPROX(((m1+m2).rowwise().sum()).sum(), m1.sum());
Gael Guennebaud05ad0832008-07-19 13:03:23 +000064 VERIFY_IS_APPROX(m1.colwise().sum(), m1.colwise().redux(ei_scalar_sum_op<Scalar>()));
Gael Guennebaud8cef5412008-06-21 17:28:07 +000065}
66
67template<typename MatrixType> void comparisons(const MatrixType& m)
68{
69 typedef typename MatrixType::Scalar Scalar;
70 typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
71
72 int rows = m.rows();
73 int cols = m.cols();
74
75 int r = ei_random<int>(0, rows-1),
76 c = ei_random<int>(0, cols-1);
77
Gael Guennebaudc10f0692008-07-21 00:34:46 +000078 MatrixType m1 = MatrixType::Random(rows, cols),
79 m2 = MatrixType::Random(rows, cols),
Gael Guennebaud8cef5412008-06-21 17:28:07 +000080 m3(rows, cols);
81
Benoit Jacobf5791ee2008-07-08 00:49:10 +000082 VERIFY(((m1.cwise() + Scalar(1)).cwise() > m1).all());
83 VERIFY(((m1.cwise() - Scalar(1)).cwise() < m1).all());
Gael Guennebaud8cef5412008-06-21 17:28:07 +000084 if (rows*cols>1)
85 {
86 m3 = m1;
87 m3(r,c) += 1;
Benoit Jacobf5791ee2008-07-08 00:49:10 +000088 VERIFY(! (m1.cwise() < m3).all() );
89 VERIFY(! (m1.cwise() > m3).all() );
Gael Guennebaud8cef5412008-06-21 17:28:07 +000090 }
Gael Guennebaud59dc1da2008-09-03 17:16:28 +000091
92 // test Select
93 VERIFY_IS_APPROX( (m1.cwise()<m2).select(m1,m2), m1.cwise().min(m2) );
94 VERIFY_IS_APPROX( (m1.cwise()>m2).select(m1,m2), m1.cwise().max(m2) );
95 Scalar mid = (m1.cwise().abs().minCoeff() + m1.cwise().abs().maxCoeff())/Scalar(2);
96 for (int j=0; j<cols; ++j)
97 for (int i=0; i<rows; ++i)
98 m3(i,j) = ei_abs(m1(i,j))<mid ? 0 : m1(i,j);
99 VERIFY_IS_APPROX( (m1.cwise().abs().cwise()<MatrixType::Constant(rows,cols,mid))
100 .select(MatrixType::Zero(rows,cols),m1), m3);
101 VERIFY_IS_APPROX( (m1.cwise().abs().cwise()<MatrixType::Constant(rows,cols,mid))
102 .select(0,m1), m3);
103 VERIFY_IS_APPROX( (m1.cwise().abs().cwise()>=MatrixType::Constant(rows,cols,mid))
104 .select(m1,0), m3);
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000105}
106
107void test_array()
108{
109 for(int i = 0; i < g_repeat; i++) {
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000110 CALL_SUBTEST( array(Matrix<float, 1, 1>()) );
111 CALL_SUBTEST( array(Matrix2f()) );
112 CALL_SUBTEST( array(Matrix4d()) );
113 CALL_SUBTEST( array(MatrixXcf(3, 3)) );
114 CALL_SUBTEST( array(MatrixXf(8, 12)) );
115 CALL_SUBTEST( array(MatrixXi(8, 12)) );
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000116 }
117 for(int i = 0; i < g_repeat; i++) {
118 CALL_SUBTEST( comparisons(Matrix<float, 1, 1>()) );
119 CALL_SUBTEST( comparisons(Matrix2f()) );
120 CALL_SUBTEST( comparisons(Matrix4d()) );
121 CALL_SUBTEST( comparisons(MatrixXf(8, 12)) );
122 CALL_SUBTEST( comparisons(MatrixXi(8, 12)) );
123 }
124}