Gael Guennebaud | 8cef541 | 2008-06-21 17:28:07 +0000 | [diff] [blame] | 1 | // 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 | |
| 28 | template<typename MatrixType> void scalarAdd(const MatrixType& m) |
| 29 | { |
| 30 | /* this test covers the following files: |
| 31 | Array.cpp |
| 32 | */ |
| 33 | |
| 34 | typedef typename MatrixType::Scalar Scalar; |
Gael Guennebaud | 2120fed | 2008-08-23 15:14:20 +0000 | [diff] [blame] | 35 | typedef typename NumTraits<Scalar>::Real RealScalar; |
Gael Guennebaud | 8cef541 | 2008-06-21 17:28:07 +0000 | [diff] [blame] | 36 | typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType; |
| 37 | |
| 38 | int rows = m.rows(); |
| 39 | int cols = m.cols(); |
| 40 | |
Gael Guennebaud | 2120fed | 2008-08-23 15:14:20 +0000 | [diff] [blame] | 41 | MatrixType m1 = test_random_matrix<MatrixType>(rows, cols), |
| 42 | m2 = test_random_matrix<MatrixType>(rows, cols), |
Gael Guennebaud | 8cef541 | 2008-06-21 17:28:07 +0000 | [diff] [blame] | 43 | m3(rows, cols); |
| 44 | |
Gael Guennebaud | 2120fed | 2008-08-23 15:14:20 +0000 | [diff] [blame] | 45 | Scalar s1 = test_random<Scalar>(), |
| 46 | s2 = test_random<Scalar>(); |
Gael Guennebaud | 8cef541 | 2008-06-21 17:28:07 +0000 | [diff] [blame] | 47 | |
Benoit Jacob | f5791ee | 2008-07-08 00:49:10 +0000 | [diff] [blame] | 48 | VERIFY_IS_APPROX(m1.cwise() + s1, s1 + m1.cwise()); |
Gael Guennebaud | c10f069 | 2008-07-21 00:34:46 +0000 | [diff] [blame] | 49 | VERIFY_IS_APPROX(m1.cwise() + s1, MatrixType::Constant(rows,cols,s1) + m1); |
| 50 | VERIFY_IS_APPROX((m1*Scalar(2)).cwise() - s2, (m1+m1) - MatrixType::Constant(rows,cols,s2) ); |
Gael Guennebaud | 8cef541 | 2008-06-21 17:28:07 +0000 | [diff] [blame] | 51 | m3 = m1; |
Benoit Jacob | f5791ee | 2008-07-08 00:49:10 +0000 | [diff] [blame] | 52 | m3.cwise() += s2; |
| 53 | VERIFY_IS_APPROX(m3, m1.cwise() + s2); |
Gael Guennebaud | 8cef541 | 2008-06-21 17:28:07 +0000 | [diff] [blame] | 54 | m3 = m1; |
Benoit Jacob | f5791ee | 2008-07-08 00:49:10 +0000 | [diff] [blame] | 55 | m3.cwise() -= s1; |
| 56 | VERIFY_IS_APPROX(m3, m1.cwise() - s1); |
Gael Guennebaud | 05ad083 | 2008-07-19 13:03:23 +0000 | [diff] [blame] | 57 | |
| 58 | VERIFY_IS_APPROX(m1.colwise().sum().sum(), m1.sum()); |
| 59 | VERIFY_IS_APPROX(m1.rowwise().sum().sum(), m1.sum()); |
Gael Guennebaud | 2120fed | 2008-08-23 15:14:20 +0000 | [diff] [blame] | 60 | if (!ei_isApprox(m1.sum(), (m1+m2).sum())) |
| 61 | VERIFY_IS_NOT_APPROX(((m1+m2).rowwise().sum()).sum(), m1.sum()); |
Gael Guennebaud | 05ad083 | 2008-07-19 13:03:23 +0000 | [diff] [blame] | 62 | VERIFY_IS_APPROX(m1.colwise().sum(), m1.colwise().redux(ei_scalar_sum_op<Scalar>())); |
Gael Guennebaud | 8cef541 | 2008-06-21 17:28:07 +0000 | [diff] [blame] | 63 | } |
| 64 | |
| 65 | template<typename MatrixType> void comparisons(const MatrixType& m) |
| 66 | { |
| 67 | typedef typename MatrixType::Scalar Scalar; |
| 68 | typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType; |
| 69 | |
| 70 | int rows = m.rows(); |
| 71 | int cols = m.cols(); |
| 72 | |
| 73 | int r = ei_random<int>(0, rows-1), |
| 74 | c = ei_random<int>(0, cols-1); |
| 75 | |
Gael Guennebaud | c10f069 | 2008-07-21 00:34:46 +0000 | [diff] [blame] | 76 | MatrixType m1 = MatrixType::Random(rows, cols), |
| 77 | m2 = MatrixType::Random(rows, cols), |
Gael Guennebaud | 8cef541 | 2008-06-21 17:28:07 +0000 | [diff] [blame] | 78 | m3(rows, cols); |
| 79 | |
Benoit Jacob | f5791ee | 2008-07-08 00:49:10 +0000 | [diff] [blame] | 80 | VERIFY(((m1.cwise() + Scalar(1)).cwise() > m1).all()); |
| 81 | VERIFY(((m1.cwise() - Scalar(1)).cwise() < m1).all()); |
Gael Guennebaud | 8cef541 | 2008-06-21 17:28:07 +0000 | [diff] [blame] | 82 | if (rows*cols>1) |
| 83 | { |
| 84 | m3 = m1; |
| 85 | m3(r,c) += 1; |
Benoit Jacob | f5791ee | 2008-07-08 00:49:10 +0000 | [diff] [blame] | 86 | VERIFY(! (m1.cwise() < m3).all() ); |
| 87 | VERIFY(! (m1.cwise() > m3).all() ); |
Gael Guennebaud | 8cef541 | 2008-06-21 17:28:07 +0000 | [diff] [blame] | 88 | } |
| 89 | } |
| 90 | |
| 91 | void test_array() |
| 92 | { |
| 93 | for(int i = 0; i < g_repeat; i++) { |
| 94 | CALL_SUBTEST( scalarAdd(Matrix<float, 1, 1>()) ); |
| 95 | CALL_SUBTEST( scalarAdd(Matrix2f()) ); |
| 96 | CALL_SUBTEST( scalarAdd(Matrix4d()) ); |
| 97 | CALL_SUBTEST( scalarAdd(MatrixXcf(3, 3)) ); |
| 98 | CALL_SUBTEST( scalarAdd(MatrixXf(8, 12)) ); |
| 99 | CALL_SUBTEST( scalarAdd(MatrixXi(8, 12)) ); |
| 100 | } |
| 101 | for(int i = 0; i < g_repeat; i++) { |
| 102 | CALL_SUBTEST( comparisons(Matrix<float, 1, 1>()) ); |
| 103 | CALL_SUBTEST( comparisons(Matrix2f()) ); |
| 104 | CALL_SUBTEST( comparisons(Matrix4d()) ); |
| 105 | CALL_SUBTEST( comparisons(MatrixXf(8, 12)) ); |
| 106 | CALL_SUBTEST( comparisons(MatrixXi(8, 12)) ); |
| 107 | } |
| 108 | } |