GteSymmetricEigensolver.h
Go to the documentation of this file.
1 // David Eberly, Geometric Tools, Redmond WA 98052
2 // Copyright (c) 1998-2017
3 // Distributed under the Boost Software License, Version 1.0.
4 // http://www.boost.org/LICENSE_1_0.txt
5 // http://www.geometrictools.com/License/Boost/LICENSE_1_0.txt
6 // File Version: 3.0.0 (2016/06/19)
7 
8 #pragma once
9 
11 #include <LowLevel/GteWrapper.h>
12 #include <algorithm>
13 #include <cmath>
14 #include <cstring>
15 #include <vector>
16 
17 // The SymmetricEigensolver class is an implementation of Algorithm 8.2.3
18 // (Symmetric QR Algorithm) described in "Matrix Computations, 2nd edition"
19 // by G. H. Golub and C. F. Van Loan, The Johns Hopkins University Press,
20 // Baltimore MD, Fourth Printing 1993. Algorithm 8.2.1 (Householder
21 // Tridiagonalization) is used to reduce matrix A to tridiagonal T.
22 // Algorithm 8.2.2 (Implicit Symmetric QR Step with Wilkinson Shift) is
23 // used for the iterative reduction from tridiagonal to diagonal. If A is
24 // the original matrix, D is the diagonal matrix of eigenvalues, and Q is
25 // the orthogonal matrix of eigenvectors, then theoretically Q^T*A*Q = D.
26 // Numerically, we have errors E = Q^T*A*Q - D. Algorithm 8.2.3 mentions
27 // that one expects |E| is approximately u*|A|, where |M| denotes the
28 // Frobenius norm of M and where u is the unit roundoff for the
29 // floating-point arithmetic: 2^{-23} for 'float', which is FLT_EPSILON
30 // = 1.192092896e-7f, and 2^{-52} for'double', which is DBL_EPSILON
31 // = 2.2204460492503131e-16.
32 //
33 // The condition |a(i,i+1)| <= epsilon*(|a(i,i) + a(i+1,i+1)|) used to
34 // determine when the reduction decouples to smaller problems is implemented
35 // as: sum = |a(i,i)| + |a(i+1,i+1)|; sum + |a(i,i+1)| == sum. The idea is
36 // that the superdiagonal term is small relative to its diagonal neighbors,
37 // and so it is effectively zero. The unit tests have shown that this
38 // interpretation of decoupling is effective.
39 //
40 // The authors suggest that once you have the tridiagonal matrix, a practical
41 // implementation will store the diagonal and superdiagonal entries in linear
42 // arrays, ignoring the theoretically zero values not in the 3-band. This is
43 // good for cache coherence. The authors also suggest storing the Householder
44 // vectors in the lower-triangular portion of the matrix to save memory. The
45 // implementation uses both suggestions.
46 //
47 // For matrices with randomly generated values in [0,1], the unit tests
48 // produce the following information for N-by-N matrices.
49 //
50 // N |A| |E| |E|/|A| iterations
51 // -------------------------------------------
52 // 2 1.2332 5.5511e-17 4.5011e-17 1
53 // 3 2.0024 1.1818e-15 5.9020e-16 5
54 // 4 2.8708 9.9287e-16 3.4585e-16 7
55 // 5 2.9040 2.5958e-15 8.9388e-16 9
56 // 6 4.0427 2.2434e-15 5.5493e-16 12
57 // 7 5.0276 3.2456e-15 6.4555e-16 15
58 // 8 5.4468 6.5626e-15 1.2048e-15 15
59 // 9 6.1510 4.0317e-15 6.5545e-16 17
60 // 10 6.7523 4.9334e-15 7.3062e-16 21
61 // 11 7.1322 7.1347e-15 1.0003e-15 22
62 // 12 7.8324 5.6106e-15 7.1633e-16 24
63 // 13 8.1073 5.1366e-15 6.3357e-16 25
64 // 14 8.6257 8.3496e-15 9.6798e-16 29
65 // 15 9.2603 6.9526e-15 7.5080e-16 31
66 // 16 9.9853 6.5807e-15 6.5904e-16 32
67 // 17 10.5388 8.0931e-15 7.6793e-16 35
68 // 18 11.2377 1.1149e-14 9.9218e-16 38
69 // 19 11.7105 1.0711e-14 9.1470e-16 42
70 // 20 12.2227 1.7723e-14 1.4500e-15 42
71 // 21 12.7411 1.2515e-14 9.8231e-16 47
72 // 22 13.4462 1.2645e-14 9.4046e-16 50
73 // 23 13.9541 1.2899e-14 9.2444e-16 47
74 // 24 14.3298 1.6337e-14 1.1400e-15 53
75 // 25 14.8050 1.4760e-14 9.9701e-16 54
76 // 26 15.3914 1.7076e-14 1.1094e-15 57
77 // 27 15.8430 1.9714e-14 1.2443e-15 60
78 // 28 16.4771 1.7148e-14 1.0407e-15 60
79 // 29 16.9909 1.7309e-14 1.0187e-15 60
80 // 30 17.4456 2.1453e-14 1.2297e-15 64
81 // 31 17.9909 2.2069e-14 1.2267e-15 68
82 //
83 // The eigenvalues and |E|/|A| values were compared to those generated by
84 // Mathematica Version 9.0; Wolfram Research, Inc., Champaign IL, 2012.
85 // The results were all comparable with eigenvalues agreeing to a large
86 // number of decimal places.
87 //
88 // Timing on an Intel (R) Core (TM) i7-3930K CPU @ 3.20 GHZ (in seconds):
89 //
90 // N |E|/|A| iters tridiag QR evecs evec[N] comperr
91 // --------------------------------------------------------------
92 // 512 4.4149e-15 1017 0.180 0.005 1.151 0.848 2.166
93 // 1024 6.1691e-15 1990 1.775 0.031 11.894 12.759 21.179
94 // 2048 8.5108e-15 3849 16.592 0.107 119.744 116.56 212.227
95 //
96 // where iters is the number of QR steps taken, tridiag is the computation
97 // of the Householder reflection vectors, evecs is the composition of
98 // Householder reflections and Givens rotations to obtain the matrix of
99 // eigenvectors, evec[N] is N calls to get the eigenvectors separately, and
100 // comperr is the computation E = Q^T*A*Q - D. The construction of the full
101 // eigenvector matrix is, of course, quite expensive. If you need only a
102 // small number of eigenvectors, use function GetEigenvector(int,Real*).
103 
104 namespace gte
105 {
106 
107 template <typename Real>
109 {
110 public:
111  // The solver processes NxN symmetric matrices, where N > 1 ('size' is N)
112  // and the matrix is stored in row-major order. The maximum number of
113  // iterations ('maxIterations') must be specified for the reduction of a
114  // tridiagonal matrix to a diagonal matrix. The goal is to compute
115  // NxN orthogonal Q and NxN diagonal D for which Q^T*A*Q = D.
116  SymmetricEigensolver(int size, unsigned int maxIterations);
117 
118  // A copy of the NxN symmetric input is made internally. The order of
119  // the eigenvalues is specified by sortType: -1 (decreasing), 0 (no
120  // sorting), or +1 (increasing). When sorted, the eigenvectors are
121  // ordered accordingly. The return value is the number of iterations
122  // consumed when convergence occurred, 0xFFFFFFFF when convergence did
123  // not occur, or 0 when N <= 1 was passed to the constructor.
124  unsigned int Solve(Real const* input, int sortType);
125 
126  // Get the eigenvalues of the matrix passed to Solve(...). The input
127  // 'eigenvalues' must have N elements.
128  void GetEigenvalues(Real* eigenvalues) const;
129 
130  // Accumulate the Householder reflections and Givens rotations to produce
131  // the orthogonal matrix Q for which Q^T*A*Q = D. The input
132  // 'eigenvectors' must be NxN and stored in row-major order.
133  void GetEigenvectors(Real* eigenvectors) const;
134 
135  // With no sorting, when N is odd the matrix returned by GetEigenvectors
136  // is a reflection and when N is even it is a rotation. With sorting
137  // enabled, the type of matrix returned depends on the permutation of
138  // columns. If the permutation has C cycles, the minimum number of column
139  // transpositions is T = N-C. Thus, when C is odd the matrix is a
140  // reflection and when C is even the matrix is a rotation.
141  bool IsRotation() const;
142 
143  // Compute a single eigenvector, which amounts to computing column c
144  // of matrix Q. The reflections and rotations are applied incrementally.
145  // This is useful when you want only a small number of the eigenvectors.
146  void GetEigenvector(int c, Real* eigenvector) const;
147  Real GetEigenvalue(int c) const;
148 
149 private:
150  // Tridiagonalize using Householder reflections. On input, mMatrix is a
151  // copy of the input matrix. On output, the upper-triangular part of
152  // mMatrix including the diagonal stores the tridiagonalization. The
153  // lower-triangular part contains 2/Dot(v,v) that are used in computing
154  // eigenvectors and the part below the subdiagonal stores the essential
155  // parts of the Householder vectors v (the elements of v after the
156  // leading 1-valued component).
157  void Tridiagonalize();
158 
159  // A helper for generating Givens rotation sine and cosine robustly.
160  void GetSinCos(Real u, Real v, Real& cs, Real& sn);
161 
162  // The QR step with implicit shift. Generally, the initial T is unreduced
163  // tridiagonal (all subdiagonal entries are nonzero). If a QR step causes
164  // a superdiagonal entry to become zero, the matrix decouples into a block
165  // diagonal matrix with two tridiagonal blocks. These blocks can be
166  // reduced independently of each other, which allows for parallelization
167  // of the algorithm. The inputs imin and imax identify the subblock of T
168  // to be processed. That block has upper-left element T(imin,imin) and
169  // lower-right element T(imax,imax).
170  void DoQRImplicitShift(int imin, int imax);
171 
172  // Sort the eigenvalues and compute the corresponding permutation of the
173  // indices of the array storing the eigenvalues. The permutation is used
174  // for reordering the eigenvalues and eigenvectors in the calls to
175  // GetEigenvalues(...) and GetEigenvectors(...).
176  void ComputePermutation(int sortType);
177 
178  // The number N of rows and columns of the matrices to be processed.
179  int mSize;
180 
181  // The maximum number of iterations for reducing the tridiagonal mtarix
182  // to a diagonal matrix.
183  unsigned int mMaxIterations;
184 
185  // The internal copy of a matrix passed to the solver. See the comments
186  // about function Tridiagonalize() about what is stored in the matrix.
187  std::vector<Real> mMatrix; // NxN elements
188 
189  // After the initial tridiagonalization by Householder reflections, we no
190  // longer need the full mMatrix. Copy the diagonal and superdiagonal
191  // entries to linear arrays in order to be cache friendly.
192  std::vector<Real> mDiagonal; // N elements
193  std::vector<Real> mSuperdiagonal; // N-1 elements
194 
195  // The Givens rotations used to reduce the initial tridiagonal matrix to
196  // a diagonal matrix. A rotation is the identity with the following
197  // replacement entries: R(index,index) = cs, R(index,index+1) = sn,
198  // R(index+1,index) = -sn, and R(index+1,index+1) = cs. If N is the
199  // matrix size and K is the maximum number of iterations, the maximum
200  // number of Givens rotations is K*(N-1). The maximum amount of memory
201  // is allocated to store these.
203  {
204  GivensRotation();
205  GivensRotation(int inIndex, Real inCs, Real inSn);
206  int index;
207  Real cs, sn;
208  };
209 
210  std::vector<GivensRotation> mGivens; // K*(N-1) elements
211 
212  // When sorting is requested, the permutation associated with the sort is
213  // stored in mPermutation. When sorting is not requested, mPermutation[0]
214  // is set to -1. mVisited is used for finding cycles in the permutation.
215  std::vector<int> mPermutation; // N elements
216  mutable std::vector<int> mVisited; // N elements
217  mutable int mIsRotation; // 1 = rotation, 0 = reflection, -1 = unknown
218 
219  // Temporary storage to compute Householder reflections and to support
220  // sorting of eigenvectors.
221  mutable std::vector<Real> mPVector; // N elements
222  mutable std::vector<Real> mVVector; // N elements
223  mutable std::vector<Real> mWVector; // N elements
224 };
225 
226 
227 template <typename Real>
229  unsigned int maxIterations)
230  :
231  mSize(0),
232  mMaxIterations(0),
233  mIsRotation(-1)
234 {
235  if (size > 1 && maxIterations > 0)
236  {
237  mSize = size;
238  mMaxIterations = maxIterations;
239  mMatrix.resize(size*size);
240  mDiagonal.resize(size);
241  mSuperdiagonal.resize(size - 1);
242  mGivens.reserve(maxIterations * (size - 1));
243  mPermutation.resize(size);
244  mVisited.resize(size);
245  mPVector.resize(size);
246  mVVector.resize(size);
247  mWVector.resize(size);
248  }
249 }
250 
251 template <typename Real>
252 unsigned int SymmetricEigensolver<Real>::Solve(Real const* input,
253  int sortType)
254 {
255  if (mSize > 0)
256  {
257  std::copy(input, input + mSize*mSize, mMatrix.begin());
258  Tridiagonalize();
259 
260  mGivens.clear();
261  for (unsigned int j = 0; j < mMaxIterations; ++j)
262  {
263  int imin = -1, imax = -1;
264  for (int i = mSize - 2; i >= 0; --i)
265  {
266  // When a01 is much smaller than its diagonal neighbors, it is
267  // effectively zero.
268  Real a00 = mDiagonal[i];
269  Real a01 = mSuperdiagonal[i];
270  Real a11 = mDiagonal[i + 1];
271  Real sum = std::abs(a00) + std::abs(a11);
272  if (sum + std::abs(a01) != sum)
273  {
274  if (imax == -1)
275  {
276  imax = i;
277  }
278  imin = i;
279  }
280  else
281  {
282  // The superdiagonal term is effectively zero compared to
283  // the neighboring diagonal terms.
284  if (imin >= 0)
285  {
286  break;
287  }
288  }
289  }
290 
291  if (imax == -1)
292  {
293  // The algorithm has converged.
294  ComputePermutation(sortType);
295  return j;
296  }
297 
298  // Process the lower-right-most unreduced tridiagonal block.
299  DoQRImplicitShift(imin, imax);
300  }
301  return 0xFFFFFFFF;
302  }
303  else
304  {
305  return 0;
306  }
307 }
308 
309 template <typename Real>
310 void SymmetricEigensolver<Real>::GetEigenvalues(Real* eigenvalues) const
311 {
312  if (eigenvalues && mSize > 0)
313  {
314  if (mPermutation[0] >= 0)
315  {
316  // Sorting was requested.
317  for (int i = 0; i < mSize; ++i)
318  {
319  int p = mPermutation[i];
320  eigenvalues[i] = mDiagonal[p];
321  }
322  }
323  else
324  {
325  // Sorting was not requested.
326  size_t numBytes = mSize*sizeof(Real);
327  Memcpy(eigenvalues, &mDiagonal[0], numBytes);
328  }
329  }
330 }
331 
332 template <typename Real>
333 void SymmetricEigensolver<Real>::GetEigenvectors(Real* eigenvectors) const
334 {
335  if (eigenvectors && mSize > 0)
336  {
337  // Start with the identity matrix.
338  std::fill(eigenvectors, eigenvectors + mSize*mSize, (Real)0);
339  for (int d = 0; d < mSize; ++d)
340  {
341  eigenvectors[d + mSize*d] = (Real)1;
342  }
343 
344  // Multiply the Householder reflections using backward accumulation.
345  int r, c;
346  for (int i = mSize - 3, rmin = i + 1; i >= 0; --i, --rmin)
347  {
348  // Copy the v vector and 2/Dot(v,v) from the matrix.
349  Real const* column = &mMatrix[i];
350  Real twoinvvdv = column[mSize*(i + 1)];
351  for (r = 0; r < i + 1; ++r)
352  {
353  mVVector[r] = (Real)0;
354  }
355  mVVector[r] = (Real)1;
356  for (++r; r < mSize; ++r)
357  {
358  mVVector[r] = column[mSize*r];
359  }
360 
361  // Compute the w vector.
362  for (r = 0; r < mSize; ++r)
363  {
364  mWVector[r] = (Real)0;
365  for (c = rmin; c < mSize; ++c)
366  {
367  mWVector[r] += mVVector[c] * eigenvectors[r + mSize*c];
368  }
369  mWVector[r] *= twoinvvdv;
370  }
371 
372  // Update the matrix, Q <- Q - v*w^T.
373  for (r = rmin; r < mSize; ++r)
374  {
375  for (c = 0; c < mSize; ++c)
376  {
377  eigenvectors[c + mSize*r] -= mVVector[r] * mWVector[c];
378  }
379  }
380  }
381 
382  // Multiply the Givens rotations.
383  for (auto const& givens : mGivens)
384  {
385  for (r = 0; r < mSize; ++r)
386  {
387  int j = givens.index + mSize*r;
388  Real& q0 = eigenvectors[j];
389  Real& q1 = eigenvectors[j + 1];
390  Real prd0 = givens.cs * q0 - givens.sn * q1;
391  Real prd1 = givens.sn * q0 + givens.cs * q1;
392  q0 = prd0;
393  q1 = prd1;
394  }
395  }
396 
397  // The number of Householder reflections is N = (mSize-2)*(mSize-1)/2.
398  // If N is even, the product of Householder reflections is a rotation;
399  // otherwise, N is odd and the product is a reflection.
400  mIsRotation = 1 - (((mSize - 2) * (mSize - 1) / 2) & 1);
401 
402  if (mPermutation[0] >= 0)
403  {
404  // Sorting was requested.
405  std::fill(mVisited.begin(), mVisited.end(), 0);
406  for (int i = 0; i < mSize; ++i)
407  {
408  if (mVisited[i] == 0 && mPermutation[i] != i)
409  {
410  // The item starts a cycle with 2 or more elements.
411  int start = i, current = i, j, next;
412  for (j = 0; j < mSize; ++j)
413  {
414  mPVector[j] = eigenvectors[i + mSize*j];
415  }
416  while ((next = mPermutation[current]) != start)
417  {
418  mIsRotation = 1 - mIsRotation;
419  mVisited[current] = 1;
420  for (j = 0; j < mSize; ++j)
421  {
422  eigenvectors[current + mSize*j] =
423  eigenvectors[next + mSize*j];
424  }
425  current = next;
426  }
427  mVisited[current] = 1;
428  for (j = 0; j < mSize; ++j)
429  {
430  eigenvectors[current + mSize*j] = mPVector[j];
431  }
432  }
433  }
434  }
435  }
436 }
437 
438 template <typename Real>
440 {
441  if (mSize > 0)
442  {
443  if (mIsRotation == -1)
444  {
445  // Without sorting, the matrix is a rotation when size is even.
446  // The number of Householder reflections is
447  // N = (mSize-2)*(mSize-1)/2. If N is even, the product of
448  // Householder reflections is a rotation; otherwise, N is odd and
449  // the product is a reflection.
450  mIsRotation = 1 - (((mSize - 2) * (mSize - 1) / 2) & 1);
451 
452  if (mPermutation[0] >= 0)
453  {
454  // With sorting, the matrix is a rotation when the number of
455  // cycles in the permutation is even.
456  std::fill(mVisited.begin(), mVisited.end(), 0);
457  for (int i = 0; i < mSize; ++i)
458  {
459  if (mVisited[i] == 0 && mPermutation[i] != i)
460  {
461  // The item starts a cycle with 2 or more elements.
462  int start = i, current = i, next;
463  while ((next = mPermutation[current]) != start)
464  {
465  mIsRotation = 1 - mIsRotation;
466  mVisited[current] = 1;
467  current = next;
468  }
469  mVisited[current] = 1;
470  }
471  }
472  }
473  }
474  return mIsRotation == 1;
475  }
476  else
477  {
478  return false;
479  }
480 }
481 
482 template <typename Real>
483 void SymmetricEigensolver<Real>::GetEigenvector(int c, Real* eigenvector)
484 const
485 {
486  if (0 <= c && c < mSize)
487  {
488  // y = H*x, then x and y are swapped for the next H
489  Real* x = eigenvector;
490  Real* y = &mPVector[0];
491 
492  // Start with the Euclidean basis vector.
493  memset(x, 0, mSize*sizeof(Real));
494  if (mPermutation[0] >= 0)
495  {
496  // Sorting was requested.
497  x[mPermutation[c]] = (Real)1;
498  }
499  else
500  {
501  x[c] = (Real)1;
502  }
503 
504  // Apply the Givens rotations.
505  for (auto const& givens : gte::reverse(mGivens))
506  {
507  Real& xr = x[givens.index];
508  Real& xrp1 = x[givens.index + 1];
509  Real tmp0 = givens.cs * xr + givens.sn * xrp1;
510  Real tmp1 = -givens.sn * xr + givens.cs * xrp1;
511  xr = tmp0;
512  xrp1 = tmp1;
513  }
514 
515  // Apply the Householder reflections.
516  for (int i = mSize - 3; i >= 0; --i)
517  {
518  // Get the Householder vector v.
519  Real const* column = &mMatrix[i];
520  Real twoinvvdv = column[mSize*(i + 1)];
521  int r;
522  for (r = 0; r < i + 1; ++r)
523  {
524  y[r] = x[r];
525  }
526 
527  // Compute s = Dot(x,v) * 2/v^T*v.
528  Real s = x[r]; // r = i+1, v[i+1] = 1
529  for (int j = r + 1; j < mSize; ++j)
530  {
531  s += x[j] * column[mSize*j];
532  }
533  s *= twoinvvdv;
534 
535  y[r] = x[r] - s; // v[i+1] = 1
536 
537  // Compute the remaining components of y.
538  for (++r; r < mSize; ++r)
539  {
540  y[r] = x[r] - s * column[mSize*r];
541  }
542 
543  std::swap(x, y);
544  }
545 
546  // The final product is stored in x.
547  if (x != eigenvector)
548  {
549  size_t numBytes = mSize*sizeof(Real);
550  Memcpy(eigenvector, x, numBytes);
551  }
552  }
553 }
554 
555 template <typename Real>
557 {
558  if (mSize > 0)
559  {
560  if (mPermutation[0] >= 0)
561  {
562  // Sorting was requested.
563  return mDiagonal[mPermutation[c]];
564  }
565  else
566  {
567  // Sorting was not requested.
568  return mDiagonal[c];
569  }
570  }
571  else
572  {
573  return std::numeric_limits<Real>::max();
574  }
575 }
576 
577 template <typename Real>
579 {
580  int r, c;
581  for (int i = 0, ip1 = 1; i < mSize - 2; ++i, ++ip1)
582  {
583  // Compute the Householder vector. Read the initial vector from the
584  // row of the matrix.
585  Real length = (Real)0;
586  for (r = 0; r < ip1; ++r)
587  {
588  mVVector[r] = (Real)0;
589  }
590  for (r = ip1; r < mSize; ++r)
591  {
592  Real vr = mMatrix[r + mSize*i];
593  mVVector[r] = vr;
594  length += vr * vr;
595  }
596  Real vdv = (Real)1;
597  length = sqrt(length);
598  if (length >(Real)0)
599  {
600  Real& v1 = mVVector[ip1];
601  Real sgn = (v1 >= (Real)0 ? (Real)1 : (Real)-1);
602  Real invDenom = ((Real)1) / (v1 + sgn * length);
603  v1 = (Real)1;
604  for (r = ip1 + 1; r < mSize; ++r)
605  {
606  Real& vr = mVVector[r];
607  vr *= invDenom;
608  vdv += vr * vr;
609  }
610  }
611 
612  // Compute the rank-1 offsets v*w^T and w*v^T.
613  Real invvdv = (Real)1 / vdv;
614  Real twoinvvdv = invvdv * (Real)2;
615  Real pdvtvdv = (Real)0;
616  for (r = i; r < mSize; ++r)
617  {
618  mPVector[r] = (Real)0;
619  for (c = i; c < r; ++c)
620  {
621  mPVector[r] += mMatrix[r + mSize*c] * mVVector[c];
622  }
623  for (; c < mSize; ++c)
624  {
625  mPVector[r] += mMatrix[c + mSize*r] * mVVector[c];
626  }
627  mPVector[r] *= twoinvvdv;
628  pdvtvdv += mPVector[r] * mVVector[r];
629  }
630 
631  pdvtvdv *= invvdv;
632  for (r = i; r < mSize; ++r)
633  {
634  mWVector[r] = mPVector[r] - pdvtvdv * mVVector[r];
635  }
636 
637  // Update the input matrix.
638  for (r = i; r < mSize; ++r)
639  {
640  Real vr = mVVector[r];
641  Real wr = mWVector[r];
642  Real offset = vr * wr * (Real)2;
643  mMatrix[r + mSize*r] -= offset;
644  for (c = r + 1; c < mSize; ++c)
645  {
646  offset = vr * mWVector[c] + wr * mVVector[c];
647  mMatrix[c + mSize*r] -= offset;
648  }
649  }
650 
651  // Copy the vector to column i of the matrix. The 0-valued components
652  // at indices 0 through i are not stored. The 1-valued component at
653  // index i+1 is also not stored; instead, the quantity 2/Dot(v,v) is
654  // stored for use in eigenvector construction. That construction must
655  // take into account the implied components that are not stored.
656  mMatrix[i + mSize*ip1] = twoinvvdv;
657  for (r = ip1 + 1; r < mSize; ++r)
658  {
659  mMatrix[i + mSize*r] = mVVector[r];
660  }
661  }
662 
663  // Copy the diagonal and subdiagonal entries for cache coherence in
664  // the QR iterations.
665  int k, ksup = mSize - 1, index = 0, delta = mSize + 1;
666  for (k = 0; k < ksup; ++k, index += delta)
667  {
668  mDiagonal[k] = mMatrix[index];
669  mSuperdiagonal[k] = mMatrix[index + 1];
670  }
671  mDiagonal[k] = mMatrix[index];
672 }
673 
674 template <typename Real>
675 void SymmetricEigensolver<Real>::GetSinCos(Real x, Real y, Real& cs, Real& sn)
676 {
677  // Solves sn*x + cs*y = 0 robustly.
678  Real tau;
679  if (y != (Real)0)
680  {
681  if (std::abs(y) > std::abs(x))
682  {
683  tau = -x / y;
684  sn = ((Real)1) / sqrt(((Real)1) + tau*tau);
685  cs = sn * tau;
686  }
687  else
688  {
689  tau = -y / x;
690  cs = ((Real)1) / sqrt(((Real)1) + tau*tau);
691  sn = cs * tau;
692  }
693  }
694  else
695  {
696  cs = (Real)1;
697  sn = (Real)0;
698  }
699 }
700 
701 template <typename Real>
703 {
704  // The implicit shift. Compute the eigenvalue u of the lower-right 2x2
705  // block that is closer to a11.
706  Real a00 = mDiagonal[imax];
707  Real a01 = mSuperdiagonal[imax];
708  Real a11 = mDiagonal[imax + 1];
709  Real dif = (a00 - a11) * (Real)0.5;
710  Real sgn = (dif >= (Real)0 ? (Real)1 : (Real)-1);
711  Real a01sqr = a01 * a01;
712  Real u = a11 - a01sqr / (dif + sgn*sqrt(dif*dif + a01sqr));
713  Real x = mDiagonal[imin] - u;
714  Real y = mSuperdiagonal[imin];
715 
716  Real a12, a22, a23, tmp11, tmp12, tmp21, tmp22, cs, sn;
717  Real a02 = (Real)0;
718  int i0 = imin - 1, i1 = imin, i2 = imin + 1;
719  for (; i1 <= imax; ++i0, ++i1, ++i2)
720  {
721  // Compute the Givens rotation and save it for use in computing the
722  // eigenvectors.
723  GetSinCos(x, y, cs, sn);
724  mGivens.push_back(GivensRotation(i1, cs, sn));
725 
726  // Update the tridiagonal matrix. This amounts to updating a 4x4
727  // subblock,
728  // b00 b01 b02 b03
729  // b01 b11 b12 b13
730  // b02 b12 b22 b23
731  // b03 b13 b23 b33
732  // The four corners (b00, b03, b33) do not change values. The
733  // The interior block {{b11,b12},{b12,b22}} is updated on each pass.
734  // For the first pass, the b0c values are out of range, so only
735  // the values (b13, b23) change. For the last pass, the br3 values
736  // are out of range, so only the values (b01, b02) change. For
737  // passes between first and last, the values (b01, b02, b13, b23)
738  // change.
739  if (i1 > imin)
740  {
741  mSuperdiagonal[i0] = cs*mSuperdiagonal[i0] - sn*a02;
742  }
743 
744  a11 = mDiagonal[i1];
745  a12 = mSuperdiagonal[i1];
746  a22 = mDiagonal[i2];
747  tmp11 = cs*a11 - sn*a12;
748  tmp12 = cs*a12 - sn*a22;
749  tmp21 = sn*a11 + cs*a12;
750  tmp22 = sn*a12 + cs*a22;
751  mDiagonal[i1] = cs*tmp11 - sn*tmp12;
752  mSuperdiagonal[i1] = sn*tmp11 + cs*tmp12;
753  mDiagonal[i2] = sn*tmp21 + cs*tmp22;
754 
755  if (i1 < imax)
756  {
757  a23 = mSuperdiagonal[i2];
758  a02 = -sn*a23;
759  mSuperdiagonal[i2] = cs*a23;
760 
761  // Update the parameters for the next Givens rotation.
762  x = mSuperdiagonal[i1];
763  y = a02;
764  }
765  }
766 }
767 
768 template <typename Real>
770 {
771  mIsRotation = -1;
772 
773  if (sortType == 0)
774  {
775  // Set a flag for GetEigenvalues() and GetEigenvectors() to know
776  // that sorted output was not requested.
777  mPermutation[0] = -1;
778  return;
779  }
780 
781  // Compute the permutation induced by sorting. Initially, we start with
782  // the identity permutation I = (0,1,...,N-1).
783  struct SortItem
784  {
785  Real eigenvalue;
786  int index;
787  };
788 
789  std::vector<SortItem> items(mSize);
790  int i;
791  for (i = 0; i < mSize; ++i)
792  {
793  items[i].eigenvalue = mDiagonal[i];
794  items[i].index = i;
795  }
796 
797  if (sortType > 0)
798  {
799  std::sort(items.begin(), items.end(),
800  [](SortItem const& item0, SortItem const& item1)
801  {
802  return item0.eigenvalue < item1.eigenvalue;
803  }
804  );
805  }
806  else
807  {
808  std::sort(items.begin(), items.end(),
809  [](SortItem const& item0, SortItem const& item1)
810  {
811  return item0.eigenvalue > item1.eigenvalue;
812  }
813  );
814  }
815 
816  i = 0;
817  for (auto const& item : items)
818  {
819  mPermutation[i++] = item.index;
820  }
821 
822  // GetEigenvectors() has nontrivial code for computing the orthogonal Q
823  // from the reflections and rotations. To avoid complicating the code
824  // further when sorting is requested, Q is computed as in the unsorted
825  // case. We then need to swap columns of Q to be consistent with the
826  // sorting of the eigenvalues. To minimize copying due to column swaps,
827  // we use permutation P. The minimum number of transpositions to obtain
828  // P from I is N minus the number of cycles of P. Each cycle is reordered
829  // with a minimum number of transpositions; that is, the eigenitems are
830  // cyclically swapped, leading to a minimum amount of copying. For
831  // example, if there is a cycle i0 -> i1 -> i2 -> i3, then the copying is
832  // save = eigenitem[i0];
833  // eigenitem[i1] = eigenitem[i2];
834  // eigenitem[i2] = eigenitem[i3];
835  // eigenitem[i3] = save;
836 }
837 
838 template <typename Real>
840 {
841  // No default initialization for fast creation of std::vector of objects
842  // of this type.
843 }
844 
845 template <typename Real>
847  Real inCs, Real inSn)
848  :
849  index(inIndex),
850  cs(inCs),
851  sn(inSn)
852 {
853 }
854 
855 
856 }
SymmetricEigensolver(int size, unsigned int maxIterations)
void GetSinCos(Real u, Real v, Real &cs, Real &sn)
ReversalObject< ReverseIterator > reverse(Iterable &&range)
gte::BSNumber< UIntegerType > abs(gte::BSNumber< UIntegerType > const &number)
Definition: GteBSNumber.h:966
GLuint start
Definition: glcorearb.h:470
void GetEigenvector(int c, Real *eigenvector) const
GLfloat GLfloat v1
Definition: glcorearb.h:812
void GetEigenvalues(Real *eigenvalues) const
void DoQRImplicitShift(int imin, int imax)
GLint GLenum GLint x
Definition: glcorearb.h:404
GLsizeiptr size
Definition: glcorearb.h:659
GLenum GLenum GLsizei void GLsizei void * column
Definition: glext.h:2725
void GetEigenvectors(Real *eigenvectors) const
const GLubyte * c
Definition: glext.h:11671
std::vector< GivensRotation > mGivens
unsigned int Solve(Real const *input, int sortType)
GLuint GLsizei GLsizei * length
Definition: glcorearb.h:790
const GLdouble * v
Definition: glcorearb.h:832
GLboolean r
Definition: glcorearb.h:1217
GLdouble s
Definition: glext.h:231
GLenum GLenum GLenum input
Definition: glext.h:9913
GLuint index
Definition: glcorearb.h:781
void Memcpy(void *target, void const *source, size_t count)
Definition: GteWrapper.cpp:16
GLintptr offset
Definition: glcorearb.h:660
GLfloat GLfloat p
Definition: glext.h:11668
GLint y
Definition: glcorearb.h:98


geometric_tools_engine
Author(s): Yijiang Huang
autogenerated on Thu Jul 18 2019 04:00:01