• Home
  • Line#
  • Scopes#
  • Navigate#
  • Raw
  • Download
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2008 Benoit Jacob <jacob.benoit.1@gmail.com>
5 //
6 // This Source Code Form is subject to the terms of the Mozilla
7 // Public License v. 2.0. If a copy of the MPL was not distributed
8 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9 
10 #include "main.h"
11 
matrixVisitor(const MatrixType & p)12 template<typename MatrixType> void matrixVisitor(const MatrixType& p)
13 {
14   typedef typename MatrixType::Scalar Scalar;
15 
16   Index rows = p.rows();
17   Index cols = p.cols();
18 
19   // construct a random matrix where all coefficients are different
20   MatrixType m;
21   m = MatrixType::Random(rows, cols);
22   for(Index i = 0; i < m.size(); i++)
23     for(Index i2 = 0; i2 < i; i2++)
24       while(m(i) == m(i2)) // yes, ==
25         m(i) = internal::random<Scalar>();
26 
27   Scalar minc = Scalar(1000), maxc = Scalar(-1000);
28   Index minrow=0,mincol=0,maxrow=0,maxcol=0;
29   for(Index j = 0; j < cols; j++)
30   for(Index i = 0; i < rows; i++)
31   {
32     if(m(i,j) < minc)
33     {
34       minc = m(i,j);
35       minrow = i;
36       mincol = j;
37     }
38     if(m(i,j) > maxc)
39     {
40       maxc = m(i,j);
41       maxrow = i;
42       maxcol = j;
43     }
44   }
45   Index eigen_minrow, eigen_mincol, eigen_maxrow, eigen_maxcol;
46   Scalar eigen_minc, eigen_maxc;
47   eigen_minc = m.minCoeff(&eigen_minrow,&eigen_mincol);
48   eigen_maxc = m.maxCoeff(&eigen_maxrow,&eigen_maxcol);
49   VERIFY(minrow == eigen_minrow);
50   VERIFY(maxrow == eigen_maxrow);
51   VERIFY(mincol == eigen_mincol);
52   VERIFY(maxcol == eigen_maxcol);
53   VERIFY_IS_APPROX(minc, eigen_minc);
54   VERIFY_IS_APPROX(maxc, eigen_maxc);
55   VERIFY_IS_APPROX(minc, m.minCoeff());
56   VERIFY_IS_APPROX(maxc, m.maxCoeff());
57 
58   eigen_maxc = (m.adjoint()*m).maxCoeff(&eigen_maxrow,&eigen_maxcol);
59   Index maxrow2=0,maxcol2=0;
60   eigen_maxc = (m.adjoint()*m).eval().maxCoeff(&maxrow2,&maxcol2);
61   VERIFY(maxrow2 == eigen_maxrow);
62   VERIFY(maxcol2 == eigen_maxcol);
63 
64   if (!NumTraits<Scalar>::IsInteger && m.size() > 2) {
65     // Test NaN propagation by replacing an element with NaN.
66     bool stop = false;
67     for (Index j = 0; j < cols && !stop; ++j) {
68       for (Index i = 0; i < rows && !stop; ++i) {
69         if (!(j == mincol && i == minrow) &&
70             !(j == maxcol && i == maxrow)) {
71           m(i,j) = NumTraits<Scalar>::quiet_NaN();
72           stop = true;
73           break;
74         }
75       }
76     }
77 
78     eigen_minc = m.template minCoeff<PropagateNumbers>(&eigen_minrow, &eigen_mincol);
79     eigen_maxc = m.template maxCoeff<PropagateNumbers>(&eigen_maxrow, &eigen_maxcol);
80     VERIFY(minrow == eigen_minrow);
81     VERIFY(maxrow == eigen_maxrow);
82     VERIFY(mincol == eigen_mincol);
83     VERIFY(maxcol == eigen_maxcol);
84     VERIFY_IS_APPROX(minc, eigen_minc);
85     VERIFY_IS_APPROX(maxc, eigen_maxc);
86     VERIFY_IS_APPROX(minc, m.template minCoeff<PropagateNumbers>());
87     VERIFY_IS_APPROX(maxc, m.template maxCoeff<PropagateNumbers>());
88 
89     eigen_minc = m.template minCoeff<PropagateNaN>(&eigen_minrow, &eigen_mincol);
90     eigen_maxc = m.template maxCoeff<PropagateNaN>(&eigen_maxrow, &eigen_maxcol);
91     VERIFY(minrow != eigen_minrow || mincol != eigen_mincol);
92     VERIFY(maxrow != eigen_maxrow || maxcol != eigen_maxcol);
93     VERIFY((numext::isnan)(eigen_minc));
94     VERIFY((numext::isnan)(eigen_maxc));
95   }
96 
97 }
98 
vectorVisitor(const VectorType & w)99 template<typename VectorType> void vectorVisitor(const VectorType& w)
100 {
101   typedef typename VectorType::Scalar Scalar;
102 
103   Index size = w.size();
104 
105   // construct a random vector where all coefficients are different
106   VectorType v;
107   v = VectorType::Random(size);
108   for(Index i = 0; i < size; i++)
109     for(Index i2 = 0; i2 < i; i2++)
110       while(v(i) == v(i2)) // yes, ==
111         v(i) = internal::random<Scalar>();
112 
113   Scalar minc = v(0), maxc = v(0);
114   Index minidx=0, maxidx=0;
115   for(Index i = 0; i < size; i++)
116   {
117     if(v(i) < minc)
118     {
119       minc = v(i);
120       minidx = i;
121     }
122     if(v(i) > maxc)
123     {
124       maxc = v(i);
125       maxidx = i;
126     }
127   }
128   Index eigen_minidx, eigen_maxidx;
129   Scalar eigen_minc, eigen_maxc;
130   eigen_minc = v.minCoeff(&eigen_minidx);
131   eigen_maxc = v.maxCoeff(&eigen_maxidx);
132   VERIFY(minidx == eigen_minidx);
133   VERIFY(maxidx == eigen_maxidx);
134   VERIFY_IS_APPROX(minc, eigen_minc);
135   VERIFY_IS_APPROX(maxc, eigen_maxc);
136   VERIFY_IS_APPROX(minc, v.minCoeff());
137   VERIFY_IS_APPROX(maxc, v.maxCoeff());
138 
139   Index idx0 = internal::random<Index>(0,size-1);
140   Index idx1 = eigen_minidx;
141   Index idx2 = eigen_maxidx;
142   VectorType v1(v), v2(v);
143   v1(idx0) = v1(idx1);
144   v2(idx0) = v2(idx2);
145   v1.minCoeff(&eigen_minidx);
146   v2.maxCoeff(&eigen_maxidx);
147   VERIFY(eigen_minidx == (std::min)(idx0,idx1));
148   VERIFY(eigen_maxidx == (std::min)(idx0,idx2));
149 
150   if (!NumTraits<Scalar>::IsInteger && size > 2) {
151     // Test NaN propagation by replacing an element with NaN.
152     for (Index i = 0; i < size; ++i) {
153       if (i != minidx && i != maxidx) {
154         v(i) = NumTraits<Scalar>::quiet_NaN();
155         break;
156       }
157     }
158     eigen_minc = v.template minCoeff<PropagateNumbers>(&eigen_minidx);
159     eigen_maxc = v.template maxCoeff<PropagateNumbers>(&eigen_maxidx);
160     VERIFY(minidx == eigen_minidx);
161     VERIFY(maxidx == eigen_maxidx);
162     VERIFY_IS_APPROX(minc, eigen_minc);
163     VERIFY_IS_APPROX(maxc, eigen_maxc);
164     VERIFY_IS_APPROX(minc, v.template minCoeff<PropagateNumbers>());
165     VERIFY_IS_APPROX(maxc, v.template maxCoeff<PropagateNumbers>());
166 
167     eigen_minc = v.template minCoeff<PropagateNaN>(&eigen_minidx);
168     eigen_maxc = v.template maxCoeff<PropagateNaN>(&eigen_maxidx);
169     VERIFY(minidx != eigen_minidx);
170     VERIFY(maxidx != eigen_maxidx);
171     VERIFY((numext::isnan)(eigen_minc));
172     VERIFY((numext::isnan)(eigen_maxc));
173   }
174 }
175 
EIGEN_DECLARE_TEST(visitor)176 EIGEN_DECLARE_TEST(visitor)
177 {
178   for(int i = 0; i < g_repeat; i++) {
179     CALL_SUBTEST_1( matrixVisitor(Matrix<float, 1, 1>()) );
180     CALL_SUBTEST_2( matrixVisitor(Matrix2f()) );
181     CALL_SUBTEST_3( matrixVisitor(Matrix4d()) );
182     CALL_SUBTEST_4( matrixVisitor(MatrixXd(8, 12)) );
183     CALL_SUBTEST_5( matrixVisitor(Matrix<double,Dynamic,Dynamic,RowMajor>(20, 20)) );
184     CALL_SUBTEST_6( matrixVisitor(MatrixXi(8, 12)) );
185   }
186   for(int i = 0; i < g_repeat; i++) {
187     CALL_SUBTEST_7( vectorVisitor(Vector4f()) );
188     CALL_SUBTEST_7( vectorVisitor(Matrix<int,12,1>()) );
189     CALL_SUBTEST_8( vectorVisitor(VectorXd(10)) );
190     CALL_SUBTEST_9( vectorVisitor(RowVectorXd(10)) );
191     CALL_SUBTEST_10( vectorVisitor(VectorXf(33)) );
192   }
193 }
194