conservative_resize.cpp 5.2 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167
  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2009 Hauke Heibel <hauke.heibel@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. #include "main.h"
  10. #include <Eigen/Core>
  11. #include "AnnoyingScalar.h"
  12. using namespace Eigen;
  13. template <typename Scalar, int Storage>
  14. void run_matrix_tests()
  15. {
  16. typedef Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic, Storage> MatrixType;
  17. MatrixType m, n;
  18. // boundary cases ...
  19. m = n = MatrixType::Random(50,50);
  20. m.conservativeResize(1,50);
  21. VERIFY_IS_APPROX(m, n.block(0,0,1,50));
  22. m = n = MatrixType::Random(50,50);
  23. m.conservativeResize(50,1);
  24. VERIFY_IS_APPROX(m, n.block(0,0,50,1));
  25. m = n = MatrixType::Random(50,50);
  26. m.conservativeResize(50,50);
  27. VERIFY_IS_APPROX(m, n.block(0,0,50,50));
  28. // random shrinking ...
  29. for (int i=0; i<25; ++i)
  30. {
  31. const Index rows = internal::random<Index>(1,50);
  32. const Index cols = internal::random<Index>(1,50);
  33. m = n = MatrixType::Random(50,50);
  34. m.conservativeResize(rows,cols);
  35. VERIFY_IS_APPROX(m, n.block(0,0,rows,cols));
  36. }
  37. // random growing with zeroing ...
  38. for (int i=0; i<25; ++i)
  39. {
  40. const Index rows = internal::random<Index>(50,75);
  41. const Index cols = internal::random<Index>(50,75);
  42. m = n = MatrixType::Random(50,50);
  43. m.conservativeResizeLike(MatrixType::Zero(rows,cols));
  44. VERIFY_IS_APPROX(m.block(0,0,n.rows(),n.cols()), n);
  45. VERIFY( rows<=50 || m.block(50,0,rows-50,cols).sum() == Scalar(0) );
  46. VERIFY( cols<=50 || m.block(0,50,rows,cols-50).sum() == Scalar(0) );
  47. }
  48. }
  49. template <typename Scalar>
  50. void run_vector_tests()
  51. {
  52. typedef Matrix<Scalar, 1, Eigen::Dynamic> VectorType;
  53. VectorType m, n;
  54. // boundary cases ...
  55. m = n = VectorType::Random(50);
  56. m.conservativeResize(1);
  57. VERIFY_IS_APPROX(m, n.segment(0,1));
  58. m = n = VectorType::Random(50);
  59. m.conservativeResize(50);
  60. VERIFY_IS_APPROX(m, n.segment(0,50));
  61. m = n = VectorType::Random(50);
  62. m.conservativeResize(m.rows(),1);
  63. VERIFY_IS_APPROX(m, n.segment(0,1));
  64. m = n = VectorType::Random(50);
  65. m.conservativeResize(m.rows(),50);
  66. VERIFY_IS_APPROX(m, n.segment(0,50));
  67. // random shrinking ...
  68. for (int i=0; i<50; ++i)
  69. {
  70. const int size = internal::random<int>(1,50);
  71. m = n = VectorType::Random(50);
  72. m.conservativeResize(size);
  73. VERIFY_IS_APPROX(m, n.segment(0,size));
  74. m = n = VectorType::Random(50);
  75. m.conservativeResize(m.rows(), size);
  76. VERIFY_IS_APPROX(m, n.segment(0,size));
  77. }
  78. // random growing with zeroing ...
  79. for (int i=0; i<50; ++i)
  80. {
  81. const int size = internal::random<int>(50,100);
  82. m = n = VectorType::Random(50);
  83. m.conservativeResizeLike(VectorType::Zero(size));
  84. VERIFY_IS_APPROX(m.segment(0,50), n);
  85. VERIFY( size<=50 || m.segment(50,size-50).sum() == Scalar(0) );
  86. m = n = VectorType::Random(50);
  87. m.conservativeResizeLike(Matrix<Scalar,Dynamic,Dynamic>::Zero(1,size));
  88. VERIFY_IS_APPROX(m.segment(0,50), n);
  89. VERIFY( size<=50 || m.segment(50,size-50).sum() == Scalar(0) );
  90. }
  91. }
  92. // Basic memory leak check with a non-copyable scalar type
  93. template<int> void noncopyable()
  94. {
  95. typedef Eigen::Matrix<AnnoyingScalar,Dynamic,1> VectorType;
  96. typedef Eigen::Matrix<AnnoyingScalar,Dynamic,Dynamic> MatrixType;
  97. {
  98. #ifndef EIGEN_TEST_ANNOYING_SCALAR_DONT_THROW
  99. AnnoyingScalar::dont_throw = true;
  100. #endif
  101. int n = 50;
  102. VectorType v0(n), v1(n);
  103. MatrixType m0(n,n), m1(n,n), m2(n,n);
  104. v0.setOnes(); v1.setOnes();
  105. m0.setOnes(); m1.setOnes(); m2.setOnes();
  106. VERIFY(m0==m1);
  107. m0.conservativeResize(2*n,2*n);
  108. VERIFY(m0.topLeftCorner(n,n) == m1);
  109. VERIFY(v0.head(n) == v1);
  110. v0.conservativeResize(2*n);
  111. VERIFY(v0.head(n) == v1);
  112. }
  113. VERIFY(AnnoyingScalar::instances==0 && "global memory leak detected in noncopyable");
  114. }
  115. EIGEN_DECLARE_TEST(conservative_resize)
  116. {
  117. for(int i=0; i<g_repeat; ++i)
  118. {
  119. CALL_SUBTEST_1((run_matrix_tests<int, Eigen::RowMajor>()));
  120. CALL_SUBTEST_1((run_matrix_tests<int, Eigen::ColMajor>()));
  121. CALL_SUBTEST_2((run_matrix_tests<float, Eigen::RowMajor>()));
  122. CALL_SUBTEST_2((run_matrix_tests<float, Eigen::ColMajor>()));
  123. CALL_SUBTEST_3((run_matrix_tests<double, Eigen::RowMajor>()));
  124. CALL_SUBTEST_3((run_matrix_tests<double, Eigen::ColMajor>()));
  125. CALL_SUBTEST_4((run_matrix_tests<std::complex<float>, Eigen::RowMajor>()));
  126. CALL_SUBTEST_4((run_matrix_tests<std::complex<float>, Eigen::ColMajor>()));
  127. CALL_SUBTEST_5((run_matrix_tests<std::complex<double>, Eigen::RowMajor>()));
  128. CALL_SUBTEST_5((run_matrix_tests<std::complex<double>, Eigen::ColMajor>()));
  129. CALL_SUBTEST_1((run_matrix_tests<int, Eigen::RowMajor | Eigen::DontAlign>()));
  130. CALL_SUBTEST_1((run_vector_tests<int>()));
  131. CALL_SUBTEST_2((run_vector_tests<float>()));
  132. CALL_SUBTEST_3((run_vector_tests<double>()));
  133. CALL_SUBTEST_4((run_vector_tests<std::complex<float> >()));
  134. CALL_SUBTEST_5((run_vector_tests<std::complex<double> >()));
  135. #ifndef EIGEN_TEST_ANNOYING_SCALAR_DONT_THROW
  136. AnnoyingScalar::dont_throw = true;
  137. #endif
  138. CALL_SUBTEST_6(( run_vector_tests<AnnoyingScalar>() ));
  139. CALL_SUBTEST_6(( noncopyable<0>() ));
  140. }
  141. }