UnsymmetricEigenvalues.h 14 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374
  1. // David Eberly, Geometric Tools, Redmond WA 98052
  2. // Copyright (c) 1998-2020
  3. // Distributed under the Boost Software License, Version 1.0.
  4. // https://www.boost.org/LICENSE_1_0.txt
  5. // https://www.geometrictools.com/License/Boost/LICENSE_1_0.txt
  6. // Version: 4.0.2019.08.13
  7. #pragma once
  8. #include <Mathematics/Math.h>
  9. #include <algorithm>
  10. #include <array>
  11. #include <cstdint>
  12. #include <cstring>
  13. #include <functional>
  14. #include <vector>
  15. // An implementation of the QR algorithm described in "Matrix Computations,
  16. // 2nd edition" by G. H. Golub and C. F. Van Loan, The Johns Hopkins
  17. // University Press, Baltimore MD, Fourth Printing 1993. In particular,
  18. // the implementation is based on Chapter 7 (The Unsymmetric Eigenvalue
  19. // Problem), Section 7.5 (The Practical QR Algorithm).
  20. namespace WwiseGTE
  21. {
  22. template <typename Real>
  23. class UnsymmetricEigenvalues
  24. {
  25. public:
  26. // The solver processes NxN matrices (not necessarily symmetric),
  27. // where N >= 3 ('size' is N) and the matrix is stored in row-major
  28. // order. The maximum number of iterations ('maxIterations') must
  29. // be specified for reducing an upper Hessenberg matrix to an upper
  30. // quasi-triangular matrix (upper triangular matrix of blocks where
  31. // the diagonal blocks are 1x1 or 2x2). The goal is to compute the
  32. // real-valued eigenvalues.
  33. UnsymmetricEigenvalues(int32_t size, uint32_t maxIterations)
  34. :
  35. mSize(0),
  36. mSizeM1(0),
  37. mMaxIterations(0),
  38. mNumEigenvalues(0)
  39. {
  40. if (size >= 3 && maxIterations > 0)
  41. {
  42. mSize = size;
  43. mSizeM1 = size - 1;
  44. mMaxIterations = maxIterations;
  45. mMatrix.resize(size * size);
  46. mX.resize(size);
  47. mV.resize(size);
  48. mScaledV.resize(size);
  49. mW.resize(size);
  50. mFlagStorage.resize(size + 1);
  51. std::fill(mFlagStorage.begin(), mFlagStorage.end(), 0);
  52. mSubdiagonalFlag = &mFlagStorage[1];
  53. mEigenvalues.resize(mSize);
  54. }
  55. }
  56. // A copy of the NxN input is made internally. The order of the
  57. // eigenvalues is specified by sortType: -1 (decreasing), 0 (no
  58. // sorting), or +1 (increasing). When sorted, the eigenvectors are
  59. // ordered accordingly. The return value is the number of iterations
  60. // consumed when convergence occurred, 0xFFFFFFFF when convergence did
  61. // not occur, or 0 when N <= 1 was passed to the constructor.
  62. uint32_t Solve(Real const* input, int32_t sortType)
  63. {
  64. if (mSize > 0)
  65. {
  66. std::copy(input, input + mSize * mSize, mMatrix.begin());
  67. ReduceToUpperHessenberg();
  68. std::array<int, 2> block;
  69. bool found = GetBlock(block);
  70. uint32_t numIterations;
  71. for (numIterations = 0; numIterations < mMaxIterations; ++numIterations)
  72. {
  73. if (found)
  74. {
  75. // Solve the current subproblem.
  76. FrancisQRStep(block[0], block[1] + 1);
  77. // Find another subproblem (if any).
  78. found = GetBlock(block);
  79. }
  80. else
  81. {
  82. break;
  83. }
  84. }
  85. // The matrix is fully uncoupled, upper Hessenberg with 1x1 or
  86. // 2x2 diagonal blocks. Golub and Van Loan call this "upper
  87. // quasi-triangular".
  88. mNumEigenvalues = 0;
  89. std::fill(mEigenvalues.begin(), mEigenvalues.end(), (Real)0);
  90. for (int i = 0; i < mSizeM1; ++i)
  91. {
  92. if (mSubdiagonalFlag[i] == 0)
  93. {
  94. if (mSubdiagonalFlag[i - 1] == 0)
  95. {
  96. // We have a 1x1 block with a real eigenvalue.
  97. mEigenvalues[mNumEigenvalues++] = A(i, i);
  98. }
  99. }
  100. else
  101. {
  102. if (mSubdiagonalFlag[i - 1] == 0 && mSubdiagonalFlag[i + 1] == 0)
  103. {
  104. // We have a 2x2 block that might have real
  105. // eigenvalues.
  106. Real a00 = A(i, i);
  107. Real a01 = A(i, i + 1);
  108. Real a10 = A(i + 1, i);
  109. Real a11 = A(i + 1, i + 1);
  110. Real tr = a00 + a11;
  111. Real det = a00 * a11 - a01 * a10;
  112. Real halfTr = tr * (Real)0.5;
  113. Real discr = halfTr * halfTr - det;
  114. if (discr >= (Real)0)
  115. {
  116. Real rootDiscr = std::sqrt(discr);
  117. mEigenvalues[mNumEigenvalues++] = halfTr - rootDiscr;
  118. mEigenvalues[mNumEigenvalues++] = halfTr + rootDiscr;
  119. }
  120. }
  121. // else:
  122. // The QR iteration failed to converge at this block.
  123. // It must also be the case that
  124. // numIterations == mMaxIterations. TODO: The caller
  125. // will be aware of this when testing the returned
  126. // numIterations. Is there a remedy for such a case?
  127. // This happened with root finding using the companion
  128. // matrix of a polynomial.)
  129. }
  130. }
  131. if (sortType != 0 && mNumEigenvalues > 1)
  132. {
  133. if (sortType > 0)
  134. {
  135. std::sort(mEigenvalues.begin(),
  136. mEigenvalues.begin() + mNumEigenvalues, std::less<Real>());
  137. }
  138. else
  139. {
  140. std::sort(mEigenvalues.begin(),
  141. mEigenvalues.begin() + mNumEigenvalues, std::greater<Real>());
  142. }
  143. }
  144. return numIterations;
  145. }
  146. return 0;
  147. }
  148. // Get the real-valued eigenvalues of the matrix passed to Solve(...).
  149. // The input 'eigenvalues' must have at least N elements.
  150. void GetEigenvalues(uint32_t& numEigenvalues, Real* eigenvalues) const
  151. {
  152. if (mSize > 0)
  153. {
  154. numEigenvalues = mNumEigenvalues;
  155. std::memcpy(eigenvalues, mEigenvalues.data(), numEigenvalues * sizeof(Real));
  156. }
  157. else
  158. {
  159. numEigenvalues = 0;
  160. }
  161. }
  162. private:
  163. // 2D accessors to elements of mMatrix[].
  164. inline Real const& A(int r, int c) const
  165. {
  166. return mMatrix[c + r * mSize];
  167. }
  168. inline Real& A(int r, int c)
  169. {
  170. return mMatrix[c + r * mSize];
  171. }
  172. // Compute the Householder vector for (X[rmin],...,x[rmax]). The
  173. // input vector is stored in mX in the index range [rmin,rmax]. The
  174. // output vector V is stored in mV in the index range [rmin,rmax].
  175. // The scaled vector is S = (-2/Dot(V,V))*V and is stored in mScaledV
  176. // in the index range [rmin,rmax].
  177. void House(int rmin, int rmax)
  178. {
  179. Real length = (Real)0;
  180. for (int r = rmin; r <= rmax; ++r)
  181. {
  182. length += mX[r] * mX[r];
  183. }
  184. length = std::sqrt(length);
  185. if (length != (Real)0)
  186. {
  187. Real sign = (mX[rmin] >= (Real)0 ? (Real)1 : (Real)-1);
  188. Real invDenom = (Real)1 / (mX[rmin] + sign * length);
  189. for (int r = rmin + 1; r <= rmax; ++r)
  190. {
  191. mV[r] = mX[r] * invDenom;
  192. }
  193. }
  194. mV[rmin] = (Real)1;
  195. Real dot = (Real)1;
  196. for (int r = rmin + 1; r <= rmax; ++r)
  197. {
  198. dot += mV[r] * mV[r];
  199. }
  200. Real scale = (Real)-2 / dot;
  201. for (int r = rmin; r <= rmax; ++r)
  202. {
  203. mScaledV[r] = scale * mV[r];
  204. }
  205. }
  206. // Support for replacing matrix A by P^T*A*P, where P is a Householder
  207. // reflection computed using House(...).
  208. void RowHouse(int rmin, int rmax, int cmin, int cmax)
  209. {
  210. for (int c = cmin; c <= cmax; ++c)
  211. {
  212. mW[c] = (Real)0;
  213. for (int r = rmin; r <= rmax; ++r)
  214. {
  215. mW[c] += mScaledV[r] * A(r, c);
  216. }
  217. }
  218. for (int r = rmin; r <= rmax; ++r)
  219. {
  220. for (int c = cmin; c <= cmax; ++c)
  221. {
  222. A(r, c) += mV[r] * mW[c];
  223. }
  224. }
  225. }
  226. void ColHouse(int rmin, int rmax, int cmin, int cmax)
  227. {
  228. for (int r = rmin; r <= rmax; ++r)
  229. {
  230. mW[r] = (Real)0;
  231. for (int c = cmin; c <= cmax; ++c)
  232. {
  233. mW[r] += mScaledV[c] * A(r, c);
  234. }
  235. }
  236. for (int r = rmin; r <= rmax; ++r)
  237. {
  238. for (int c = cmin; c <= cmax; ++c)
  239. {
  240. A(r, c) += mW[r] * mV[c];
  241. }
  242. }
  243. }
  244. void ReduceToUpperHessenberg()
  245. {
  246. for (int c = 0, cp1 = 1; c <= mSize - 3; ++c, ++cp1)
  247. {
  248. for (int r = cp1; r <= mSizeM1; ++r)
  249. {
  250. mX[r] = A(r, c);
  251. }
  252. House(cp1, mSizeM1);
  253. RowHouse(cp1, mSizeM1, c, mSizeM1);
  254. ColHouse(0, mSizeM1, cp1, mSizeM1);
  255. }
  256. }
  257. void FrancisQRStep(int rmin, int rmax)
  258. {
  259. // Apply the double implicit shift step.
  260. int const i0 = rmax - 1, i1 = rmax;
  261. Real a00 = A(i0, i0);
  262. Real a01 = A(i0, i1);
  263. Real a10 = A(i1, i0);
  264. Real a11 = A(i1, i1);
  265. Real tr = a00 + a11;
  266. Real det = a00 * a11 - a01 * a10;
  267. int const j0 = rmin, j1 = j0 + 1, j2 = j1 + 1;
  268. Real b00 = A(j0, j0);
  269. Real b01 = A(j0, j1);
  270. Real b10 = A(j1, j0);
  271. Real b11 = A(j1, j1);
  272. Real b21 = A(j2, j1);
  273. mX[rmin] = b00 * (b00 - tr) + b01 * b10 + det;
  274. mX[rmin + 1] = b10 * (b00 + b11 - tr);
  275. mX[rmin + 2] = b10 * b21;
  276. House(rmin, rmin + 2);
  277. RowHouse(rmin, rmin + 2, rmin, rmax);
  278. ColHouse(rmin, std::min(rmax, rmin + 3), rmin, rmin + 2);
  279. // Apply Householder reflections to restore the matrix to upper
  280. // Hessenberg form.
  281. for (int c = 0, cp1 = 1; c <= mSize - 3; ++c, ++cp1)
  282. {
  283. int kmax = std::min(cp1 + 2, mSizeM1);
  284. for (int r = cp1; r <= kmax; ++r)
  285. {
  286. mX[r] = A(r, c);
  287. }
  288. House(cp1, kmax);
  289. RowHouse(cp1, kmax, c, mSizeM1);
  290. ColHouse(0, mSizeM1, cp1, kmax);
  291. }
  292. }
  293. bool GetBlock(std::array<int, 2>& block)
  294. {
  295. for (int i = 0; i < mSizeM1; ++i)
  296. {
  297. Real a00 = A(i, i);
  298. Real a11 = A(i + 1, i + 1);
  299. Real a21 = A(i + 1, i);
  300. Real sum0 = a00 + a11;
  301. Real sum1 = sum0 + a21;
  302. mSubdiagonalFlag[i] = (sum1 != sum0 ? 1 : 0);
  303. }
  304. for (int i = 0; i < mSizeM1; ++i)
  305. {
  306. if (mSubdiagonalFlag[i] == 1)
  307. {
  308. block = { i, -1 };
  309. while (i < mSizeM1 && mSubdiagonalFlag[i] == 1)
  310. {
  311. block[1] = i++;
  312. }
  313. if (block[1] != block[0])
  314. {
  315. return true;
  316. }
  317. }
  318. }
  319. return false;
  320. }
  321. // The number N of rows and columns of the matrices to be processed.
  322. int32_t mSize, mSizeM1;
  323. // The maximum number of iterations for reducing the tridiagonal
  324. // matrix to a diagonal matrix.
  325. uint32_t mMaxIterations;
  326. // The internal copy of a matrix passed to the solver.
  327. std::vector<Real> mMatrix; // NxN elements
  328. // Temporary storage to compute Householder reflections.
  329. std::vector<Real> mX, mV, mScaledV, mW; // N elements
  330. // Flags about the zeroness of the subdiagonal entries. This is used
  331. // to detect uncoupled submatrices and apply the QR algorithm to the
  332. // corresponding subproblems. The storage is padded on both ends with
  333. // zeros to avoid additional code logic when packing the eigenvalues
  334. // for access by the caller.
  335. std::vector<int> mFlagStorage;
  336. int* mSubdiagonalFlag;
  337. int mNumEigenvalues;
  338. std::vector<Real> mEigenvalues;
  339. };
  340. }