linearstructure.cpp 5.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148
  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2006-2008 Benoit Jacob <jacob.benoit.1@gmail.com>
  5. // Copyright (C) 2014 Gael Guennebaud <gael.guennebaud@inria.fr>
  6. //
  7. // This Source Code Form is subject to the terms of the Mozilla
  8. // Public License v. 2.0. If a copy of the MPL was not distributed
  9. // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
  10. static bool g_called;
  11. #define EIGEN_SCALAR_BINARY_OP_PLUGIN { g_called |= (!internal::is_same<LhsScalar,RhsScalar>::value); }
  12. #include "main.h"
  13. template<typename MatrixType> void linearStructure(const MatrixType& m)
  14. {
  15. using std::abs;
  16. /* this test covers the following files:
  17. CwiseUnaryOp.h, CwiseBinaryOp.h, SelfCwiseBinaryOp.h
  18. */
  19. typedef typename MatrixType::Scalar Scalar;
  20. typedef typename MatrixType::RealScalar RealScalar;
  21. Index rows = m.rows();
  22. Index cols = m.cols();
  23. // this test relies a lot on Random.h, and there's not much more that we can do
  24. // to test it, hence I consider that we will have tested Random.h
  25. MatrixType m1 = MatrixType::Random(rows, cols),
  26. m2 = MatrixType::Random(rows, cols),
  27. m3(rows, cols);
  28. Scalar s1 = internal::random<Scalar>();
  29. while (abs(s1)<RealScalar(1e-3)) s1 = internal::random<Scalar>();
  30. Index r = internal::random<Index>(0, rows-1),
  31. c = internal::random<Index>(0, cols-1);
  32. VERIFY_IS_APPROX(-(-m1), m1);
  33. VERIFY_IS_APPROX(m1+m1, 2*m1);
  34. VERIFY_IS_APPROX(m1+m2-m1, m2);
  35. VERIFY_IS_APPROX(-m2+m1+m2, m1);
  36. VERIFY_IS_APPROX(m1*s1, s1*m1);
  37. VERIFY_IS_APPROX((m1+m2)*s1, s1*m1+s1*m2);
  38. VERIFY_IS_APPROX((-m1+m2)*s1, -s1*m1+s1*m2);
  39. m3 = m2; m3 += m1;
  40. VERIFY_IS_APPROX(m3, m1+m2);
  41. m3 = m2; m3 -= m1;
  42. VERIFY_IS_APPROX(m3, m2-m1);
  43. m3 = m2; m3 *= s1;
  44. VERIFY_IS_APPROX(m3, s1*m2);
  45. if(!NumTraits<Scalar>::IsInteger)
  46. {
  47. m3 = m2; m3 /= s1;
  48. VERIFY_IS_APPROX(m3, m2/s1);
  49. }
  50. // again, test operator() to check const-qualification
  51. VERIFY_IS_APPROX((-m1)(r,c), -(m1(r,c)));
  52. VERIFY_IS_APPROX((m1-m2)(r,c), (m1(r,c))-(m2(r,c)));
  53. VERIFY_IS_APPROX((m1+m2)(r,c), (m1(r,c))+(m2(r,c)));
  54. VERIFY_IS_APPROX((s1*m1)(r,c), s1*(m1(r,c)));
  55. VERIFY_IS_APPROX((m1*s1)(r,c), (m1(r,c))*s1);
  56. if(!NumTraits<Scalar>::IsInteger)
  57. VERIFY_IS_APPROX((m1/s1)(r,c), (m1(r,c))/s1);
  58. // use .block to disable vectorization and compare to the vectorized version
  59. VERIFY_IS_APPROX(m1+m1.block(0,0,rows,cols), m1+m1);
  60. VERIFY_IS_APPROX(m1.cwiseProduct(m1.block(0,0,rows,cols)), m1.cwiseProduct(m1));
  61. VERIFY_IS_APPROX(m1 - m1.block(0,0,rows,cols), m1 - m1);
  62. VERIFY_IS_APPROX(m1.block(0,0,rows,cols) * s1, m1 * s1);
  63. }
  64. // Make sure that complex * real and real * complex are properly optimized
  65. template<typename MatrixType> void real_complex(DenseIndex rows = MatrixType::RowsAtCompileTime, DenseIndex cols = MatrixType::ColsAtCompileTime)
  66. {
  67. typedef typename MatrixType::Scalar Scalar;
  68. typedef typename MatrixType::RealScalar RealScalar;
  69. RealScalar s = internal::random<RealScalar>();
  70. MatrixType m1 = MatrixType::Random(rows, cols);
  71. g_called = false;
  72. VERIFY_IS_APPROX(s*m1, Scalar(s)*m1);
  73. VERIFY(g_called && "real * matrix<complex> not properly optimized");
  74. g_called = false;
  75. VERIFY_IS_APPROX(m1*s, m1*Scalar(s));
  76. VERIFY(g_called && "matrix<complex> * real not properly optimized");
  77. g_called = false;
  78. VERIFY_IS_APPROX(m1/s, m1/Scalar(s));
  79. VERIFY(g_called && "matrix<complex> / real not properly optimized");
  80. g_called = false;
  81. VERIFY_IS_APPROX(s+m1.array(), Scalar(s)+m1.array());
  82. VERIFY(g_called && "real + matrix<complex> not properly optimized");
  83. g_called = false;
  84. VERIFY_IS_APPROX(m1.array()+s, m1.array()+Scalar(s));
  85. VERIFY(g_called && "matrix<complex> + real not properly optimized");
  86. g_called = false;
  87. VERIFY_IS_APPROX(s-m1.array(), Scalar(s)-m1.array());
  88. VERIFY(g_called && "real - matrix<complex> not properly optimized");
  89. g_called = false;
  90. VERIFY_IS_APPROX(m1.array()-s, m1.array()-Scalar(s));
  91. VERIFY(g_called && "matrix<complex> - real not properly optimized");
  92. }
  93. void test_linearstructure()
  94. {
  95. g_called = true;
  96. VERIFY(g_called); // avoid `unneeded-internal-declaration` warning.
  97. for(int i = 0; i < g_repeat; i++) {
  98. CALL_SUBTEST_1( linearStructure(Matrix<float, 1, 1>()) );
  99. CALL_SUBTEST_2( linearStructure(Matrix2f()) );
  100. CALL_SUBTEST_3( linearStructure(Vector3d()) );
  101. CALL_SUBTEST_4( linearStructure(Matrix4d()) );
  102. CALL_SUBTEST_5( linearStructure(MatrixXcf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE/2), internal::random<int>(1,EIGEN_TEST_MAX_SIZE/2))) );
  103. CALL_SUBTEST_6( linearStructure(MatrixXf (internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  104. CALL_SUBTEST_7( linearStructure(MatrixXi (internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  105. CALL_SUBTEST_8( linearStructure(MatrixXcd(internal::random<int>(1,EIGEN_TEST_MAX_SIZE/2), internal::random<int>(1,EIGEN_TEST_MAX_SIZE/2))) );
  106. CALL_SUBTEST_9( linearStructure(ArrayXXf (internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  107. CALL_SUBTEST_10( linearStructure(ArrayXXcf (internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  108. CALL_SUBTEST_11( real_complex<Matrix4cd>() );
  109. CALL_SUBTEST_11( real_complex<MatrixXcf>(10,10) );
  110. CALL_SUBTEST_11( real_complex<ArrayXXcf>(10,10) );
  111. }
  112. #ifdef EIGEN_TEST_PART_4
  113. {
  114. // make sure that /=scalar and /scalar do not overflow
  115. // rational: 1.0/4.94e-320 overflow, but m/4.94e-320 should not
  116. Matrix4d m2, m3;
  117. m3 = m2 = Matrix4d::Random()*1e-20;
  118. m2 = m2 / 4.9e-320;
  119. VERIFY_IS_APPROX(m2.cwiseQuotient(m2), Matrix4d::Ones());
  120. m3 /= 4.9e-320;
  121. VERIFY_IS_APPROX(m3.cwiseQuotient(m3), Matrix4d::Ones());
  122. }
  123. #endif
  124. }