triangular.cpp 12 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292
  1. // This file is triangularView of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2008-2009 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. #ifdef EIGEN_TEST_PART_100
  10. # define EIGEN_NO_DEPRECATED_WARNING
  11. #endif
  12. #include "main.h"
  13. template<typename MatrixType> void triangular_deprecated(const MatrixType &m)
  14. {
  15. Index rows = m.rows();
  16. Index cols = m.cols();
  17. MatrixType m1, m2, m3, m4;
  18. m1.setRandom(rows,cols);
  19. m2.setRandom(rows,cols);
  20. m3 = m1; m4 = m2;
  21. // deprecated method:
  22. m1.template triangularView<Eigen::Upper>().swap(m2);
  23. // use this method instead:
  24. m3.template triangularView<Eigen::Upper>().swap(m4.template triangularView<Eigen::Upper>());
  25. VERIFY_IS_APPROX(m1,m3);
  26. VERIFY_IS_APPROX(m2,m4);
  27. // deprecated method:
  28. m1.template triangularView<Eigen::Lower>().swap(m4);
  29. // use this method instead:
  30. m3.template triangularView<Eigen::Lower>().swap(m2.template triangularView<Eigen::Lower>());
  31. VERIFY_IS_APPROX(m1,m3);
  32. VERIFY_IS_APPROX(m2,m4);
  33. }
  34. template<typename MatrixType> void triangular_square(const MatrixType& m)
  35. {
  36. typedef typename MatrixType::Scalar Scalar;
  37. typedef typename NumTraits<Scalar>::Real RealScalar;
  38. typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
  39. RealScalar largerEps = 10*test_precision<RealScalar>();
  40. Index rows = m.rows();
  41. Index cols = m.cols();
  42. MatrixType m1 = MatrixType::Random(rows, cols),
  43. m2 = MatrixType::Random(rows, cols),
  44. m3(rows, cols),
  45. m4(rows, cols),
  46. r1(rows, cols),
  47. r2(rows, cols);
  48. VectorType v2 = VectorType::Random(rows);
  49. MatrixType m1up = m1.template triangularView<Upper>();
  50. MatrixType m2up = m2.template triangularView<Upper>();
  51. if (rows*cols>1)
  52. {
  53. VERIFY(m1up.isUpperTriangular());
  54. VERIFY(m2up.transpose().isLowerTriangular());
  55. VERIFY(!m2.isLowerTriangular());
  56. }
  57. // VERIFY_IS_APPROX(m1up.transpose() * m2, m1.upper().transpose().lower() * m2);
  58. // test overloaded operator+=
  59. r1.setZero();
  60. r2.setZero();
  61. r1.template triangularView<Upper>() += m1;
  62. r2 += m1up;
  63. VERIFY_IS_APPROX(r1,r2);
  64. // test overloaded operator=
  65. m1.setZero();
  66. m1.template triangularView<Upper>() = m2.transpose() + m2;
  67. m3 = m2.transpose() + m2;
  68. VERIFY_IS_APPROX(m3.template triangularView<Lower>().transpose().toDenseMatrix(), m1);
  69. // test overloaded operator=
  70. m1.setZero();
  71. m1.template triangularView<Lower>() = m2.transpose() + m2;
  72. VERIFY_IS_APPROX(m3.template triangularView<Lower>().toDenseMatrix(), m1);
  73. VERIFY_IS_APPROX(m3.template triangularView<Lower>().conjugate().toDenseMatrix(),
  74. m3.conjugate().template triangularView<Lower>().toDenseMatrix());
  75. m1 = MatrixType::Random(rows, cols);
  76. for (int i=0; i<rows; ++i)
  77. while (numext::abs2(m1(i,i))<RealScalar(1e-1)) m1(i,i) = internal::random<Scalar>();
  78. Transpose<MatrixType> trm4(m4);
  79. // test back and forward substitution with a vector as the rhs
  80. m3 = m1.template triangularView<Upper>();
  81. VERIFY(v2.isApprox(m3.adjoint() * (m1.adjoint().template triangularView<Lower>().solve(v2)), largerEps));
  82. m3 = m1.template triangularView<Lower>();
  83. VERIFY(v2.isApprox(m3.transpose() * (m1.transpose().template triangularView<Upper>().solve(v2)), largerEps));
  84. m3 = m1.template triangularView<Upper>();
  85. VERIFY(v2.isApprox(m3 * (m1.template triangularView<Upper>().solve(v2)), largerEps));
  86. m3 = m1.template triangularView<Lower>();
  87. VERIFY(v2.isApprox(m3.conjugate() * (m1.conjugate().template triangularView<Lower>().solve(v2)), largerEps));
  88. // test back and forward substitution with a matrix as the rhs
  89. m3 = m1.template triangularView<Upper>();
  90. VERIFY(m2.isApprox(m3.adjoint() * (m1.adjoint().template triangularView<Lower>().solve(m2)), largerEps));
  91. m3 = m1.template triangularView<Lower>();
  92. VERIFY(m2.isApprox(m3.transpose() * (m1.transpose().template triangularView<Upper>().solve(m2)), largerEps));
  93. m3 = m1.template triangularView<Upper>();
  94. VERIFY(m2.isApprox(m3 * (m1.template triangularView<Upper>().solve(m2)), largerEps));
  95. m3 = m1.template triangularView<Lower>();
  96. VERIFY(m2.isApprox(m3.conjugate() * (m1.conjugate().template triangularView<Lower>().solve(m2)), largerEps));
  97. // check M * inv(L) using in place API
  98. m4 = m3;
  99. m1.transpose().template triangularView<Eigen::Upper>().solveInPlace(trm4);
  100. VERIFY_IS_APPROX(m4 * m1.template triangularView<Eigen::Lower>(), m3);
  101. // check M * inv(U) using in place API
  102. m3 = m1.template triangularView<Upper>();
  103. m4 = m3;
  104. m3.transpose().template triangularView<Eigen::Lower>().solveInPlace(trm4);
  105. VERIFY_IS_APPROX(m4 * m1.template triangularView<Eigen::Upper>(), m3);
  106. // check solve with unit diagonal
  107. m3 = m1.template triangularView<UnitUpper>();
  108. VERIFY(m2.isApprox(m3 * (m1.template triangularView<UnitUpper>().solve(m2)), largerEps));
  109. // VERIFY(( m1.template triangularView<Upper>()
  110. // * m2.template triangularView<Upper>()).isUpperTriangular());
  111. // test swap
  112. m1.setOnes();
  113. m2.setZero();
  114. m2.template triangularView<Upper>().swap(m1.template triangularView<Eigen::Upper>());
  115. m3.setZero();
  116. m3.template triangularView<Upper>().setOnes();
  117. VERIFY_IS_APPROX(m2,m3);
  118. VERIFY_RAISES_STATIC_ASSERT(m1.template triangularView<Eigen::Lower>().swap(m2.template triangularView<Eigen::Upper>()));
  119. m1.setRandom();
  120. m3 = m1.template triangularView<Upper>();
  121. Matrix<Scalar, MatrixType::ColsAtCompileTime, Dynamic> m5(cols, internal::random<int>(1,20)); m5.setRandom();
  122. Matrix<Scalar, Dynamic, MatrixType::RowsAtCompileTime> m6(internal::random<int>(1,20), rows); m6.setRandom();
  123. VERIFY_IS_APPROX(m1.template triangularView<Upper>() * m5, m3*m5);
  124. VERIFY_IS_APPROX(m6*m1.template triangularView<Upper>(), m6*m3);
  125. m1up = m1.template triangularView<Upper>();
  126. VERIFY_IS_APPROX(m1.template selfadjointView<Upper>().template triangularView<Upper>().toDenseMatrix(), m1up);
  127. VERIFY_IS_APPROX(m1up.template selfadjointView<Upper>().template triangularView<Upper>().toDenseMatrix(), m1up);
  128. VERIFY_IS_APPROX(m1.template selfadjointView<Upper>().template triangularView<Lower>().toDenseMatrix(), m1up.adjoint());
  129. VERIFY_IS_APPROX(m1up.template selfadjointView<Upper>().template triangularView<Lower>().toDenseMatrix(), m1up.adjoint());
  130. VERIFY_IS_APPROX(m1.template selfadjointView<Upper>().diagonal(), m1.diagonal());
  131. m3.setRandom();
  132. const MatrixType& m3c(m3);
  133. VERIFY( is_same_type(m3c.template triangularView<Lower>(),m3.template triangularView<Lower>().template conjugateIf<false>()) );
  134. VERIFY( is_same_type(m3c.template triangularView<Lower>().conjugate(),m3.template triangularView<Lower>().template conjugateIf<true>()) );
  135. VERIFY_IS_APPROX(m3.template triangularView<Lower>().template conjugateIf<true>().toDenseMatrix(),
  136. m3.conjugate().template triangularView<Lower>().toDenseMatrix());
  137. VERIFY_IS_APPROX(m3.template triangularView<Lower>().template conjugateIf<false>().toDenseMatrix(),
  138. m3.template triangularView<Lower>().toDenseMatrix());
  139. VERIFY( is_same_type(m3c.template selfadjointView<Lower>(),m3.template selfadjointView<Lower>().template conjugateIf<false>()) );
  140. VERIFY( is_same_type(m3c.template selfadjointView<Lower>().conjugate(),m3.template selfadjointView<Lower>().template conjugateIf<true>()) );
  141. VERIFY_IS_APPROX(m3.template selfadjointView<Lower>().template conjugateIf<true>().toDenseMatrix(),
  142. m3.conjugate().template selfadjointView<Lower>().toDenseMatrix());
  143. VERIFY_IS_APPROX(m3.template selfadjointView<Lower>().template conjugateIf<false>().toDenseMatrix(),
  144. m3.template selfadjointView<Lower>().toDenseMatrix());
  145. }
  146. template<typename MatrixType> void triangular_rect(const MatrixType& m)
  147. {
  148. typedef typename MatrixType::Scalar Scalar;
  149. typedef typename NumTraits<Scalar>::Real RealScalar;
  150. enum { Rows = MatrixType::RowsAtCompileTime, Cols = MatrixType::ColsAtCompileTime };
  151. Index rows = m.rows();
  152. Index cols = m.cols();
  153. MatrixType m1 = MatrixType::Random(rows, cols),
  154. m2 = MatrixType::Random(rows, cols),
  155. m3(rows, cols),
  156. m4(rows, cols),
  157. r1(rows, cols),
  158. r2(rows, cols);
  159. MatrixType m1up = m1.template triangularView<Upper>();
  160. MatrixType m2up = m2.template triangularView<Upper>();
  161. if (rows>1 && cols>1)
  162. {
  163. VERIFY(m1up.isUpperTriangular());
  164. VERIFY(m2up.transpose().isLowerTriangular());
  165. VERIFY(!m2.isLowerTriangular());
  166. }
  167. // test overloaded operator+=
  168. r1.setZero();
  169. r2.setZero();
  170. r1.template triangularView<Upper>() += m1;
  171. r2 += m1up;
  172. VERIFY_IS_APPROX(r1,r2);
  173. // test overloaded operator=
  174. m1.setZero();
  175. m1.template triangularView<Upper>() = 3 * m2;
  176. m3 = 3 * m2;
  177. VERIFY_IS_APPROX(m3.template triangularView<Upper>().toDenseMatrix(), m1);
  178. m1.setZero();
  179. m1.template triangularView<Lower>() = 3 * m2;
  180. VERIFY_IS_APPROX(m3.template triangularView<Lower>().toDenseMatrix(), m1);
  181. m1.setZero();
  182. m1.template triangularView<StrictlyUpper>() = 3 * m2;
  183. VERIFY_IS_APPROX(m3.template triangularView<StrictlyUpper>().toDenseMatrix(), m1);
  184. m1.setZero();
  185. m1.template triangularView<StrictlyLower>() = 3 * m2;
  186. VERIFY_IS_APPROX(m3.template triangularView<StrictlyLower>().toDenseMatrix(), m1);
  187. m1.setRandom();
  188. m2 = m1.template triangularView<Upper>();
  189. VERIFY(m2.isUpperTriangular());
  190. VERIFY(!m2.isLowerTriangular());
  191. m2 = m1.template triangularView<StrictlyUpper>();
  192. VERIFY(m2.isUpperTriangular());
  193. VERIFY(m2.diagonal().isMuchSmallerThan(RealScalar(1)));
  194. m2 = m1.template triangularView<UnitUpper>();
  195. VERIFY(m2.isUpperTriangular());
  196. m2.diagonal().array() -= Scalar(1);
  197. VERIFY(m2.diagonal().isMuchSmallerThan(RealScalar(1)));
  198. m2 = m1.template triangularView<Lower>();
  199. VERIFY(m2.isLowerTriangular());
  200. VERIFY(!m2.isUpperTriangular());
  201. m2 = m1.template triangularView<StrictlyLower>();
  202. VERIFY(m2.isLowerTriangular());
  203. VERIFY(m2.diagonal().isMuchSmallerThan(RealScalar(1)));
  204. m2 = m1.template triangularView<UnitLower>();
  205. VERIFY(m2.isLowerTriangular());
  206. m2.diagonal().array() -= Scalar(1);
  207. VERIFY(m2.diagonal().isMuchSmallerThan(RealScalar(1)));
  208. // test swap
  209. m1.setOnes();
  210. m2.setZero();
  211. m2.template triangularView<Upper>().swap(m1.template triangularView<Eigen::Upper>());
  212. m3.setZero();
  213. m3.template triangularView<Upper>().setOnes();
  214. VERIFY_IS_APPROX(m2,m3);
  215. }
  216. void bug_159()
  217. {
  218. Matrix3d m = Matrix3d::Random().triangularView<Lower>();
  219. EIGEN_UNUSED_VARIABLE(m)
  220. }
  221. EIGEN_DECLARE_TEST(triangular)
  222. {
  223. int maxsize = (std::min)(EIGEN_TEST_MAX_SIZE,20);
  224. for(int i = 0; i < g_repeat ; i++)
  225. {
  226. int r = internal::random<int>(2,maxsize); TEST_SET_BUT_UNUSED_VARIABLE(r)
  227. int c = internal::random<int>(2,maxsize); TEST_SET_BUT_UNUSED_VARIABLE(c)
  228. CALL_SUBTEST_1( triangular_square(Matrix<float, 1, 1>()) );
  229. CALL_SUBTEST_2( triangular_square(Matrix<float, 2, 2>()) );
  230. CALL_SUBTEST_3( triangular_square(Matrix3d()) );
  231. CALL_SUBTEST_4( triangular_square(Matrix<std::complex<float>,8, 8>()) );
  232. CALL_SUBTEST_5( triangular_square(MatrixXcd(r,r)) );
  233. CALL_SUBTEST_6( triangular_square(Matrix<float,Dynamic,Dynamic,RowMajor>(r, r)) );
  234. CALL_SUBTEST_7( triangular_rect(Matrix<float, 4, 5>()) );
  235. CALL_SUBTEST_8( triangular_rect(Matrix<double, 6, 2>()) );
  236. CALL_SUBTEST_9( triangular_rect(MatrixXcf(r, c)) );
  237. CALL_SUBTEST_5( triangular_rect(MatrixXcd(r, c)) );
  238. CALL_SUBTEST_6( triangular_rect(Matrix<float,Dynamic,Dynamic,RowMajor>(r, c)) );
  239. CALL_SUBTEST_100( triangular_deprecated(Matrix<float, 5, 7>()) );
  240. CALL_SUBTEST_100( triangular_deprecated(MatrixXd(r,c)) );
  241. }
  242. CALL_SUBTEST_1( bug_159() );
  243. }