xref: /aosp_15_r20/external/eigen/test/basicstuff.cpp (revision bf2c37156dfe67e5dfebd6d394bad8b2ab5804d4)
1*bf2c3715SXin Li // This file is part of Eigen, a lightweight C++ template library
2*bf2c3715SXin Li // for linear algebra.
3*bf2c3715SXin Li //
4*bf2c3715SXin Li // Copyright (C) 2006-2008 Benoit Jacob <[email protected]>
5*bf2c3715SXin Li //
6*bf2c3715SXin Li // This Source Code Form is subject to the terms of the Mozilla
7*bf2c3715SXin Li // Public License v. 2.0. If a copy of the MPL was not distributed
8*bf2c3715SXin Li // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9*bf2c3715SXin Li 
10*bf2c3715SXin Li #define EIGEN_NO_STATIC_ASSERT
11*bf2c3715SXin Li 
12*bf2c3715SXin Li #include "main.h"
13*bf2c3715SXin Li #include "random_without_cast_overflow.h"
14*bf2c3715SXin Li 
basicStuff(const MatrixType & m)15*bf2c3715SXin Li template<typename MatrixType> void basicStuff(const MatrixType& m)
16*bf2c3715SXin Li {
17*bf2c3715SXin Li   typedef typename MatrixType::Scalar Scalar;
18*bf2c3715SXin Li   typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
19*bf2c3715SXin Li   typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> SquareMatrixType;
20*bf2c3715SXin Li 
21*bf2c3715SXin Li   Index rows = m.rows();
22*bf2c3715SXin Li   Index cols = m.cols();
23*bf2c3715SXin Li 
24*bf2c3715SXin Li   // this test relies a lot on Random.h, and there's not much more that we can do
25*bf2c3715SXin Li   // to test it, hence I consider that we will have tested Random.h
26*bf2c3715SXin Li   MatrixType m1 = MatrixType::Random(rows, cols),
27*bf2c3715SXin Li              m2 = MatrixType::Random(rows, cols),
28*bf2c3715SXin Li              m3(rows, cols),
29*bf2c3715SXin Li              mzero = MatrixType::Zero(rows, cols),
30*bf2c3715SXin Li              square = Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime>::Random(rows, rows);
31*bf2c3715SXin Li   VectorType v1 = VectorType::Random(rows),
32*bf2c3715SXin Li              vzero = VectorType::Zero(rows);
33*bf2c3715SXin Li   SquareMatrixType sm1 = SquareMatrixType::Random(rows,rows), sm2(rows,rows);
34*bf2c3715SXin Li 
35*bf2c3715SXin Li   Scalar x = 0;
36*bf2c3715SXin Li   while(x == Scalar(0)) x = internal::random<Scalar>();
37*bf2c3715SXin Li 
38*bf2c3715SXin Li   Index r = internal::random<Index>(0, rows-1),
39*bf2c3715SXin Li         c = internal::random<Index>(0, cols-1);
40*bf2c3715SXin Li 
41*bf2c3715SXin Li   m1.coeffRef(r,c) = x;
42*bf2c3715SXin Li   VERIFY_IS_APPROX(x, m1.coeff(r,c));
43*bf2c3715SXin Li   m1(r,c) = x;
44*bf2c3715SXin Li   VERIFY_IS_APPROX(x, m1(r,c));
45*bf2c3715SXin Li   v1.coeffRef(r) = x;
46*bf2c3715SXin Li   VERIFY_IS_APPROX(x, v1.coeff(r));
47*bf2c3715SXin Li   v1(r) = x;
48*bf2c3715SXin Li   VERIFY_IS_APPROX(x, v1(r));
49*bf2c3715SXin Li   v1[r] = x;
50*bf2c3715SXin Li   VERIFY_IS_APPROX(x, v1[r]);
51*bf2c3715SXin Li 
52*bf2c3715SXin Li   // test fetching with various index types.
53*bf2c3715SXin Li   Index r1 = internal::random<Index>(0, numext::mini(Index(127),rows-1));
54*bf2c3715SXin Li   x = v1(static_cast<char>(r1));
55*bf2c3715SXin Li   x = v1(static_cast<signed char>(r1));
56*bf2c3715SXin Li   x = v1(static_cast<unsigned char>(r1));
57*bf2c3715SXin Li   x = v1(static_cast<signed short>(r1));
58*bf2c3715SXin Li   x = v1(static_cast<unsigned short>(r1));
59*bf2c3715SXin Li   x = v1(static_cast<signed int>(r1));
60*bf2c3715SXin Li   x = v1(static_cast<unsigned int>(r1));
61*bf2c3715SXin Li   x = v1(static_cast<signed long>(r1));
62*bf2c3715SXin Li   x = v1(static_cast<unsigned long>(r1));
63*bf2c3715SXin Li #if EIGEN_HAS_CXX11
64*bf2c3715SXin Li   x = v1(static_cast<long long int>(r1));
65*bf2c3715SXin Li   x = v1(static_cast<unsigned long long int>(r1));
66*bf2c3715SXin Li #endif
67*bf2c3715SXin Li 
68*bf2c3715SXin Li   VERIFY_IS_APPROX(               v1,    v1);
69*bf2c3715SXin Li   VERIFY_IS_NOT_APPROX(           v1,    2*v1);
70*bf2c3715SXin Li   VERIFY_IS_MUCH_SMALLER_THAN(    vzero, v1);
71*bf2c3715SXin Li   VERIFY_IS_MUCH_SMALLER_THAN(  vzero, v1.squaredNorm());
72*bf2c3715SXin Li   VERIFY_IS_NOT_MUCH_SMALLER_THAN(v1,    v1);
73*bf2c3715SXin Li   VERIFY_IS_APPROX(               vzero, v1-v1);
74*bf2c3715SXin Li   VERIFY_IS_APPROX(               m1,    m1);
75*bf2c3715SXin Li   VERIFY_IS_NOT_APPROX(           m1,    2*m1);
76*bf2c3715SXin Li   VERIFY_IS_MUCH_SMALLER_THAN(    mzero, m1);
77*bf2c3715SXin Li   VERIFY_IS_NOT_MUCH_SMALLER_THAN(m1,    m1);
78*bf2c3715SXin Li   VERIFY_IS_APPROX(               mzero, m1-m1);
79*bf2c3715SXin Li 
80*bf2c3715SXin Li   // always test operator() on each read-only expression class,
81*bf2c3715SXin Li   // in order to check const-qualifiers.
82*bf2c3715SXin Li   // indeed, if an expression class (here Zero) is meant to be read-only,
83*bf2c3715SXin Li   // hence has no _write() method, the corresponding MatrixBase method (here zero())
84*bf2c3715SXin Li   // should return a const-qualified object so that it is the const-qualified
85*bf2c3715SXin Li   // operator() that gets called, which in turn calls _read().
86*bf2c3715SXin Li   VERIFY_IS_MUCH_SMALLER_THAN(MatrixType::Zero(rows,cols)(r,c), static_cast<Scalar>(1));
87*bf2c3715SXin Li 
88*bf2c3715SXin Li   // now test copying a row-vector into a (column-)vector and conversely.
89*bf2c3715SXin Li   square.col(r) = square.row(r).eval();
90*bf2c3715SXin Li   Matrix<Scalar, 1, MatrixType::RowsAtCompileTime> rv(rows);
91*bf2c3715SXin Li   Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> cv(rows);
92*bf2c3715SXin Li   rv = square.row(r);
93*bf2c3715SXin Li   cv = square.col(r);
94*bf2c3715SXin Li 
95*bf2c3715SXin Li   VERIFY_IS_APPROX(rv, cv.transpose());
96*bf2c3715SXin Li 
97*bf2c3715SXin Li   if(cols!=1 && rows!=1 && MatrixType::SizeAtCompileTime!=Dynamic)
98*bf2c3715SXin Li   {
99*bf2c3715SXin Li     VERIFY_RAISES_ASSERT(m1 = (m2.block(0,0, rows-1, cols-1)));
100*bf2c3715SXin Li   }
101*bf2c3715SXin Li 
102*bf2c3715SXin Li   if(cols!=1 && rows!=1)
103*bf2c3715SXin Li   {
104*bf2c3715SXin Li     VERIFY_RAISES_ASSERT(m1[0]);
105*bf2c3715SXin Li     VERIFY_RAISES_ASSERT((m1+m1)[0]);
106*bf2c3715SXin Li   }
107*bf2c3715SXin Li 
108*bf2c3715SXin Li   VERIFY_IS_APPROX(m3 = m1,m1);
109*bf2c3715SXin Li   MatrixType m4;
110*bf2c3715SXin Li   VERIFY_IS_APPROX(m4 = m1,m1);
111*bf2c3715SXin Li 
112*bf2c3715SXin Li   m3.real() = m1.real();
113*bf2c3715SXin Li   VERIFY_IS_APPROX(static_cast<const MatrixType&>(m3).real(), static_cast<const MatrixType&>(m1).real());
114*bf2c3715SXin Li   VERIFY_IS_APPROX(static_cast<const MatrixType&>(m3).real(), m1.real());
115*bf2c3715SXin Li 
116*bf2c3715SXin Li   // check == / != operators
117*bf2c3715SXin Li   VERIFY(m1==m1);
118*bf2c3715SXin Li   VERIFY(m1!=m2);
119*bf2c3715SXin Li   VERIFY(!(m1==m2));
120*bf2c3715SXin Li   VERIFY(!(m1!=m1));
121*bf2c3715SXin Li   m1 = m2;
122*bf2c3715SXin Li   VERIFY(m1==m2);
123*bf2c3715SXin Li   VERIFY(!(m1!=m2));
124*bf2c3715SXin Li 
125*bf2c3715SXin Li   // check automatic transposition
126*bf2c3715SXin Li   sm2.setZero();
127*bf2c3715SXin Li   for(Index i=0;i<rows;++i)
128*bf2c3715SXin Li     sm2.col(i) = sm1.row(i);
129*bf2c3715SXin Li   VERIFY_IS_APPROX(sm2,sm1.transpose());
130*bf2c3715SXin Li 
131*bf2c3715SXin Li   sm2.setZero();
132*bf2c3715SXin Li   for(Index i=0;i<rows;++i)
133*bf2c3715SXin Li     sm2.col(i).noalias() = sm1.row(i);
134*bf2c3715SXin Li   VERIFY_IS_APPROX(sm2,sm1.transpose());
135*bf2c3715SXin Li 
136*bf2c3715SXin Li   sm2.setZero();
137*bf2c3715SXin Li   for(Index i=0;i<rows;++i)
138*bf2c3715SXin Li     sm2.col(i).noalias() += sm1.row(i);
139*bf2c3715SXin Li   VERIFY_IS_APPROX(sm2,sm1.transpose());
140*bf2c3715SXin Li 
141*bf2c3715SXin Li   sm2.setZero();
142*bf2c3715SXin Li   for(Index i=0;i<rows;++i)
143*bf2c3715SXin Li     sm2.col(i).noalias() -= sm1.row(i);
144*bf2c3715SXin Li   VERIFY_IS_APPROX(sm2,-sm1.transpose());
145*bf2c3715SXin Li 
146*bf2c3715SXin Li   // check ternary usage
147*bf2c3715SXin Li   {
148*bf2c3715SXin Li     bool b = internal::random<int>(0,10)>5;
149*bf2c3715SXin Li     m3 = b ? m1 : m2;
150*bf2c3715SXin Li     if(b) VERIFY_IS_APPROX(m3,m1);
151*bf2c3715SXin Li     else  VERIFY_IS_APPROX(m3,m2);
152*bf2c3715SXin Li     m3 = b ? -m1 : m2;
153*bf2c3715SXin Li     if(b) VERIFY_IS_APPROX(m3,-m1);
154*bf2c3715SXin Li     else  VERIFY_IS_APPROX(m3,m2);
155*bf2c3715SXin Li     m3 = b ? m1 : -m2;
156*bf2c3715SXin Li     if(b) VERIFY_IS_APPROX(m3,m1);
157*bf2c3715SXin Li     else  VERIFY_IS_APPROX(m3,-m2);
158*bf2c3715SXin Li   }
159*bf2c3715SXin Li }
160*bf2c3715SXin Li 
basicStuffComplex(const MatrixType & m)161*bf2c3715SXin Li template<typename MatrixType> void basicStuffComplex(const MatrixType& m)
162*bf2c3715SXin Li {
163*bf2c3715SXin Li   typedef typename MatrixType::Scalar Scalar;
164*bf2c3715SXin Li   typedef typename NumTraits<Scalar>::Real RealScalar;
165*bf2c3715SXin Li   typedef Matrix<RealScalar, MatrixType::RowsAtCompileTime, MatrixType::ColsAtCompileTime> RealMatrixType;
166*bf2c3715SXin Li 
167*bf2c3715SXin Li   Index rows = m.rows();
168*bf2c3715SXin Li   Index cols = m.cols();
169*bf2c3715SXin Li 
170*bf2c3715SXin Li   Scalar s1 = internal::random<Scalar>(),
171*bf2c3715SXin Li          s2 = internal::random<Scalar>();
172*bf2c3715SXin Li 
173*bf2c3715SXin Li   VERIFY(numext::real(s1)==numext::real_ref(s1));
174*bf2c3715SXin Li   VERIFY(numext::imag(s1)==numext::imag_ref(s1));
175*bf2c3715SXin Li   numext::real_ref(s1) = numext::real(s2);
176*bf2c3715SXin Li   numext::imag_ref(s1) = numext::imag(s2);
177*bf2c3715SXin Li   VERIFY(internal::isApprox(s1, s2, NumTraits<RealScalar>::epsilon()));
178*bf2c3715SXin Li   // extended precision in Intel FPUs means that s1 == s2 in the line above is not guaranteed.
179*bf2c3715SXin Li 
180*bf2c3715SXin Li   RealMatrixType rm1 = RealMatrixType::Random(rows,cols),
181*bf2c3715SXin Li                  rm2 = RealMatrixType::Random(rows,cols);
182*bf2c3715SXin Li   MatrixType cm(rows,cols);
183*bf2c3715SXin Li   cm.real() = rm1;
184*bf2c3715SXin Li   cm.imag() = rm2;
185*bf2c3715SXin Li   VERIFY_IS_APPROX(static_cast<const MatrixType&>(cm).real(), rm1);
186*bf2c3715SXin Li   VERIFY_IS_APPROX(static_cast<const MatrixType&>(cm).imag(), rm2);
187*bf2c3715SXin Li   rm1.setZero();
188*bf2c3715SXin Li   rm2.setZero();
189*bf2c3715SXin Li   rm1 = cm.real();
190*bf2c3715SXin Li   rm2 = cm.imag();
191*bf2c3715SXin Li   VERIFY_IS_APPROX(static_cast<const MatrixType&>(cm).real(), rm1);
192*bf2c3715SXin Li   VERIFY_IS_APPROX(static_cast<const MatrixType&>(cm).imag(), rm2);
193*bf2c3715SXin Li   cm.real().setZero();
194*bf2c3715SXin Li   VERIFY(static_cast<const MatrixType&>(cm).real().isZero());
195*bf2c3715SXin Li   VERIFY(!static_cast<const MatrixType&>(cm).imag().isZero());
196*bf2c3715SXin Li }
197*bf2c3715SXin Li 
198*bf2c3715SXin Li template<typename SrcScalar, typename TgtScalar>
199*bf2c3715SXin Li struct casting_test {
runcasting_test200*bf2c3715SXin Li   static void run() {
201*bf2c3715SXin Li     Matrix<SrcScalar,4,4> m;
202*bf2c3715SXin Li     for (int i=0; i<m.rows(); ++i) {
203*bf2c3715SXin Li       for (int j=0; j<m.cols(); ++j) {
204*bf2c3715SXin Li         m(i, j) = internal::random_without_cast_overflow<SrcScalar,TgtScalar>::value();
205*bf2c3715SXin Li       }
206*bf2c3715SXin Li     }
207*bf2c3715SXin Li     Matrix<TgtScalar,4,4> n = m.template cast<TgtScalar>();
208*bf2c3715SXin Li     for (int i=0; i<m.rows(); ++i) {
209*bf2c3715SXin Li       for (int j=0; j<m.cols(); ++j) {
210*bf2c3715SXin Li         VERIFY_IS_APPROX(n(i, j), (internal::cast<SrcScalar,TgtScalar>(m(i, j))));
211*bf2c3715SXin Li       }
212*bf2c3715SXin Li     }
213*bf2c3715SXin Li   }
214*bf2c3715SXin Li };
215*bf2c3715SXin Li 
216*bf2c3715SXin Li template<typename SrcScalar, typename EnableIf = void>
217*bf2c3715SXin Li struct casting_test_runner {
runcasting_test_runner218*bf2c3715SXin Li   static void run() {
219*bf2c3715SXin Li     casting_test<SrcScalar, bool>::run();
220*bf2c3715SXin Li     casting_test<SrcScalar, int8_t>::run();
221*bf2c3715SXin Li     casting_test<SrcScalar, uint8_t>::run();
222*bf2c3715SXin Li     casting_test<SrcScalar, int16_t>::run();
223*bf2c3715SXin Li     casting_test<SrcScalar, uint16_t>::run();
224*bf2c3715SXin Li     casting_test<SrcScalar, int32_t>::run();
225*bf2c3715SXin Li     casting_test<SrcScalar, uint32_t>::run();
226*bf2c3715SXin Li #if EIGEN_HAS_CXX11
227*bf2c3715SXin Li     casting_test<SrcScalar, int64_t>::run();
228*bf2c3715SXin Li     casting_test<SrcScalar, uint64_t>::run();
229*bf2c3715SXin Li #endif
230*bf2c3715SXin Li     casting_test<SrcScalar, half>::run();
231*bf2c3715SXin Li     casting_test<SrcScalar, bfloat16>::run();
232*bf2c3715SXin Li     casting_test<SrcScalar, float>::run();
233*bf2c3715SXin Li     casting_test<SrcScalar, double>::run();
234*bf2c3715SXin Li     casting_test<SrcScalar, std::complex<float> >::run();
235*bf2c3715SXin Li     casting_test<SrcScalar, std::complex<double> >::run();
236*bf2c3715SXin Li   }
237*bf2c3715SXin Li };
238*bf2c3715SXin Li 
239*bf2c3715SXin Li template<typename SrcScalar>
240*bf2c3715SXin Li struct casting_test_runner<SrcScalar, typename internal::enable_if<(NumTraits<SrcScalar>::IsComplex)>::type>
241*bf2c3715SXin Li {
runcasting_test_runner242*bf2c3715SXin Li   static void run() {
243*bf2c3715SXin Li     // Only a few casts from std::complex<T> are defined.
244*bf2c3715SXin Li     casting_test<SrcScalar, half>::run();
245*bf2c3715SXin Li     casting_test<SrcScalar, bfloat16>::run();
246*bf2c3715SXin Li     casting_test<SrcScalar, std::complex<float> >::run();
247*bf2c3715SXin Li     casting_test<SrcScalar, std::complex<double> >::run();
248*bf2c3715SXin Li   }
249*bf2c3715SXin Li };
250*bf2c3715SXin Li 
casting_all()251*bf2c3715SXin Li void casting_all() {
252*bf2c3715SXin Li   casting_test_runner<bool>::run();
253*bf2c3715SXin Li   casting_test_runner<int8_t>::run();
254*bf2c3715SXin Li   casting_test_runner<uint8_t>::run();
255*bf2c3715SXin Li   casting_test_runner<int16_t>::run();
256*bf2c3715SXin Li   casting_test_runner<uint16_t>::run();
257*bf2c3715SXin Li   casting_test_runner<int32_t>::run();
258*bf2c3715SXin Li   casting_test_runner<uint32_t>::run();
259*bf2c3715SXin Li #if EIGEN_HAS_CXX11
260*bf2c3715SXin Li   casting_test_runner<int64_t>::run();
261*bf2c3715SXin Li   casting_test_runner<uint64_t>::run();
262*bf2c3715SXin Li #endif
263*bf2c3715SXin Li   casting_test_runner<half>::run();
264*bf2c3715SXin Li   casting_test_runner<bfloat16>::run();
265*bf2c3715SXin Li   casting_test_runner<float>::run();
266*bf2c3715SXin Li   casting_test_runner<double>::run();
267*bf2c3715SXin Li   casting_test_runner<std::complex<float> >::run();
268*bf2c3715SXin Li   casting_test_runner<std::complex<double> >::run();
269*bf2c3715SXin Li }
270*bf2c3715SXin Li 
271*bf2c3715SXin Li template <typename Scalar>
fixedSizeMatrixConstruction()272*bf2c3715SXin Li void fixedSizeMatrixConstruction()
273*bf2c3715SXin Li {
274*bf2c3715SXin Li   Scalar raw[4];
275*bf2c3715SXin Li   for(int k=0; k<4; ++k)
276*bf2c3715SXin Li     raw[k] = internal::random<Scalar>();
277*bf2c3715SXin Li 
278*bf2c3715SXin Li   {
279*bf2c3715SXin Li     Matrix<Scalar,4,1> m(raw);
280*bf2c3715SXin Li     Array<Scalar,4,1> a(raw);
281*bf2c3715SXin Li     for(int k=0; k<4; ++k) VERIFY(m(k) == raw[k]);
282*bf2c3715SXin Li     for(int k=0; k<4; ++k) VERIFY(a(k) == raw[k]);
283*bf2c3715SXin Li     VERIFY_IS_EQUAL(m,(Matrix<Scalar,4,1>(raw[0],raw[1],raw[2],raw[3])));
284*bf2c3715SXin Li     VERIFY((a==(Array<Scalar,4,1>(raw[0],raw[1],raw[2],raw[3]))).all());
285*bf2c3715SXin Li   }
286*bf2c3715SXin Li   {
287*bf2c3715SXin Li     Matrix<Scalar,3,1> m(raw);
288*bf2c3715SXin Li     Array<Scalar,3,1> a(raw);
289*bf2c3715SXin Li     for(int k=0; k<3; ++k) VERIFY(m(k) == raw[k]);
290*bf2c3715SXin Li     for(int k=0; k<3; ++k) VERIFY(a(k) == raw[k]);
291*bf2c3715SXin Li     VERIFY_IS_EQUAL(m,(Matrix<Scalar,3,1>(raw[0],raw[1],raw[2])));
292*bf2c3715SXin Li     VERIFY((a==Array<Scalar,3,1>(raw[0],raw[1],raw[2])).all());
293*bf2c3715SXin Li   }
294*bf2c3715SXin Li   {
295*bf2c3715SXin Li     Matrix<Scalar,2,1> m(raw), m2( (DenseIndex(raw[0])), (DenseIndex(raw[1])) );
296*bf2c3715SXin Li     Array<Scalar,2,1> a(raw),  a2( (DenseIndex(raw[0])), (DenseIndex(raw[1])) );
297*bf2c3715SXin Li     for(int k=0; k<2; ++k) VERIFY(m(k) == raw[k]);
298*bf2c3715SXin Li     for(int k=0; k<2; ++k) VERIFY(a(k) == raw[k]);
299*bf2c3715SXin Li     VERIFY_IS_EQUAL(m,(Matrix<Scalar,2,1>(raw[0],raw[1])));
300*bf2c3715SXin Li     VERIFY((a==Array<Scalar,2,1>(raw[0],raw[1])).all());
301*bf2c3715SXin Li     for(int k=0; k<2; ++k) VERIFY(m2(k) == DenseIndex(raw[k]));
302*bf2c3715SXin Li     for(int k=0; k<2; ++k) VERIFY(a2(k) == DenseIndex(raw[k]));
303*bf2c3715SXin Li   }
304*bf2c3715SXin Li   {
305*bf2c3715SXin Li     Matrix<Scalar,1,2> m(raw),
306*bf2c3715SXin Li                        m2( (DenseIndex(raw[0])), (DenseIndex(raw[1])) ),
307*bf2c3715SXin Li                        m3( (int(raw[0])), (int(raw[1])) ),
308*bf2c3715SXin Li                        m4( (float(raw[0])), (float(raw[1])) );
309*bf2c3715SXin Li     Array<Scalar,1,2> a(raw),  a2( (DenseIndex(raw[0])), (DenseIndex(raw[1])) );
310*bf2c3715SXin Li     for(int k=0; k<2; ++k) VERIFY(m(k) == raw[k]);
311*bf2c3715SXin Li     for(int k=0; k<2; ++k) VERIFY(a(k) == raw[k]);
312*bf2c3715SXin Li     VERIFY_IS_EQUAL(m,(Matrix<Scalar,1,2>(raw[0],raw[1])));
313*bf2c3715SXin Li     VERIFY((a==Array<Scalar,1,2>(raw[0],raw[1])).all());
314*bf2c3715SXin Li     for(int k=0; k<2; ++k) VERIFY(m2(k) == DenseIndex(raw[k]));
315*bf2c3715SXin Li     for(int k=0; k<2; ++k) VERIFY(a2(k) == DenseIndex(raw[k]));
316*bf2c3715SXin Li     for(int k=0; k<2; ++k) VERIFY(m3(k) == int(raw[k]));
317*bf2c3715SXin Li     for(int k=0; k<2; ++k) VERIFY((m4(k)) == Scalar(float(raw[k])));
318*bf2c3715SXin Li   }
319*bf2c3715SXin Li   {
320*bf2c3715SXin Li     Matrix<Scalar,1,1> m(raw), m1(raw[0]), m2( (DenseIndex(raw[0])) ), m3( (int(raw[0])) );
321*bf2c3715SXin Li     Array<Scalar,1,1> a(raw), a1(raw[0]), a2( (DenseIndex(raw[0])) );
322*bf2c3715SXin Li     VERIFY(m(0) == raw[0]);
323*bf2c3715SXin Li     VERIFY(a(0) == raw[0]);
324*bf2c3715SXin Li     VERIFY(m1(0) == raw[0]);
325*bf2c3715SXin Li     VERIFY(a1(0) == raw[0]);
326*bf2c3715SXin Li     VERIFY(m2(0) == DenseIndex(raw[0]));
327*bf2c3715SXin Li     VERIFY(a2(0) == DenseIndex(raw[0]));
328*bf2c3715SXin Li     VERIFY(m3(0) == int(raw[0]));
329*bf2c3715SXin Li     VERIFY_IS_EQUAL(m,(Matrix<Scalar,1,1>(raw[0])));
330*bf2c3715SXin Li     VERIFY((a==Array<Scalar,1,1>(raw[0])).all());
331*bf2c3715SXin Li   }
332*bf2c3715SXin Li }
333*bf2c3715SXin Li 
EIGEN_DECLARE_TEST(basicstuff)334*bf2c3715SXin Li EIGEN_DECLARE_TEST(basicstuff)
335*bf2c3715SXin Li {
336*bf2c3715SXin Li   for(int i = 0; i < g_repeat; i++) {
337*bf2c3715SXin Li     CALL_SUBTEST_1( basicStuff(Matrix<float, 1, 1>()) );
338*bf2c3715SXin Li     CALL_SUBTEST_2( basicStuff(Matrix4d()) );
339*bf2c3715SXin Li     CALL_SUBTEST_3( basicStuff(MatrixXcf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
340*bf2c3715SXin Li     CALL_SUBTEST_4( basicStuff(MatrixXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
341*bf2c3715SXin Li     CALL_SUBTEST_5( basicStuff(MatrixXcd(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
342*bf2c3715SXin Li     CALL_SUBTEST_6( basicStuff(Matrix<float, 100, 100>()) );
343*bf2c3715SXin Li     CALL_SUBTEST_7( basicStuff(Matrix<long double,Dynamic,Dynamic>(internal::random<int>(1,EIGEN_TEST_MAX_SIZE),internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
344*bf2c3715SXin Li     CALL_SUBTEST_8( casting_all() );
345*bf2c3715SXin Li 
346*bf2c3715SXin Li     CALL_SUBTEST_3( basicStuffComplex(MatrixXcf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
347*bf2c3715SXin Li     CALL_SUBTEST_5( basicStuffComplex(MatrixXcd(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
348*bf2c3715SXin Li   }
349*bf2c3715SXin Li 
350*bf2c3715SXin Li   CALL_SUBTEST_1(fixedSizeMatrixConstruction<unsigned char>());
351*bf2c3715SXin Li   CALL_SUBTEST_1(fixedSizeMatrixConstruction<float>());
352*bf2c3715SXin Li   CALL_SUBTEST_1(fixedSizeMatrixConstruction<double>());
353*bf2c3715SXin Li   CALL_SUBTEST_1(fixedSizeMatrixConstruction<int>());
354*bf2c3715SXin Li   CALL_SUBTEST_1(fixedSizeMatrixConstruction<long int>());
355*bf2c3715SXin Li   CALL_SUBTEST_1(fixedSizeMatrixConstruction<std::ptrdiff_t>());
356*bf2c3715SXin Li }
357