ref.cpp 11 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294
  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 20013 Gael Guennebaud <gael.guennebaud@inria.fr>
  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. // This unit test cannot be easily written to work with EIGEN_DEFAULT_TO_ROW_MAJOR
  10. #ifdef EIGEN_DEFAULT_TO_ROW_MAJOR
  11. #undef EIGEN_DEFAULT_TO_ROW_MAJOR
  12. #endif
  13. #define TEST_ENABLE_TEMPORARY_TRACKING
  14. #define TEST_CHECK_STATIC_ASSERTIONS
  15. #include "main.h"
  16. // test Ref.h
  17. // Deal with i387 extended precision
  18. #if EIGEN_ARCH_i386 && !(EIGEN_ARCH_x86_64)
  19. #if EIGEN_COMP_GNUC_STRICT && EIGEN_GNUC_AT_LEAST(4,4)
  20. #pragma GCC optimize ("-ffloat-store")
  21. #else
  22. #undef VERIFY_IS_EQUAL
  23. #define VERIFY_IS_EQUAL(X,Y) VERIFY_IS_APPROX(X,Y)
  24. #endif
  25. #endif
  26. template<typename MatrixType> void ref_matrix(const MatrixType& m)
  27. {
  28. typedef typename MatrixType::Scalar Scalar;
  29. typedef typename MatrixType::RealScalar RealScalar;
  30. typedef Matrix<Scalar,Dynamic,Dynamic,MatrixType::Options> DynMatrixType;
  31. typedef Matrix<RealScalar,Dynamic,Dynamic,MatrixType::Options> RealDynMatrixType;
  32. typedef Ref<MatrixType> RefMat;
  33. typedef Ref<DynMatrixType> RefDynMat;
  34. typedef Ref<const DynMatrixType> ConstRefDynMat;
  35. typedef Ref<RealDynMatrixType , 0, Stride<Dynamic,Dynamic> > RefRealMatWithStride;
  36. Index rows = m.rows(), cols = m.cols();
  37. MatrixType m1 = MatrixType::Random(rows, cols),
  38. m2 = m1;
  39. Index i = internal::random<Index>(0,rows-1);
  40. Index j = internal::random<Index>(0,cols-1);
  41. Index brows = internal::random<Index>(1,rows-i);
  42. Index bcols = internal::random<Index>(1,cols-j);
  43. RefMat rm0 = m1;
  44. VERIFY_IS_EQUAL(rm0, m1);
  45. RefDynMat rm1 = m1;
  46. VERIFY_IS_EQUAL(rm1, m1);
  47. RefDynMat rm2 = m1.block(i,j,brows,bcols);
  48. VERIFY_IS_EQUAL(rm2, m1.block(i,j,brows,bcols));
  49. rm2.setOnes();
  50. m2.block(i,j,brows,bcols).setOnes();
  51. VERIFY_IS_EQUAL(m1, m2);
  52. m2.block(i,j,brows,bcols).setRandom();
  53. rm2 = m2.block(i,j,brows,bcols);
  54. VERIFY_IS_EQUAL(m1, m2);
  55. ConstRefDynMat rm3 = m1.block(i,j,brows,bcols);
  56. m1.block(i,j,brows,bcols) *= 2;
  57. m2.block(i,j,brows,bcols) *= 2;
  58. VERIFY_IS_EQUAL(rm3, m2.block(i,j,brows,bcols));
  59. RefRealMatWithStride rm4 = m1.real();
  60. VERIFY_IS_EQUAL(rm4, m2.real());
  61. rm4.array() += 1;
  62. m2.real().array() += 1;
  63. VERIFY_IS_EQUAL(m1, m2);
  64. }
  65. template<typename VectorType> void ref_vector(const VectorType& m)
  66. {
  67. typedef typename VectorType::Scalar Scalar;
  68. typedef typename VectorType::RealScalar RealScalar;
  69. typedef Matrix<Scalar,Dynamic,1,VectorType::Options> DynMatrixType;
  70. typedef Matrix<Scalar,Dynamic,Dynamic,ColMajor> MatrixType;
  71. typedef Matrix<RealScalar,Dynamic,1,VectorType::Options> RealDynMatrixType;
  72. typedef Ref<VectorType> RefMat;
  73. typedef Ref<DynMatrixType> RefDynMat;
  74. typedef Ref<const DynMatrixType> ConstRefDynMat;
  75. typedef Ref<RealDynMatrixType , 0, InnerStride<> > RefRealMatWithStride;
  76. typedef Ref<DynMatrixType , 0, InnerStride<> > RefMatWithStride;
  77. Index size = m.size();
  78. VectorType v1 = VectorType::Random(size),
  79. v2 = v1;
  80. MatrixType mat1 = MatrixType::Random(size,size),
  81. mat2 = mat1,
  82. mat3 = MatrixType::Random(size,size);
  83. Index i = internal::random<Index>(0,size-1);
  84. Index bsize = internal::random<Index>(1,size-i);
  85. { RefMat rm0 = v1; VERIFY_IS_EQUAL(rm0, v1); }
  86. { RefMat rm0 = v1.block(0,0,size,1); VERIFY_IS_EQUAL(rm0, v1); }
  87. { RefDynMat rv1 = v1; VERIFY_IS_EQUAL(rv1, v1); }
  88. { RefDynMat rv1 = v1.block(0,0,size,1); VERIFY_IS_EQUAL(rv1, v1); }
  89. { VERIFY_RAISES_ASSERT( RefMat rm0 = v1.block(0, 0, size, 0); EIGEN_UNUSED_VARIABLE(rm0); ); }
  90. if(VectorType::SizeAtCompileTime!=1)
  91. { VERIFY_RAISES_ASSERT( RefDynMat rv1 = v1.block(0, 0, size, 0); EIGEN_UNUSED_VARIABLE(rv1); ); }
  92. RefDynMat rv2 = v1.segment(i,bsize);
  93. VERIFY_IS_EQUAL(rv2, v1.segment(i,bsize));
  94. rv2.setOnes();
  95. v2.segment(i,bsize).setOnes();
  96. VERIFY_IS_EQUAL(v1, v2);
  97. v2.segment(i,bsize).setRandom();
  98. rv2 = v2.segment(i,bsize);
  99. VERIFY_IS_EQUAL(v1, v2);
  100. ConstRefDynMat rm3 = v1.segment(i,bsize);
  101. v1.segment(i,bsize) *= 2;
  102. v2.segment(i,bsize) *= 2;
  103. VERIFY_IS_EQUAL(rm3, v2.segment(i,bsize));
  104. RefRealMatWithStride rm4 = v1.real();
  105. VERIFY_IS_EQUAL(rm4, v2.real());
  106. rm4.array() += 1;
  107. v2.real().array() += 1;
  108. VERIFY_IS_EQUAL(v1, v2);
  109. RefMatWithStride rm5 = mat1.row(i).transpose();
  110. VERIFY_IS_EQUAL(rm5, mat1.row(i).transpose());
  111. rm5.array() += 1;
  112. mat2.row(i).array() += 1;
  113. VERIFY_IS_EQUAL(mat1, mat2);
  114. rm5.noalias() = rm4.transpose() * mat3;
  115. mat2.row(i) = v2.real().transpose() * mat3;
  116. VERIFY_IS_APPROX(mat1, mat2);
  117. }
  118. template<typename PlainObjectType> void check_const_correctness(const PlainObjectType&)
  119. {
  120. // verify that ref-to-const don't have LvalueBit
  121. typedef typename internal::add_const<PlainObjectType>::type ConstPlainObjectType;
  122. VERIFY( !(internal::traits<Ref<ConstPlainObjectType> >::Flags & LvalueBit) );
  123. VERIFY( !(internal::traits<Ref<ConstPlainObjectType, Aligned> >::Flags & LvalueBit) );
  124. VERIFY( !(Ref<ConstPlainObjectType>::Flags & LvalueBit) );
  125. VERIFY( !(Ref<ConstPlainObjectType, Aligned>::Flags & LvalueBit) );
  126. }
  127. template<typename B>
  128. EIGEN_DONT_INLINE void call_ref_1(Ref<VectorXf> a, const B &b) { VERIFY_IS_EQUAL(a,b); }
  129. template<typename B>
  130. EIGEN_DONT_INLINE void call_ref_2(const Ref<const VectorXf>& a, const B &b) { VERIFY_IS_EQUAL(a,b); }
  131. template<typename B>
  132. EIGEN_DONT_INLINE void call_ref_3(Ref<VectorXf,0,InnerStride<> > a, const B &b) { VERIFY_IS_EQUAL(a,b); }
  133. template<typename B>
  134. EIGEN_DONT_INLINE void call_ref_4(const Ref<const VectorXf,0,InnerStride<> >& a, const B &b) { VERIFY_IS_EQUAL(a,b); }
  135. template<typename B>
  136. EIGEN_DONT_INLINE void call_ref_5(Ref<MatrixXf,0,OuterStride<> > a, const B &b) { VERIFY_IS_EQUAL(a,b); }
  137. template<typename B>
  138. EIGEN_DONT_INLINE void call_ref_6(const Ref<const MatrixXf,0,OuterStride<> >& a, const B &b) { VERIFY_IS_EQUAL(a,b); }
  139. template<typename B>
  140. EIGEN_DONT_INLINE void call_ref_7(Ref<Matrix<float,Dynamic,3> > a, const B &b) { VERIFY_IS_EQUAL(a,b); }
  141. void call_ref()
  142. {
  143. VectorXcf ca = VectorXcf::Random(10);
  144. VectorXf a = VectorXf::Random(10);
  145. RowVectorXf b = RowVectorXf::Random(10);
  146. MatrixXf A = MatrixXf::Random(10,10);
  147. RowVector3f c = RowVector3f::Random();
  148. const VectorXf& ac(a);
  149. VectorBlock<VectorXf> ab(a,0,3);
  150. const VectorBlock<VectorXf> abc(a,0,3);
  151. VERIFY_EVALUATION_COUNT( call_ref_1(a,a), 0);
  152. VERIFY_EVALUATION_COUNT( call_ref_1(b,b.transpose()), 0);
  153. // call_ref_1(ac,a<c); // does not compile because ac is const
  154. VERIFY_EVALUATION_COUNT( call_ref_1(ab,ab), 0);
  155. VERIFY_EVALUATION_COUNT( call_ref_1(a.head(4),a.head(4)), 0);
  156. VERIFY_EVALUATION_COUNT( call_ref_1(abc,abc), 0);
  157. VERIFY_EVALUATION_COUNT( call_ref_1(A.col(3),A.col(3)), 0);
  158. // call_ref_1(A.row(3),A.row(3)); // does not compile because innerstride!=1
  159. VERIFY_EVALUATION_COUNT( call_ref_3(A.row(3),A.row(3).transpose()), 0);
  160. VERIFY_EVALUATION_COUNT( call_ref_4(A.row(3),A.row(3).transpose()), 0);
  161. // call_ref_1(a+a, a+a); // does not compile for obvious reason
  162. MatrixXf tmp = A*A.col(1);
  163. VERIFY_EVALUATION_COUNT( call_ref_2(A*A.col(1), tmp), 1); // evaluated into a temp
  164. VERIFY_EVALUATION_COUNT( call_ref_2(ac.head(5),ac.head(5)), 0);
  165. VERIFY_EVALUATION_COUNT( call_ref_2(ac,ac), 0);
  166. VERIFY_EVALUATION_COUNT( call_ref_2(a,a), 0);
  167. VERIFY_EVALUATION_COUNT( call_ref_2(ab,ab), 0);
  168. VERIFY_EVALUATION_COUNT( call_ref_2(a.head(4),a.head(4)), 0);
  169. tmp = a+a;
  170. VERIFY_EVALUATION_COUNT( call_ref_2(a+a,tmp), 1); // evaluated into a temp
  171. VERIFY_EVALUATION_COUNT( call_ref_2(ca.imag(),ca.imag()), 1); // evaluated into a temp
  172. VERIFY_EVALUATION_COUNT( call_ref_4(ac.head(5),ac.head(5)), 0);
  173. tmp = a+a;
  174. VERIFY_EVALUATION_COUNT( call_ref_4(a+a,tmp), 1); // evaluated into a temp
  175. VERIFY_EVALUATION_COUNT( call_ref_4(ca.imag(),ca.imag()), 0);
  176. VERIFY_EVALUATION_COUNT( call_ref_5(a,a), 0);
  177. VERIFY_EVALUATION_COUNT( call_ref_5(a.head(3),a.head(3)), 0);
  178. VERIFY_EVALUATION_COUNT( call_ref_5(A,A), 0);
  179. // call_ref_5(A.transpose(),A.transpose()); // does not compile because storage order does not match
  180. VERIFY_EVALUATION_COUNT( call_ref_5(A.block(1,1,2,2),A.block(1,1,2,2)), 0);
  181. VERIFY_EVALUATION_COUNT( call_ref_5(b,b), 0); // storage order do not match, but this is a degenerate case that should work
  182. VERIFY_EVALUATION_COUNT( call_ref_5(a.row(3),a.row(3)), 0);
  183. VERIFY_EVALUATION_COUNT( call_ref_6(a,a), 0);
  184. VERIFY_EVALUATION_COUNT( call_ref_6(a.head(3),a.head(3)), 0);
  185. VERIFY_EVALUATION_COUNT( call_ref_6(A.row(3),A.row(3)), 1); // evaluated into a temp thouth it could be avoided by viewing it as a 1xn matrix
  186. tmp = A+A;
  187. VERIFY_EVALUATION_COUNT( call_ref_6(A+A,tmp), 1); // evaluated into a temp
  188. VERIFY_EVALUATION_COUNT( call_ref_6(A,A), 0);
  189. VERIFY_EVALUATION_COUNT( call_ref_6(A.transpose(),A.transpose()), 1); // evaluated into a temp because the storage orders do not match
  190. VERIFY_EVALUATION_COUNT( call_ref_6(A.block(1,1,2,2),A.block(1,1,2,2)), 0);
  191. VERIFY_EVALUATION_COUNT( call_ref_7(c,c), 0);
  192. }
  193. typedef Matrix<double,Dynamic,Dynamic,RowMajor> RowMatrixXd;
  194. int test_ref_overload_fun1(Ref<MatrixXd> ) { return 1; }
  195. int test_ref_overload_fun1(Ref<RowMatrixXd> ) { return 2; }
  196. int test_ref_overload_fun1(Ref<MatrixXf> ) { return 3; }
  197. int test_ref_overload_fun2(Ref<const MatrixXd> ) { return 4; }
  198. int test_ref_overload_fun2(Ref<const MatrixXf> ) { return 5; }
  199. void test_ref_ambiguous(const Ref<const ArrayXd> &A, Ref<ArrayXd> B)
  200. {
  201. B = A;
  202. B = A - A;
  203. }
  204. // See also bug 969
  205. void test_ref_overloads()
  206. {
  207. MatrixXd Ad, Bd;
  208. RowMatrixXd rAd, rBd;
  209. VERIFY( test_ref_overload_fun1(Ad)==1 );
  210. VERIFY( test_ref_overload_fun1(rAd)==2 );
  211. MatrixXf Af, Bf;
  212. VERIFY( test_ref_overload_fun2(Ad)==4 );
  213. VERIFY( test_ref_overload_fun2(Ad+Bd)==4 );
  214. VERIFY( test_ref_overload_fun2(Af+Bf)==5 );
  215. ArrayXd A, B;
  216. test_ref_ambiguous(A, B);
  217. }
  218. void test_ref_fixed_size_assert()
  219. {
  220. Vector4f v4;
  221. VectorXf vx(10);
  222. VERIFY_RAISES_STATIC_ASSERT( Ref<Vector3f> y = v4; (void)y; );
  223. VERIFY_RAISES_STATIC_ASSERT( Ref<Vector3f> y = vx.head<4>(); (void)y; );
  224. VERIFY_RAISES_STATIC_ASSERT( Ref<const Vector3f> y = v4; (void)y; );
  225. VERIFY_RAISES_STATIC_ASSERT( Ref<const Vector3f> y = vx.head<4>(); (void)y; );
  226. VERIFY_RAISES_STATIC_ASSERT( Ref<const Vector3f> y = 2*v4; (void)y; );
  227. }
  228. void test_ref()
  229. {
  230. for(int i = 0; i < g_repeat; i++) {
  231. CALL_SUBTEST_1( ref_vector(Matrix<float, 1, 1>()) );
  232. CALL_SUBTEST_1( check_const_correctness(Matrix<float, 1, 1>()) );
  233. CALL_SUBTEST_2( ref_vector(Vector4d()) );
  234. CALL_SUBTEST_2( check_const_correctness(Matrix4d()) );
  235. CALL_SUBTEST_3( ref_vector(Vector4cf()) );
  236. CALL_SUBTEST_4( ref_vector(VectorXcf(8)) );
  237. CALL_SUBTEST_5( ref_vector(VectorXi(12)) );
  238. CALL_SUBTEST_5( check_const_correctness(VectorXi(12)) );
  239. CALL_SUBTEST_1( ref_matrix(Matrix<float, 1, 1>()) );
  240. CALL_SUBTEST_2( ref_matrix(Matrix4d()) );
  241. CALL_SUBTEST_1( ref_matrix(Matrix<float,3,5>()) );
  242. CALL_SUBTEST_4( ref_matrix(MatrixXcf(internal::random<int>(1,10),internal::random<int>(1,10))) );
  243. CALL_SUBTEST_4( ref_matrix(Matrix<std::complex<double>,10,15>()) );
  244. CALL_SUBTEST_5( ref_matrix(MatrixXi(internal::random<int>(1,10),internal::random<int>(1,10))) );
  245. CALL_SUBTEST_6( call_ref() );
  246. }
  247. CALL_SUBTEST_7( test_ref_overloads() );
  248. CALL_SUBTEST_7( test_ref_fixed_size_assert() );
  249. }