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.

195 lines
7.1 KiB

  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2010 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. // The computeRoots function included in this is based on materials
  10. // covered by the following copyright and license:
  11. //
  12. // Geometric Tools, LLC
  13. // Copyright (c) 1998-2010
  14. // Distributed under the Boost Software License, Version 1.0.
  15. //
  16. // Permission is hereby granted, free of charge, to any person or organization
  17. // obtaining a copy of the software and accompanying documentation covered by
  18. // this license (the "Software") to use, reproduce, display, distribute,
  19. // execute, and transmit the Software, and to prepare derivative works of the
  20. // Software, and to permit third-parties to whom the Software is furnished to
  21. // do so, all subject to the following:
  22. //
  23. // The copyright notices in the Software and this entire statement, including
  24. // the above license grant, this restriction and the following disclaimer,
  25. // must be included in all copies of the Software, in whole or in part, and
  26. // all derivative works of the Software, unless such copies or derivative
  27. // works are solely in the form of machine-executable object code generated by
  28. // a source language processor.
  29. //
  30. // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
  31. // IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
  32. // FITNESS FOR A PARTICULAR PURPOSE, TITLE AND NON-INFRINGEMENT. IN NO EVENT
  33. // SHALL THE COPYRIGHT HOLDERS OR ANYONE DISTRIBUTING THE SOFTWARE BE LIABLE
  34. // FOR ANY DAMAGES OR OTHER LIABILITY, WHETHER IN CONTRACT, TORT OR OTHERWISE,
  35. // ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
  36. // DEALINGS IN THE SOFTWARE.
  37. #include <iostream>
  38. #include <Eigen/Core>
  39. #include <Eigen/Eigenvalues>
  40. #include <Eigen/Geometry>
  41. #include <bench/BenchTimer.h>
  42. using namespace Eigen;
  43. using namespace std;
  44. template<typename Matrix, typename Roots>
  45. inline void computeRoots(const Matrix& m, Roots& roots)
  46. {
  47. typedef typename Matrix::Scalar Scalar;
  48. const Scalar s_inv3 = 1.0/3.0;
  49. const Scalar s_sqrt3 = std::sqrt(Scalar(3.0));
  50. // The characteristic equation is x^3 - c2*x^2 + c1*x - c0 = 0. The
  51. // eigenvalues are the roots to this equation, all guaranteed to be
  52. // real-valued, because the matrix is symmetric.
  53. Scalar c0 = m(0,0)*m(1,1)*m(2,2) + Scalar(2)*m(0,1)*m(0,2)*m(1,2) - m(0,0)*m(1,2)*m(1,2) - m(1,1)*m(0,2)*m(0,2) - m(2,2)*m(0,1)*m(0,1);
  54. Scalar c1 = m(0,0)*m(1,1) - m(0,1)*m(0,1) + m(0,0)*m(2,2) - m(0,2)*m(0,2) + m(1,1)*m(2,2) - m(1,2)*m(1,2);
  55. Scalar c2 = m(0,0) + m(1,1) + m(2,2);
  56. // Construct the parameters used in classifying the roots of the equation
  57. // and in solving the equation for the roots in closed form.
  58. Scalar c2_over_3 = c2*s_inv3;
  59. Scalar a_over_3 = (c1 - c2*c2_over_3)*s_inv3;
  60. if (a_over_3 > Scalar(0))
  61. a_over_3 = Scalar(0);
  62. Scalar half_b = Scalar(0.5)*(c0 + c2_over_3*(Scalar(2)*c2_over_3*c2_over_3 - c1));
  63. Scalar q = half_b*half_b + a_over_3*a_over_3*a_over_3;
  64. if (q > Scalar(0))
  65. q = Scalar(0);
  66. // Compute the eigenvalues by solving for the roots of the polynomial.
  67. Scalar rho = std::sqrt(-a_over_3);
  68. Scalar theta = std::atan2(std::sqrt(-q),half_b)*s_inv3;
  69. Scalar cos_theta = std::cos(theta);
  70. Scalar sin_theta = std::sin(theta);
  71. roots(2) = c2_over_3 + Scalar(2)*rho*cos_theta;
  72. roots(0) = c2_over_3 - rho*(cos_theta + s_sqrt3*sin_theta);
  73. roots(1) = c2_over_3 - rho*(cos_theta - s_sqrt3*sin_theta);
  74. }
  75. template<typename Matrix, typename Vector>
  76. void eigen33(const Matrix& mat, Matrix& evecs, Vector& evals)
  77. {
  78. typedef typename Matrix::Scalar Scalar;
  79. // Scale the matrix so its entries are in [-1,1]. The scaling is applied
  80. // only when at least one matrix entry has magnitude larger than 1.
  81. Scalar shift = mat.trace()/3;
  82. Matrix scaledMat = mat;
  83. scaledMat.diagonal().array() -= shift;
  84. Scalar scale = scaledMat.cwiseAbs()/*.template triangularView<Lower>()*/.maxCoeff();
  85. scale = std::max(scale,Scalar(1));
  86. scaledMat/=scale;
  87. // Compute the eigenvalues
  88. // scaledMat.setZero();
  89. computeRoots(scaledMat,evals);
  90. // compute the eigen vectors
  91. // **here we assume 3 differents eigenvalues**
  92. // "optimized version" which appears to be slower with gcc!
  93. // Vector base;
  94. // Scalar alpha, beta;
  95. // base << scaledMat(1,0) * scaledMat(2,1),
  96. // scaledMat(1,0) * scaledMat(2,0),
  97. // -scaledMat(1,0) * scaledMat(1,0);
  98. // for(int k=0; k<2; ++k)
  99. // {
  100. // alpha = scaledMat(0,0) - evals(k);
  101. // beta = scaledMat(1,1) - evals(k);
  102. // evecs.col(k) = (base + Vector(-beta*scaledMat(2,0), -alpha*scaledMat(2,1), alpha*beta)).normalized();
  103. // }
  104. // evecs.col(2) = evecs.col(0).cross(evecs.col(1)).normalized();
  105. // // naive version
  106. // Matrix tmp;
  107. // tmp = scaledMat;
  108. // tmp.diagonal().array() -= evals(0);
  109. // evecs.col(0) = tmp.row(0).cross(tmp.row(1)).normalized();
  110. //
  111. // tmp = scaledMat;
  112. // tmp.diagonal().array() -= evals(1);
  113. // evecs.col(1) = tmp.row(0).cross(tmp.row(1)).normalized();
  114. //
  115. // tmp = scaledMat;
  116. // tmp.diagonal().array() -= evals(2);
  117. // evecs.col(2) = tmp.row(0).cross(tmp.row(1)).normalized();
  118. // a more stable version:
  119. if((evals(2)-evals(0))<=Eigen::NumTraits<Scalar>::epsilon())
  120. {
  121. evecs.setIdentity();
  122. }
  123. else
  124. {
  125. Matrix tmp;
  126. tmp = scaledMat;
  127. tmp.diagonal ().array () -= evals (2);
  128. evecs.col (2) = tmp.row (0).cross (tmp.row (1)).normalized ();
  129. tmp = scaledMat;
  130. tmp.diagonal ().array () -= evals (1);
  131. evecs.col(1) = tmp.row (0).cross(tmp.row (1));
  132. Scalar n1 = evecs.col(1).norm();
  133. if(n1<=Eigen::NumTraits<Scalar>::epsilon())
  134. evecs.col(1) = evecs.col(2).unitOrthogonal();
  135. else
  136. evecs.col(1) /= n1;
  137. // make sure that evecs[1] is orthogonal to evecs[2]
  138. evecs.col(1) = evecs.col(2).cross(evecs.col(1).cross(evecs.col(2))).normalized();
  139. evecs.col(0) = evecs.col(2).cross(evecs.col(1));
  140. }
  141. // Rescale back to the original size.
  142. evals *= scale;
  143. evals.array()+=shift;
  144. }
  145. int main()
  146. {
  147. BenchTimer t;
  148. int tries = 10;
  149. int rep = 400000;
  150. typedef Matrix3d Mat;
  151. typedef Vector3d Vec;
  152. Mat A = Mat::Random(3,3);
  153. A = A.adjoint() * A;
  154. // Mat Q = A.householderQr().householderQ();
  155. // A = Q * Vec(2.2424567,2.2424566,7.454353).asDiagonal() * Q.transpose();
  156. SelfAdjointEigenSolver<Mat> eig(A);
  157. BENCH(t, tries, rep, eig.compute(A));
  158. std::cout << "Eigen iterative: " << t.best() << "s\n";
  159. BENCH(t, tries, rep, eig.computeDirect(A));
  160. std::cout << "Eigen direct : " << t.best() << "s\n";
  161. Mat evecs;
  162. Vec evals;
  163. BENCH(t, tries, rep, eigen33(A,evecs,evals));
  164. std::cout << "Direct: " << t.best() << "s\n\n";
  165. // std::cerr << "Eigenvalue/eigenvector diffs:\n";
  166. // std::cerr << (evals - eig.eigenvalues()).transpose() << "\n";
  167. // for(int k=0;k<3;++k)
  168. // if(evecs.col(k).dot(eig.eigenvectors().col(k))<0)
  169. // evecs.col(k) = -evecs.col(k);
  170. // std::cerr << evecs - eig.eigenvectors() << "\n\n";
  171. }