1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
|
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra. Eigen itself is part of the KDE project.
//
// Copyright (C) 2008 Gael Guennebaud <g.gael@free.fr>
//
// Eigen is free software; you can redistribute it and/or
// modify it under the terms of the GNU Lesser General Public
// License as published by the Free Software Foundation; either
// version 3 of the License, or (at your option) any later version.
//
// Alternatively, you can redistribute it and/or
// modify it under the terms of the GNU General Public License as
// published by the Free Software Foundation; either version 2 of
// the License, or (at your option) any later version.
//
// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
// GNU General Public License for more details.
//
// You should have received a copy of the GNU Lesser General Public
// License and a copy of the GNU General Public License along with
// Eigen. If not, see <http://www.gnu.org/licenses/>.
#include "main.h"
#include <Eigen/Array>
template<typename MatrixType> void array(const MatrixType& m)
{
/* this test covers the following files:
Array.cpp
*/
typedef typename MatrixType::Scalar Scalar;
typedef typename NumTraits<Scalar>::Real RealScalar;
typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
int rows = m.rows();
int cols = m.cols();
MatrixType m1 = MatrixType::Random(rows, cols),
m2 = MatrixType::Random(rows, cols),
m3(rows, cols);
Scalar s1 = ei_random<Scalar>(),
s2 = ei_random<Scalar>();
// scalar addition
VERIFY_IS_APPROX(m1.cwise() + s1, s1 + m1.cwise());
VERIFY_IS_APPROX(m1.cwise() + s1, MatrixType::Constant(rows,cols,s1) + m1);
VERIFY_IS_APPROX((m1*Scalar(2)).cwise() - s2, (m1+m1) - MatrixType::Constant(rows,cols,s2) );
m3 = m1;
m3.cwise() += s2;
VERIFY_IS_APPROX(m3, m1.cwise() + s2);
m3 = m1;
m3.cwise() -= s1;
VERIFY_IS_APPROX(m3, m1.cwise() - s1);
// reductions
VERIFY_IS_APPROX(m1.colwise().sum().sum(), m1.sum());
VERIFY_IS_APPROX(m1.rowwise().sum().sum(), m1.sum());
if (!ei_isApprox(m1.sum(), (m1+m2).sum()))
VERIFY_IS_NOT_APPROX(((m1+m2).rowwise().sum()).sum(), m1.sum());
VERIFY_IS_APPROX(m1.colwise().sum(), m1.colwise().redux(ei_scalar_sum_op<Scalar>()));
}
template<typename MatrixType> void comparisons(const MatrixType& m)
{
typedef typename MatrixType::Scalar Scalar;
typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
int rows = m.rows();
int cols = m.cols();
int r = ei_random<int>(0, rows-1),
c = ei_random<int>(0, cols-1);
MatrixType m1 = MatrixType::Random(rows, cols),
m2 = MatrixType::Random(rows, cols),
m3(rows, cols);
VERIFY(((m1.cwise() + Scalar(1)).cwise() > m1).all());
VERIFY(((m1.cwise() - Scalar(1)).cwise() < m1).all());
if (rows*cols>1)
{
m3 = m1;
m3(r,c) += 1;
VERIFY(! (m1.cwise() < m3).all() );
VERIFY(! (m1.cwise() > m3).all() );
}
// comparisons to scalar
VERIFY( (m1.cwise() != (m1(r,c)+1) ).any() );
VERIFY( (m1.cwise() > (m1(r,c)-1) ).any() );
VERIFY( (m1.cwise() < (m1(r,c)+1) ).any() );
VERIFY( (m1.cwise() == m1(r,c) ).any() );
// test Select
VERIFY_IS_APPROX( (m1.cwise()<m2).select(m1,m2), m1.cwise().min(m2) );
VERIFY_IS_APPROX( (m1.cwise()>m2).select(m1,m2), m1.cwise().max(m2) );
Scalar mid = (m1.cwise().abs().minCoeff() + m1.cwise().abs().maxCoeff())/Scalar(2);
for (int j=0; j<cols; ++j)
for (int i=0; i<rows; ++i)
m3(i,j) = ei_abs(m1(i,j))<mid ? 0 : m1(i,j);
VERIFY_IS_APPROX( (m1.cwise().abs().cwise()<MatrixType::Constant(rows,cols,mid))
.select(MatrixType::Zero(rows,cols),m1), m3);
// shorter versions:
VERIFY_IS_APPROX( (m1.cwise().abs().cwise()<MatrixType::Constant(rows,cols,mid))
.select(0,m1), m3);
VERIFY_IS_APPROX( (m1.cwise().abs().cwise()>=MatrixType::Constant(rows,cols,mid))
.select(m1,0), m3);
// even shorter version:
VERIFY_IS_APPROX( (m1.cwise().abs().cwise()<mid).select(0,m1), m3);
}
template<typename VectorType> void lpNorm(const VectorType& v)
{
VectorType u = VectorType::Random(v.size());
VERIFY_IS_APPROX(u.template lpNorm<Infinity>(), u.cwise().abs().maxCoeff());
VERIFY_IS_APPROX(u.template lpNorm<1>(), u.cwise().abs().sum());
VERIFY_IS_APPROX(u.template lpNorm<2>(), ei_sqrt(u.cwise().abs().cwise().square().sum()));
VERIFY_IS_APPROX(ei_pow(u.template lpNorm<5>(), typename VectorType::RealScalar(5)), u.cwise().abs().cwise().pow(5).sum());
}
void test_array()
{
for(int i = 0; i < g_repeat; i++) {
CALL_SUBTEST( array(Matrix<float, 1, 1>()) );
CALL_SUBTEST( array(Matrix2f()) );
CALL_SUBTEST( array(Matrix4d()) );
CALL_SUBTEST( array(MatrixXcf(3, 3)) );
CALL_SUBTEST( array(MatrixXf(8, 12)) );
CALL_SUBTEST( array(MatrixXi(8, 12)) );
}
for(int i = 0; i < g_repeat; i++) {
CALL_SUBTEST( comparisons(Matrix<float, 1, 1>()) );
CALL_SUBTEST( comparisons(Matrix2f()) );
CALL_SUBTEST( comparisons(Matrix4d()) );
CALL_SUBTEST( comparisons(MatrixXf(8, 12)) );
CALL_SUBTEST( comparisons(MatrixXi(8, 12)) );
}
for(int i = 0; i < g_repeat; i++) {
CALL_SUBTEST( lpNorm(Matrix<float, 1, 1>()) );
CALL_SUBTEST( lpNorm(Vector2f()) );
CALL_SUBTEST( lpNorm(Vector3d()) );
CALL_SUBTEST( lpNorm(Vector4f()) );
CALL_SUBTEST( lpNorm(VectorXf(16)) );
CALL_SUBTEST( lpNorm(VectorXcd(10)) );
}
}
|