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.

789 lines
30 KiB

  1. /* -*- c++ -*- (enables emacs c++ mode) */
  2. /*===========================================================================
  3. Copyright (C) 2003-2015 Yves Renard
  4. This file is a part of GETFEM++
  5. Getfem++ is free software; you can redistribute it and/or modify it
  6. under the terms of the GNU Lesser General Public License as published
  7. by the Free Software Foundation; either version 3 of the License, or
  8. (at your option) any later version along with the GCC Runtime Library
  9. Exception either version 3.1 or (at your option) any later version.
  10. This program is distributed in the hope that it will be useful, but
  11. WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
  12. or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
  13. License and GCC Runtime Library Exception for more details.
  14. You should have received a copy of the GNU Lesser General Public License
  15. along with this program; if not, write to the Free Software Foundation,
  16. Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301, USA.
  17. As a special exception, you may use this file as it is a part of a free
  18. software library without restriction. Specifically, if other files
  19. instantiate templates or use macros or inline functions from this file,
  20. or you compile this file and link it with other files to produce an
  21. executable, this file does not by itself cause the resulting executable
  22. to be covered by the GNU Lesser General Public License. This exception
  23. does not however invalidate any other reasons why the executable file
  24. might be covered by the GNU Lesser General Public License.
  25. ===========================================================================*/
  26. /**@file gmm_dense_qr.h
  27. @author Caroline Lecalvez, Caroline.Lecalvez@gmm.insa-tlse.fr, Yves Renard <Yves.Renard@insa-lyon.fr>
  28. @date September 12, 2003.
  29. @brief Dense QR factorization.
  30. */
  31. #ifndef GMM_DENSE_QR_H
  32. #define GMM_DENSE_QR_H
  33. #include "gmm_dense_Householder.h"
  34. namespace gmm {
  35. /**
  36. QR factorization using Householder method (complex and real version).
  37. */
  38. template <typename MAT1>
  39. void qr_factor(const MAT1 &A_) {
  40. MAT1 &A = const_cast<MAT1 &>(A_);
  41. typedef typename linalg_traits<MAT1>::value_type value_type;
  42. size_type m = mat_nrows(A), n = mat_ncols(A);
  43. GMM_ASSERT2(m >= n, "dimensions mismatch");
  44. std::vector<value_type> W(m), V(m);
  45. for (size_type j = 0; j < n; ++j) {
  46. sub_interval SUBI(j, m-j), SUBJ(j, n-j);
  47. V.resize(m-j); W.resize(n-j);
  48. for (size_type i = j; i < m; ++i) V[i-j] = A(i, j);
  49. house_vector(V);
  50. row_house_update(sub_matrix(A, SUBI, SUBJ), V, W);
  51. for (size_type i = j+1; i < m; ++i) A(i, j) = V[i-j];
  52. }
  53. }
  54. // QR comes from QR_factor(QR) where the upper triangular part stands for R
  55. // and the lower part contains the Householder reflectors.
  56. // A <- AQ
  57. template <typename MAT1, typename MAT2>
  58. void apply_house_right(const MAT1 &QR, const MAT2 &A_) {
  59. MAT2 &A = const_cast<MAT2 &>(A_);
  60. typedef typename linalg_traits<MAT1>::value_type T;
  61. size_type m = mat_nrows(QR), n = mat_ncols(QR);
  62. GMM_ASSERT2(m == mat_ncols(A), "dimensions mismatch");
  63. if (m == 0) return;
  64. std::vector<T> V(m), W(mat_nrows(A));
  65. V[0] = T(1);
  66. for (size_type j = 0; j < n; ++j) {
  67. V.resize(m-j);
  68. for (size_type i = j+1; i < m; ++i) V[i-j] = QR(i, j);
  69. col_house_update(sub_matrix(A, sub_interval(0, mat_nrows(A)),
  70. sub_interval(j, m-j)), V, W);
  71. }
  72. }
  73. // QR comes from QR_factor(QR) where the upper triangular part stands for R
  74. // and the lower part contains the Householder reflectors.
  75. // A <- Q*A
  76. template <typename MAT1, typename MAT2>
  77. void apply_house_left(const MAT1 &QR, const MAT2 &A_) {
  78. MAT2 &A = const_cast<MAT2 &>(A_);
  79. typedef typename linalg_traits<MAT1>::value_type T;
  80. size_type m = mat_nrows(QR), n = mat_ncols(QR);
  81. GMM_ASSERT2(m == mat_nrows(A), "dimensions mismatch");
  82. if (m == 0) return;
  83. std::vector<T> V(m), W(mat_ncols(A));
  84. V[0] = T(1);
  85. for (size_type j = 0; j < n; ++j) {
  86. V.resize(m-j);
  87. for (size_type i = j+1; i < m; ++i) V[i-j] = QR(i, j);
  88. row_house_update(sub_matrix(A, sub_interval(j, m-j),
  89. sub_interval(0, mat_ncols(A))), V, W);
  90. }
  91. }
  92. /** Compute the QR factorization, where Q is assembled. */
  93. template <typename MAT1, typename MAT2, typename MAT3>
  94. void qr_factor(const MAT1 &A, const MAT2 &QQ, const MAT3 &RR) {
  95. MAT2 &Q = const_cast<MAT2 &>(QQ); MAT3 &R = const_cast<MAT3 &>(RR);
  96. typedef typename linalg_traits<MAT1>::value_type value_type;
  97. size_type m = mat_nrows(A), n = mat_ncols(A);
  98. GMM_ASSERT2(m >= n, "dimensions mismatch");
  99. gmm::copy(A, Q);
  100. std::vector<value_type> W(m);
  101. dense_matrix<value_type> VV(m, n);
  102. for (size_type j = 0; j < n; ++j) {
  103. sub_interval SUBI(j, m-j), SUBJ(j, n-j);
  104. for (size_type i = j; i < m; ++i) VV(i,j) = Q(i, j);
  105. house_vector(sub_vector(mat_col(VV,j), SUBI));
  106. row_house_update(sub_matrix(Q, SUBI, SUBJ),
  107. sub_vector(mat_col(VV,j), SUBI), sub_vector(W, SUBJ));
  108. }
  109. gmm::copy(sub_matrix(Q, sub_interval(0, n), sub_interval(0, n)), R);
  110. gmm::copy(identity_matrix(), Q);
  111. for (size_type j = n-1; j != size_type(-1); --j) {
  112. sub_interval SUBI(j, m-j), SUBJ(j, n-j);
  113. row_house_update(sub_matrix(Q, SUBI, SUBJ),
  114. sub_vector(mat_col(VV,j), SUBI), sub_vector(W, SUBJ));
  115. }
  116. }
  117. ///@cond DOXY_SHOW_ALL_FUNCTIONS
  118. template <typename TA, typename TV, typename Ttol,
  119. typename MAT, typename VECT>
  120. void extract_eig(const MAT &A, VECT &V, Ttol tol, TA, TV) {
  121. size_type n = mat_nrows(A);
  122. if (n == 0) return;
  123. tol *= Ttol(2);
  124. Ttol tol_i = tol * gmm::abs(A(0,0)), tol_cplx = tol_i;
  125. for (size_type i = 0; i < n; ++i) {
  126. if (i < n-1) {
  127. tol_i = (gmm::abs(A(i,i))+gmm::abs(A(i+1,i+1)))*tol;
  128. tol_cplx = std::max(tol_cplx, tol_i);
  129. }
  130. if ((i < n-1) && gmm::abs(A(i+1,i)) >= tol_i) {
  131. TA tr = A(i,i) + A(i+1, i+1);
  132. TA det = A(i,i)*A(i+1, i+1) - A(i,i+1)*A(i+1, i);
  133. TA delta = tr*tr - TA(4) * det;
  134. if (delta < -tol_cplx) {
  135. GMM_WARNING1("A complex eigenvalue has been detected : "
  136. << std::complex<TA>(tr/TA(2), gmm::sqrt(-delta)/TA(2)));
  137. V[i] = V[i+1] = tr / TA(2);
  138. }
  139. else {
  140. delta = std::max(TA(0), delta);
  141. V[i ] = TA(tr + gmm::sqrt(delta))/ TA(2);
  142. V[i+1] = TA(tr - gmm::sqrt(delta))/ TA(2);
  143. }
  144. ++i;
  145. }
  146. else
  147. V[i] = TV(A(i,i));
  148. }
  149. }
  150. template <typename TA, typename TV, typename Ttol,
  151. typename MAT, typename VECT>
  152. void extract_eig(const MAT &A, VECT &V, Ttol tol, TA, std::complex<TV>) {
  153. size_type n = mat_nrows(A);
  154. tol *= Ttol(2);
  155. for (size_type i = 0; i < n; ++i)
  156. if ((i == n-1) ||
  157. gmm::abs(A(i+1,i)) < (gmm::abs(A(i,i))+gmm::abs(A(i+1,i+1)))*tol)
  158. V[i] = std::complex<TV>(A(i,i));
  159. else {
  160. TA tr = A(i,i) + A(i+1, i+1);
  161. TA det = A(i,i)*A(i+1, i+1) - A(i,i+1)*A(i+1, i);
  162. TA delta = tr*tr - TA(4) * det;
  163. if (delta < TA(0)) {
  164. V[i] = std::complex<TV>(tr / TA(2), gmm::sqrt(-delta) / TA(2));
  165. V[i+1] = std::complex<TV>(tr / TA(2), -gmm::sqrt(-delta)/ TA(2));
  166. }
  167. else {
  168. V[i ] = TA(tr + gmm::sqrt(delta)) / TA(2);
  169. V[i+1] = TA(tr - gmm::sqrt(delta)) / TA(2);
  170. }
  171. ++i;
  172. }
  173. }
  174. template <typename TA, typename TV, typename Ttol,
  175. typename MAT, typename VECT>
  176. void extract_eig(const MAT &A, VECT &V, Ttol tol, std::complex<TA>, TV) {
  177. typedef std::complex<TA> T;
  178. size_type n = mat_nrows(A);
  179. if (n == 0) return;
  180. tol *= Ttol(2);
  181. Ttol tol_i = tol * gmm::abs(A(0,0)), tol_cplx = tol_i;
  182. for (size_type i = 0; i < n; ++i) {
  183. if (i < n-1) {
  184. tol_i = (gmm::abs(A(i,i))+gmm::abs(A(i+1,i+1)))*tol;
  185. tol_cplx = std::max(tol_cplx, tol_i);
  186. }
  187. if ((i == n-1) || gmm::abs(A(i+1,i)) < tol_i) {
  188. if (gmm::abs(std::imag(A(i,i))) > tol_cplx)
  189. GMM_WARNING1("A complex eigenvalue has been detected : "
  190. << T(A(i,i)) << " : " << gmm::abs(std::imag(A(i,i)))
  191. / gmm::abs(std::real(A(i,i))) << " : " << tol_cplx);
  192. V[i] = std::real(A(i,i));
  193. }
  194. else {
  195. T tr = A(i,i) + A(i+1, i+1);
  196. T det = A(i,i)*A(i+1, i+1) - A(i,i+1)*A(i+1, i);
  197. T delta = tr*tr - TA(4) * det;
  198. T a1 = (tr + gmm::sqrt(delta)) / TA(2);
  199. T a2 = (tr - gmm::sqrt(delta)) / TA(2);
  200. if (gmm::abs(std::imag(a1)) > tol_cplx)
  201. GMM_WARNING1("A complex eigenvalue has been detected : " << a1);
  202. if (gmm::abs(std::imag(a2)) > tol_cplx)
  203. GMM_WARNING1("A complex eigenvalue has been detected : " << a2);
  204. V[i] = std::real(a1); V[i+1] = std::real(a2);
  205. ++i;
  206. }
  207. }
  208. }
  209. template <typename TA, typename TV, typename Ttol,
  210. typename MAT, typename VECT>
  211. void extract_eig(const MAT &A, VECT &V, Ttol tol,
  212. std::complex<TA>, std::complex<TV>) {
  213. size_type n = mat_nrows(A);
  214. tol *= Ttol(2);
  215. for (size_type i = 0; i < n; ++i)
  216. if ((i == n-1) ||
  217. gmm::abs(A(i+1,i)) < (gmm::abs(A(i,i))+gmm::abs(A(i+1,i+1)))*tol)
  218. V[i] = std::complex<TV>(A(i,i));
  219. else {
  220. std::complex<TA> tr = A(i,i) + A(i+1, i+1);
  221. std::complex<TA> det = A(i,i)*A(i+1, i+1) - A(i,i+1)*A(i+1, i);
  222. std::complex<TA> delta = tr*tr - TA(4) * det;
  223. V[i] = (tr + gmm::sqrt(delta)) / TA(2);
  224. V[i+1] = (tr - gmm::sqrt(delta)) / TA(2);
  225. ++i;
  226. }
  227. }
  228. ///@endcond
  229. /**
  230. Compute eigenvalue vector.
  231. */
  232. template <typename MAT, typename Ttol, typename VECT> inline
  233. void extract_eig(const MAT &A, const VECT &V, Ttol tol) {
  234. extract_eig(A, const_cast<VECT&>(V), tol,
  235. typename linalg_traits<MAT>::value_type(),
  236. typename linalg_traits<VECT>::value_type());
  237. }
  238. /* ********************************************************************* */
  239. /* Stop criterion for QR algorithms */
  240. /* ********************************************************************* */
  241. template <typename MAT, typename Ttol>
  242. void qr_stop_criterion(MAT &A, size_type &p, size_type &q, Ttol tol) {
  243. typedef typename linalg_traits<MAT>::value_type T;
  244. typedef typename number_traits<T>::magnitude_type R;
  245. R rmin = default_min(R()) * R(2);
  246. size_type n = mat_nrows(A);
  247. if (n <= 2) { q = n; p = 0; }
  248. else {
  249. for (size_type i = 1; i < n-q; ++i)
  250. if (gmm::abs(A(i,i-1)) < (gmm::abs(A(i,i))+ gmm::abs(A(i-1,i-1)))*tol
  251. || gmm::abs(A(i,i-1)) < rmin)
  252. A(i,i-1) = T(0);
  253. while ((q < n-1 && A(n-1-q, n-2-q) == T(0)) ||
  254. (q < n-2 && A(n-2-q, n-3-q) == T(0))) ++q;
  255. if (q >= n-2) q = n;
  256. p = n-q; if (p) --p; if (p) --p;
  257. while (p > 0 && A(p,p-1) != T(0)) --p;
  258. }
  259. }
  260. template <typename MAT, typename Ttol> inline
  261. void symmetric_qr_stop_criterion(const MAT &AA, size_type &p, size_type &q,
  262. Ttol tol) {
  263. typedef typename linalg_traits<MAT>::value_type T;
  264. typedef typename number_traits<T>::magnitude_type R;
  265. R rmin = default_min(R()) * R(2);
  266. MAT& A = const_cast<MAT&>(AA);
  267. size_type n = mat_nrows(A);
  268. if (n <= 1) { q = n; p = 0; }
  269. else {
  270. for (size_type i = 1; i < n-q; ++i)
  271. if (gmm::abs(A(i,i-1)) < (gmm::abs(A(i,i))+ gmm::abs(A(i-1,i-1)))*tol
  272. || gmm::abs(A(i,i-1)) < rmin)
  273. A(i,i-1) = T(0);
  274. while (q < n-1 && A(n-1-q, n-2-q) == T(0)) ++q;
  275. if (q >= n-1) q = n;
  276. p = n-q; if (p) --p; if (p) --p;
  277. while (p > 0 && A(p,p-1) != T(0)) --p;
  278. }
  279. }
  280. template <typename VECT1, typename VECT2, typename Ttol> inline
  281. void symmetric_qr_stop_criterion(const VECT1 &diag, const VECT2 &sdiag_,
  282. size_type &p, size_type &q, Ttol tol) {
  283. typedef typename linalg_traits<VECT2>::value_type T;
  284. typedef typename number_traits<T>::magnitude_type R;
  285. R rmin = default_min(R()) * R(2);
  286. VECT2 &sdiag = const_cast<VECT2 &>(sdiag_);
  287. size_type n = vect_size(diag);
  288. if (n <= 1) { q = n; p = 0; return; }
  289. for (size_type i = 1; i < n-q; ++i)
  290. if (gmm::abs(sdiag[i-1]) < (gmm::abs(diag[i])+ gmm::abs(diag[i-1]))*tol
  291. || gmm::abs(sdiag[i-1]) < rmin)
  292. sdiag[i-1] = T(0);
  293. while (q < n-1 && sdiag[n-2-q] == T(0)) ++q;
  294. if (q >= n-1) q = n;
  295. p = n-q; if (p) --p; if (p) --p;
  296. while (p > 0 && sdiag[p-1] != T(0)) --p;
  297. }
  298. /* ********************************************************************* */
  299. /* 2x2 blocks reduction for Schur vectors */
  300. /* ********************************************************************* */
  301. template <typename MATH, typename MATQ, typename Ttol>
  302. void block2x2_reduction(MATH &H, MATQ &Q, Ttol tol) {
  303. typedef typename linalg_traits<MATH>::value_type T;
  304. typedef typename number_traits<T>::magnitude_type R;
  305. size_type n = mat_nrows(H), nq = mat_nrows(Q);
  306. if (n < 2) return;
  307. sub_interval SUBQ(0, nq), SUBL(0, 2);
  308. std::vector<T> v(2), w(std::max(n, nq)); v[0] = T(1);
  309. tol *= Ttol(2);
  310. Ttol tol_i = tol * gmm::abs(H(0,0)), tol_cplx = tol_i;
  311. for (size_type i = 0; i < n-1; ++i) {
  312. tol_i = (gmm::abs(H(i,i))+gmm::abs(H(i+1,i+1)))*tol;
  313. tol_cplx = std::max(tol_cplx, tol_i);
  314. if (gmm::abs(H(i+1,i)) > tol_i) { // 2x2 block detected
  315. T tr = (H(i+1, i+1) - H(i,i)) / T(2);
  316. T delta = tr*tr + H(i,i+1)*H(i+1, i);
  317. if (is_complex(T()) || gmm::real(delta) >= R(0)) {
  318. sub_interval SUBI(i, 2);
  319. T theta = (tr - gmm::sqrt(delta)) / H(i+1,i);
  320. R a = gmm::abs(theta);
  321. v[1] = (a == R(0)) ? T(-1)
  322. : gmm::conj(theta) * (R(1) - gmm::sqrt(a*a + R(1)) / a);
  323. row_house_update(sub_matrix(H, SUBI), v, sub_vector(w, SUBL));
  324. col_house_update(sub_matrix(H, SUBI), v, sub_vector(w, SUBL));
  325. col_house_update(sub_matrix(Q, SUBQ, SUBI), v, sub_vector(w, SUBQ));
  326. }
  327. ++i;
  328. }
  329. }
  330. }
  331. /* ********************************************************************* */
  332. /* Basic qr algorithm. */
  333. /* ********************************************************************* */
  334. #define tol_type_for_qr typename number_traits<typename \
  335. linalg_traits<MAT1>::value_type>::magnitude_type
  336. #define default_tol_for_qr \
  337. (gmm::default_tol(tol_type_for_qr()) * tol_type_for_qr(3))
  338. // QR method for real or complex square matrices based on QR factorisation.
  339. // eigval has to be a complex vector if A has complex eigeinvalues.
  340. // Very slow method. Use implicit_qr_method instead.
  341. template <typename MAT1, typename VECT, typename MAT2>
  342. void rudimentary_qr_algorithm(const MAT1 &A, const VECT &eigval_,
  343. const MAT2 &eigvect_,
  344. tol_type_for_qr tol = default_tol_for_qr,
  345. bool compvect = true) {
  346. VECT &eigval = const_cast<VECT &>(eigval_);
  347. MAT2 &eigvect = const_cast<MAT2 &>(eigvect_);
  348. typedef typename linalg_traits<MAT1>::value_type value_type;
  349. size_type n = mat_nrows(A), p, q = 0, ite = 0;
  350. dense_matrix<value_type> Q(n, n), R(n,n), A1(n,n);
  351. gmm::copy(A, A1);
  352. Hessenberg_reduction(A1, eigvect, compvect);
  353. qr_stop_criterion(A1, p, q, tol);
  354. while (q < n) {
  355. qr_factor(A1, Q, R);
  356. gmm::mult(R, Q, A1);
  357. if (compvect) { gmm::mult(eigvect, Q, R); gmm::copy(R, eigvect); }
  358. qr_stop_criterion(A1, p, q, tol);
  359. ++ite;
  360. GMM_ASSERT1(ite < n*1000, "QR algorithm failed");
  361. }
  362. if (compvect) block2x2_reduction(A1, Q, tol);
  363. extract_eig(A1, eigval, tol);
  364. }
  365. template <typename MAT1, typename VECT>
  366. void rudimentary_qr_algorithm(const MAT1 &a, VECT &eigval,
  367. tol_type_for_qr tol = default_tol_for_qr) {
  368. dense_matrix<typename linalg_traits<MAT1>::value_type> m(0,0);
  369. rudimentary_qr_algorithm(a, eigval, m, tol, false);
  370. }
  371. /* ********************************************************************* */
  372. /* Francis QR step. */
  373. /* ********************************************************************* */
  374. template <typename MAT1, typename MAT2>
  375. void Francis_qr_step(const MAT1& HH, const MAT2 &QQ, bool compute_Q) {
  376. MAT1& H = const_cast<MAT1&>(HH); MAT2& Q = const_cast<MAT2&>(QQ);
  377. typedef typename linalg_traits<MAT1>::value_type value_type;
  378. size_type n = mat_nrows(H), nq = mat_nrows(Q);
  379. std::vector<value_type> v(3), w(std::max(n, nq));
  380. value_type s = H(n-2, n-2) + H(n-1, n-1);
  381. value_type t = H(n-2, n-2) * H(n-1, n-1) - H(n-2, n-1) * H(n-1, n-2);
  382. value_type x = H(0, 0) * H(0, 0) + H(0,1) * H(1, 0) - s * H(0,0) + t;
  383. value_type y = H(1, 0) * (H(0,0) + H(1,1) - s);
  384. value_type z = H(1, 0) * H(2, 1);
  385. sub_interval SUBQ(0, nq);
  386. for (size_type k = 0; k < n - 2; ++k) {
  387. v[0] = x; v[1] = y; v[2] = z;
  388. house_vector(v);
  389. size_type r = std::min(k+4, n), q = (k==0) ? 0 : k-1;
  390. sub_interval SUBI(k, 3), SUBJ(0, r), SUBK(q, n-q);
  391. row_house_update(sub_matrix(H, SUBI, SUBK), v, sub_vector(w, SUBK));
  392. col_house_update(sub_matrix(H, SUBJ, SUBI), v, sub_vector(w, SUBJ));
  393. if (compute_Q)
  394. col_house_update(sub_matrix(Q, SUBQ, SUBI), v, sub_vector(w, SUBQ));
  395. x = H(k+1, k); y = H(k+2, k);
  396. if (k < n-3) z = H(k+3, k);
  397. }
  398. sub_interval SUBI(n-2,2), SUBJ(0, n), SUBK(n-3,3), SUBL(0, 3);
  399. v.resize(2);
  400. v[0] = x; v[1] = y;
  401. house_vector(v);
  402. row_house_update(sub_matrix(H, SUBI, SUBK), v, sub_vector(w, SUBL));
  403. col_house_update(sub_matrix(H, SUBJ, SUBI), v, sub_vector(w, SUBJ));
  404. if (compute_Q)
  405. col_house_update(sub_matrix(Q, SUBQ, SUBI), v, sub_vector(w, SUBQ));
  406. }
  407. /* ********************************************************************* */
  408. /* Wilkinson Double shift QR step (from Lapack). */
  409. /* ********************************************************************* */
  410. template <typename MAT1, typename MAT2, typename Ttol>
  411. void Wilkinson_double_shift_qr_step(const MAT1& HH, const MAT2 &QQ,
  412. Ttol tol, bool exc, bool compute_Q) {
  413. MAT1& H = const_cast<MAT1&>(HH); MAT2& Q = const_cast<MAT2&>(QQ);
  414. typedef typename linalg_traits<MAT1>::value_type T;
  415. typedef typename number_traits<T>::magnitude_type R;
  416. size_type n = mat_nrows(H), nq = mat_nrows(Q), m;
  417. std::vector<T> v(3), w(std::max(n, nq));
  418. const R dat1(0.75), dat2(-0.4375);
  419. T h33, h44, h43h34, v1(0), v2(0), v3(0);
  420. if (exc) { /* Exceptional shift. */
  421. R s = gmm::abs(H(n-1, n-2)) + gmm::abs(H(n-2, n-3));
  422. h33 = h44 = dat1 * s;
  423. h43h34 = dat2*s*s;
  424. }
  425. else { /* Wilkinson double shift. */
  426. h44 = H(n-1,n-1); h33 = H(n-2, n-2);
  427. h43h34 = H(n-1, n-2) * H(n-2, n-1);
  428. }
  429. /* Look for two consecutive small subdiagonal elements. */
  430. /* Determine the effect of starting the double-shift QR iteration at */
  431. /* row m, and see if this would make H(m-1, m-2) negligible. */
  432. for (m = n-2; m != 0; --m) {
  433. T h11 = H(m-1, m-1), h22 = H(m, m);
  434. T h21 = H(m, m-1), h12 = H(m-1, m);
  435. T h44s = h44 - h11, h33s = h33 - h11;
  436. v1 = (h33s*h44s-h43h34) / h21 + h12;
  437. v2 = h22 - h11 - h33s - h44s;
  438. v3 = H(m+1, m);
  439. R s = gmm::abs(v1) + gmm::abs(v2) + gmm::abs(v3);
  440. v1 /= s; v2 /= s; v3 /= s;
  441. if (m == 1) break;
  442. T h00 = H(m-2, m-2);
  443. T h10 = H(m-1, m-2);
  444. R tst1 = gmm::abs(v1)*(gmm::abs(h00)+gmm::abs(h11)+gmm::abs(h22));
  445. if (gmm::abs(h10)*(gmm::abs(v2)+gmm::abs(v3)) <= tol * tst1) break;
  446. }
  447. /* Double shift QR step. */
  448. sub_interval SUBQ(0, nq);
  449. for (size_type k = (m == 0) ? 0 : m-1; k < n-2; ++k) {
  450. v[0] = v1; v[1] = v2; v[2] = v3;
  451. house_vector(v);
  452. size_type r = std::min(k+4, n), q = (k==0) ? 0 : k-1;
  453. sub_interval SUBI(k, 3), SUBJ(0, r), SUBK(q, n-q);
  454. row_house_update(sub_matrix(H, SUBI, SUBK), v, sub_vector(w, SUBK));
  455. col_house_update(sub_matrix(H, SUBJ, SUBI), v, sub_vector(w, SUBJ));
  456. if (k > m-1) { H(k+1, k-1) = T(0); if (k < n-3) H(k+2, k-1) = T(0); }
  457. if (compute_Q)
  458. col_house_update(sub_matrix(Q, SUBQ, SUBI), v, sub_vector(w, SUBQ));
  459. v1 = H(k+1, k); v2 = H(k+2, k);
  460. if (k < n-3) v3 = H(k+3, k);
  461. }
  462. sub_interval SUBI(n-2,2), SUBJ(0, n), SUBK(n-3,3), SUBL(0, 3);
  463. v.resize(2); v[0] = v1; v[1] = v2;
  464. house_vector(v);
  465. row_house_update(sub_matrix(H, SUBI, SUBK), v, sub_vector(w, SUBL));
  466. col_house_update(sub_matrix(H, SUBJ, SUBI), v, sub_vector(w, SUBJ));
  467. if (compute_Q)
  468. col_house_update(sub_matrix(Q, SUBQ, SUBI), v, sub_vector(w, SUBQ));
  469. }
  470. /* ********************************************************************* */
  471. /* Implicit QR algorithm. */
  472. /* ********************************************************************* */
  473. // QR method for real or complex square matrices based on an
  474. // implicit QR factorisation. eigval has to be a complex vector
  475. // if A has complex eigenvalues. Complexity about 10n^3, 25n^3 if
  476. // eigenvectors are computed
  477. template <typename MAT1, typename VECT, typename MAT2>
  478. void implicit_qr_algorithm(const MAT1 &A, const VECT &eigval_,
  479. const MAT2 &Q_,
  480. tol_type_for_qr tol = default_tol_for_qr,
  481. bool compvect = true) {
  482. VECT &eigval = const_cast<VECT &>(eigval_);
  483. MAT2 &Q = const_cast<MAT2 &>(Q_);
  484. typedef typename linalg_traits<MAT1>::value_type value_type;
  485. size_type n(mat_nrows(A)), q(0), q_old, p(0), ite(0), its(0);
  486. dense_matrix<value_type> H(n,n);
  487. sub_interval SUBK(0,0);
  488. gmm::copy(A, H);
  489. Hessenberg_reduction(H, Q, compvect);
  490. qr_stop_criterion(H, p, q, tol);
  491. while (q < n) {
  492. sub_interval SUBI(p, n-p-q), SUBJ(0, mat_ncols(Q));
  493. if (compvect) SUBK = SUBI;
  494. // Francis_qr_step(sub_matrix(H, SUBI),
  495. // sub_matrix(Q, SUBJ, SUBK), compvect);
  496. Wilkinson_double_shift_qr_step(sub_matrix(H, SUBI),
  497. sub_matrix(Q, SUBJ, SUBK),
  498. tol, (its == 10 || its == 20), compvect);
  499. q_old = q;
  500. qr_stop_criterion(H, p, q, tol*2);
  501. if (q != q_old) its = 0;
  502. ++its; ++ite;
  503. GMM_ASSERT1(ite < n*100, "QR algorithm failed");
  504. }
  505. if (compvect) block2x2_reduction(H, Q, tol);
  506. extract_eig(H, eigval, tol);
  507. }
  508. template <typename MAT1, typename VECT>
  509. void implicit_qr_algorithm(const MAT1 &a, VECT &eigval,
  510. tol_type_for_qr tol = default_tol_for_qr) {
  511. dense_matrix<typename linalg_traits<MAT1>::value_type> m(1,1);
  512. implicit_qr_algorithm(a, eigval, m, tol, false);
  513. }
  514. /* ********************************************************************* */
  515. /* Implicit symmetric QR step with Wilkinson Shift. */
  516. /* ********************************************************************* */
  517. template <typename MAT1, typename MAT2>
  518. void symmetric_Wilkinson_qr_step(const MAT1& MM, const MAT2 &ZZ,
  519. bool compute_z) {
  520. MAT1& M = const_cast<MAT1&>(MM); MAT2& Z = const_cast<MAT2&>(ZZ);
  521. typedef typename linalg_traits<MAT1>::value_type T;
  522. typedef typename number_traits<T>::magnitude_type R;
  523. size_type n = mat_nrows(M);
  524. for (size_type i = 0; i < n; ++i) {
  525. M(i, i) = T(gmm::real(M(i, i)));
  526. if (i > 0) {
  527. T a = (M(i, i-1) + gmm::conj(M(i-1, i)))/R(2);
  528. M(i, i-1) = a; M(i-1, i) = gmm::conj(a);
  529. }
  530. }
  531. R d = gmm::real(M(n-2, n-2) - M(n-1, n-1)) / R(2);
  532. R e = gmm::abs_sqr(M(n-1, n-2));
  533. R nu = d + gmm::sgn(d)*gmm::sqrt(d*d+e);
  534. if (nu == R(0)) { M(n-1, n-2) = T(0); return; }
  535. R mu = gmm::real(M(n-1, n-1)) - e / nu;
  536. T x = M(0,0) - T(mu), z = M(1, 0), c, s;
  537. for (size_type k = 1; k < n; ++k) {
  538. Givens_rotation(x, z, c, s);
  539. if (k > 1) Apply_Givens_rotation_left(M(k-1,k-2), M(k,k-2), c, s);
  540. Apply_Givens_rotation_left(M(k-1,k-1), M(k,k-1), c, s);
  541. Apply_Givens_rotation_left(M(k-1,k ), M(k,k ), c, s);
  542. if (k < n-1) Apply_Givens_rotation_left(M(k-1,k+1), M(k,k+1), c, s);
  543. if (k > 1) Apply_Givens_rotation_right(M(k-2,k-1), M(k-2,k), c, s);
  544. Apply_Givens_rotation_right(M(k-1,k-1), M(k-1,k), c, s);
  545. Apply_Givens_rotation_right(M(k ,k-1), M(k,k) , c, s);
  546. if (k < n-1) Apply_Givens_rotation_right(M(k+1,k-1), M(k+1,k), c, s);
  547. if (compute_z) col_rot(Z, c, s, k-1, k);
  548. if (k < n-1) { x = M(k, k-1); z = M(k+1, k-1); }
  549. }
  550. }
  551. template <typename VECT1, typename VECT2, typename MAT>
  552. void symmetric_Wilkinson_qr_step(const VECT1& diag_, const VECT2& sdiag_,
  553. const MAT &ZZ, bool compute_z) {
  554. VECT1& diag = const_cast<VECT1&>(diag_);
  555. VECT2& sdiag = const_cast<VECT2&>(sdiag_);
  556. MAT& Z = const_cast<MAT&>(ZZ);
  557. typedef typename linalg_traits<VECT2>::value_type T;
  558. typedef typename number_traits<T>::magnitude_type R;
  559. size_type n = vect_size(diag);
  560. R d = (diag[n-2] - diag[n-1]) / R(2);
  561. R e = gmm::abs_sqr(sdiag[n-2]);
  562. R nu = d + gmm::sgn(d)*gmm::sqrt(d*d+e);
  563. if (nu == R(0)) { sdiag[n-2] = T(0); return; }
  564. R mu = diag[n-1] - e / nu;
  565. T x = diag[0] - T(mu), z = sdiag[0], c, s;
  566. T a01(0), a02(0);
  567. T a10(0), a11(diag[0]), a12(gmm::conj(sdiag[0])), a13(0);
  568. T a20(0), a21(sdiag[0]), a22(diag[1]), a23(gmm::conj(sdiag[1]));
  569. T a31(0), a32(sdiag[1]);
  570. for (size_type k = 1; k < n; ++k) {
  571. Givens_rotation(x, z, c, s);
  572. if (k > 1) Apply_Givens_rotation_left(a10, a20, c, s);
  573. Apply_Givens_rotation_left(a11, a21, c, s);
  574. Apply_Givens_rotation_left(a12, a22, c, s);
  575. if (k < n-1) Apply_Givens_rotation_left(a13, a23, c, s);
  576. if (k > 1) Apply_Givens_rotation_right(a01, a02, c, s);
  577. Apply_Givens_rotation_right(a11, a12, c, s);
  578. Apply_Givens_rotation_right(a21, a22, c, s);
  579. if (k < n-1) Apply_Givens_rotation_right(a31, a32, c, s);
  580. if (compute_z) col_rot(Z, c, s, k-1, k);
  581. diag[k-1] = gmm::real(a11);
  582. diag[k] = gmm::real(a22);
  583. if (k > 1) sdiag[k-2] = (gmm::conj(a01) + a10) / R(2);
  584. sdiag[k-1] = (gmm::conj(a12) + a21) / R(2);
  585. x = sdiag[k-1]; z = (gmm::conj(a13) + a31) / R(2);
  586. a01 = a12; a02 = a13;
  587. a10 = a21; a11 = a22; a12 = a23; a13 = T(0);
  588. a20 = a31; a21 = a32; a31 = T(0);
  589. if (k < n-1) {
  590. sdiag[k] = (gmm::conj(a23) + a32) / R(2);
  591. a22 = T(diag[k+1]); a32 = sdiag[k+1]; a23 = gmm::conj(a32);
  592. }
  593. }
  594. }
  595. /* ********************************************************************* */
  596. /* Implicit QR algorithm for symmetric or hermitian matrices. */
  597. /* ********************************************************************* */
  598. // implicit QR method for real square symmetric matrices or complex
  599. // hermitian matrices.
  600. // eigval has to be a complex vector if A has complex eigeinvalues.
  601. // complexity about 4n^3/3, 9n^3 if eigenvectors are computed
  602. template <typename MAT1, typename VECT, typename MAT2>
  603. void symmetric_qr_algorithm_old(const MAT1 &A, const VECT &eigval_,
  604. const MAT2 &eigvect_,
  605. tol_type_for_qr tol = default_tol_for_qr,
  606. bool compvect = true) {
  607. VECT &eigval = const_cast<VECT &>(eigval_);
  608. MAT2 &eigvect = const_cast<MAT2 &>(eigvect_);
  609. typedef typename linalg_traits<MAT1>::value_type T;
  610. typedef typename number_traits<T>::magnitude_type R;
  611. if (compvect) gmm::copy(identity_matrix(), eigvect);
  612. size_type n = mat_nrows(A), q = 0, p, ite = 0;
  613. dense_matrix<T> Tri(n, n);
  614. gmm::copy(A, Tri);
  615. Householder_tridiagonalization(Tri, eigvect, compvect);
  616. symmetric_qr_stop_criterion(Tri, p, q, tol);
  617. while (q < n) {
  618. sub_interval SUBI(p, n-p-q), SUBJ(0, mat_ncols(eigvect)), SUBK(p, n-p-q);
  619. if (!compvect) SUBK = sub_interval(0,0);
  620. symmetric_Wilkinson_qr_step(sub_matrix(Tri, SUBI),
  621. sub_matrix(eigvect, SUBJ, SUBK), compvect);
  622. symmetric_qr_stop_criterion(Tri, p, q, tol*R(2));
  623. ++ite;
  624. GMM_ASSERT1(ite < n*100, "QR algorithm failed. Probably, your matrix"
  625. " is not real symmetric or complex hermitian");
  626. }
  627. extract_eig(Tri, eigval, tol);
  628. }
  629. template <typename MAT1, typename VECT, typename MAT2>
  630. void symmetric_qr_algorithm(const MAT1 &A, const VECT &eigval_,
  631. const MAT2 &eigvect_,
  632. tol_type_for_qr tol = default_tol_for_qr,
  633. bool compvect = true) {
  634. VECT &eigval = const_cast<VECT &>(eigval_);
  635. MAT2 &eigvect = const_cast<MAT2 &>(eigvect_);
  636. typedef typename linalg_traits<MAT1>::value_type T;
  637. typedef typename number_traits<T>::magnitude_type R;
  638. size_type n = mat_nrows(A), q = 0, p, ite = 0;
  639. if (compvect) gmm::copy(identity_matrix(), eigvect);
  640. if (n == 0) return;
  641. if (n == 1) { eigval[0]=gmm::real(A(0,0)); return; }
  642. dense_matrix<T> Tri(n, n);
  643. gmm::copy(A, Tri);
  644. Householder_tridiagonalization(Tri, eigvect, compvect);
  645. std::vector<R> diag(n);
  646. std::vector<T> sdiag(n);
  647. for (size_type i = 0; i < n; ++i)
  648. { diag[i] = gmm::real(Tri(i, i)); if (i+1 < n) sdiag[i] = Tri(i+1, i); }
  649. symmetric_qr_stop_criterion(diag, sdiag, p, q, tol);
  650. while (q < n) {
  651. sub_interval SUBI(p, n-p-q), SUBJ(0, mat_ncols(eigvect)), SUBK(p, n-p-q);
  652. if (!compvect) SUBK = sub_interval(0,0);
  653. symmetric_Wilkinson_qr_step(sub_vector(diag, SUBI),
  654. sub_vector(sdiag, SUBI),
  655. sub_matrix(eigvect, SUBJ, SUBK), compvect);
  656. symmetric_qr_stop_criterion(diag, sdiag, p, q, tol*R(3));
  657. ++ite;
  658. GMM_ASSERT1(ite < n*100, "QR algorithm failed.");
  659. }
  660. gmm::copy(diag, eigval);
  661. }
  662. template <typename MAT1, typename VECT>
  663. void symmetric_qr_algorithm(const MAT1 &a, VECT &eigval,
  664. tol_type_for_qr tol = default_tol_for_qr) {
  665. dense_matrix<typename linalg_traits<MAT1>::value_type> m(0,0);
  666. symmetric_qr_algorithm(a, eigval, m, tol, false);
  667. }
  668. }
  669. #endif