You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.

107 lines
4.3 KiB

  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2010 Hauke Heibel <hauke.heibel@gmail.com>
  5. // Copyright (C) 2015 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. #define TEST_ENABLE_TEMPORARY_TRACKING
  11. #include "main.h"
  12. template <int N, typename XprType>
  13. void use_n_times(const XprType &xpr)
  14. {
  15. typename internal::nested_eval<XprType,N>::type mat(xpr);
  16. typename XprType::PlainObject res(mat.rows(), mat.cols());
  17. nb_temporaries--; // remove res
  18. res.setZero();
  19. for(int i=0; i<N; ++i)
  20. res += mat;
  21. }
  22. template <int N, typename ReferenceType, typename XprType>
  23. bool verify_eval_type(const XprType &, const ReferenceType&)
  24. {
  25. typedef typename internal::nested_eval<XprType,N>::type EvalType;
  26. return internal::is_same<typename internal::remove_all<EvalType>::type, typename internal::remove_all<ReferenceType>::type>::value;
  27. }
  28. template <typename MatrixType> void run_nesting_ops_1(const MatrixType& _m)
  29. {
  30. typename internal::nested_eval<MatrixType,2>::type m(_m);
  31. // Make really sure that we are in debug mode!
  32. VERIFY_RAISES_ASSERT(eigen_assert(false));
  33. // The only intention of these tests is to ensure that this code does
  34. // not trigger any asserts or segmentation faults... more to come.
  35. VERIFY_IS_APPROX( (m.transpose() * m).diagonal().sum(), (m.transpose() * m).diagonal().sum() );
  36. VERIFY_IS_APPROX( (m.transpose() * m).diagonal().array().abs().sum(), (m.transpose() * m).diagonal().array().abs().sum() );
  37. VERIFY_IS_APPROX( (m.transpose() * m).array().abs().sum(), (m.transpose() * m).array().abs().sum() );
  38. }
  39. template <typename MatrixType> void run_nesting_ops_2(const MatrixType& _m)
  40. {
  41. typedef typename MatrixType::Scalar Scalar;
  42. Index rows = _m.rows();
  43. Index cols = _m.cols();
  44. MatrixType m1 = MatrixType::Random(rows,cols);
  45. Matrix<Scalar,MatrixType::RowsAtCompileTime,MatrixType::ColsAtCompileTime,ColMajor> m2;
  46. if((MatrixType::SizeAtCompileTime==Dynamic))
  47. {
  48. VERIFY_EVALUATION_COUNT( use_n_times<1>(m1 + m1*m1), 1 );
  49. VERIFY_EVALUATION_COUNT( use_n_times<10>(m1 + m1*m1), 1 );
  50. VERIFY_EVALUATION_COUNT( use_n_times<1>(m1.template triangularView<Lower>().solve(m1.col(0))), 1 );
  51. VERIFY_EVALUATION_COUNT( use_n_times<10>(m1.template triangularView<Lower>().solve(m1.col(0))), 1 );
  52. VERIFY_EVALUATION_COUNT( use_n_times<1>(Scalar(2)*m1.template triangularView<Lower>().solve(m1.col(0))), 2 ); // FIXME could be one by applying the scaling in-place on the solve result
  53. VERIFY_EVALUATION_COUNT( use_n_times<1>(m1.col(0)+m1.template triangularView<Lower>().solve(m1.col(0))), 2 ); // FIXME could be one by adding m1.col() inplace
  54. VERIFY_EVALUATION_COUNT( use_n_times<10>(m1.col(0)+m1.template triangularView<Lower>().solve(m1.col(0))), 2 );
  55. }
  56. {
  57. VERIFY( verify_eval_type<10>(m1, m1) );
  58. if(!NumTraits<Scalar>::IsComplex)
  59. {
  60. VERIFY( verify_eval_type<3>(2*m1, 2*m1) );
  61. VERIFY( verify_eval_type<4>(2*m1, m1) );
  62. }
  63. else
  64. {
  65. VERIFY( verify_eval_type<1>(2*m1, 2*m1) );
  66. VERIFY( verify_eval_type<2>(2*m1, m1) );
  67. }
  68. VERIFY( verify_eval_type<2>(m1+m1, m1+m1) );
  69. VERIFY( verify_eval_type<3>(m1+m1, m1) );
  70. VERIFY( verify_eval_type<1>(m1*m1.transpose(), m2) );
  71. VERIFY( verify_eval_type<1>(m1*(m1+m1).transpose(), m2) );
  72. VERIFY( verify_eval_type<2>(m1*m1.transpose(), m2) );
  73. VERIFY( verify_eval_type<1>(m1+m1*m1, m1) );
  74. VERIFY( verify_eval_type<1>(m1.template triangularView<Lower>().solve(m1), m1) );
  75. VERIFY( verify_eval_type<1>(m1+m1.template triangularView<Lower>().solve(m1), m1) );
  76. }
  77. }
  78. void test_nesting_ops()
  79. {
  80. CALL_SUBTEST_1(run_nesting_ops_1(MatrixXf::Random(25,25)));
  81. CALL_SUBTEST_2(run_nesting_ops_1(MatrixXcd::Random(25,25)));
  82. CALL_SUBTEST_3(run_nesting_ops_1(Matrix4f::Random()));
  83. CALL_SUBTEST_4(run_nesting_ops_1(Matrix2d::Random()));
  84. Index s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE);
  85. CALL_SUBTEST_1( run_nesting_ops_2(MatrixXf(s,s)) );
  86. CALL_SUBTEST_2( run_nesting_ops_2(MatrixXcd(s,s)) );
  87. CALL_SUBTEST_3( run_nesting_ops_2(Matrix4f()) );
  88. CALL_SUBTEST_4( run_nesting_ops_2(Matrix2d()) );
  89. TEST_SET_BUT_UNUSED_VARIABLE(s)
  90. }