123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374 |
- #pragma once
- #include <Mathematics/Math.h>
- #include <algorithm>
- #include <array>
- #include <cstdint>
- #include <cstring>
- #include <functional>
- #include <vector>
- namespace WwiseGTE
- {
- template <typename Real>
- class UnsymmetricEigenvalues
- {
- public:
-
-
-
-
-
-
-
- UnsymmetricEigenvalues(int32_t size, uint32_t maxIterations)
- :
- mSize(0),
- mSizeM1(0),
- mMaxIterations(0),
- mNumEigenvalues(0)
- {
- if (size >= 3 && maxIterations > 0)
- {
- mSize = size;
- mSizeM1 = size - 1;
- mMaxIterations = maxIterations;
- mMatrix.resize(size * size);
- mX.resize(size);
- mV.resize(size);
- mScaledV.resize(size);
- mW.resize(size);
- mFlagStorage.resize(size + 1);
- std::fill(mFlagStorage.begin(), mFlagStorage.end(), 0);
- mSubdiagonalFlag = &mFlagStorage[1];
- mEigenvalues.resize(mSize);
- }
- }
-
-
-
-
-
-
- uint32_t Solve(Real const* input, int32_t sortType)
- {
- if (mSize > 0)
- {
- std::copy(input, input + mSize * mSize, mMatrix.begin());
- ReduceToUpperHessenberg();
- std::array<int, 2> block;
- bool found = GetBlock(block);
- uint32_t numIterations;
- for (numIterations = 0; numIterations < mMaxIterations; ++numIterations)
- {
- if (found)
- {
-
- FrancisQRStep(block[0], block[1] + 1);
-
- found = GetBlock(block);
- }
- else
- {
- break;
- }
- }
-
-
-
- mNumEigenvalues = 0;
- std::fill(mEigenvalues.begin(), mEigenvalues.end(), (Real)0);
- for (int i = 0; i < mSizeM1; ++i)
- {
- if (mSubdiagonalFlag[i] == 0)
- {
- if (mSubdiagonalFlag[i - 1] == 0)
- {
-
- mEigenvalues[mNumEigenvalues++] = A(i, i);
- }
- }
- else
- {
- if (mSubdiagonalFlag[i - 1] == 0 && mSubdiagonalFlag[i + 1] == 0)
- {
-
-
- Real a00 = A(i, i);
- Real a01 = A(i, i + 1);
- Real a10 = A(i + 1, i);
- Real a11 = A(i + 1, i + 1);
- Real tr = a00 + a11;
- Real det = a00 * a11 - a01 * a10;
- Real halfTr = tr * (Real)0.5;
- Real discr = halfTr * halfTr - det;
- if (discr >= (Real)0)
- {
- Real rootDiscr = std::sqrt(discr);
- mEigenvalues[mNumEigenvalues++] = halfTr - rootDiscr;
- mEigenvalues[mNumEigenvalues++] = halfTr + rootDiscr;
- }
- }
-
-
-
-
-
-
-
-
- }
- }
- if (sortType != 0 && mNumEigenvalues > 1)
- {
- if (sortType > 0)
- {
- std::sort(mEigenvalues.begin(),
- mEigenvalues.begin() + mNumEigenvalues, std::less<Real>());
- }
- else
- {
- std::sort(mEigenvalues.begin(),
- mEigenvalues.begin() + mNumEigenvalues, std::greater<Real>());
- }
- }
- return numIterations;
- }
- return 0;
- }
-
-
- void GetEigenvalues(uint32_t& numEigenvalues, Real* eigenvalues) const
- {
- if (mSize > 0)
- {
- numEigenvalues = mNumEigenvalues;
- std::memcpy(eigenvalues, mEigenvalues.data(), numEigenvalues * sizeof(Real));
- }
- else
- {
- numEigenvalues = 0;
- }
- }
- private:
-
- inline Real const& A(int r, int c) const
- {
- return mMatrix[c + r * mSize];
- }
- inline Real& A(int r, int c)
- {
- return mMatrix[c + r * mSize];
- }
-
-
-
-
-
- void House(int rmin, int rmax)
- {
- Real length = (Real)0;
- for (int r = rmin; r <= rmax; ++r)
- {
- length += mX[r] * mX[r];
- }
- length = std::sqrt(length);
- if (length != (Real)0)
- {
- Real sign = (mX[rmin] >= (Real)0 ? (Real)1 : (Real)-1);
- Real invDenom = (Real)1 / (mX[rmin] + sign * length);
- for (int r = rmin + 1; r <= rmax; ++r)
- {
- mV[r] = mX[r] * invDenom;
- }
- }
- mV[rmin] = (Real)1;
- Real dot = (Real)1;
- for (int r = rmin + 1; r <= rmax; ++r)
- {
- dot += mV[r] * mV[r];
- }
- Real scale = (Real)-2 / dot;
- for (int r = rmin; r <= rmax; ++r)
- {
- mScaledV[r] = scale * mV[r];
- }
- }
-
-
- void RowHouse(int rmin, int rmax, int cmin, int cmax)
- {
- for (int c = cmin; c <= cmax; ++c)
- {
- mW[c] = (Real)0;
- for (int r = rmin; r <= rmax; ++r)
- {
- mW[c] += mScaledV[r] * A(r, c);
- }
- }
- for (int r = rmin; r <= rmax; ++r)
- {
- for (int c = cmin; c <= cmax; ++c)
- {
- A(r, c) += mV[r] * mW[c];
- }
- }
- }
- void ColHouse(int rmin, int rmax, int cmin, int cmax)
- {
- for (int r = rmin; r <= rmax; ++r)
- {
- mW[r] = (Real)0;
- for (int c = cmin; c <= cmax; ++c)
- {
- mW[r] += mScaledV[c] * A(r, c);
- }
- }
- for (int r = rmin; r <= rmax; ++r)
- {
- for (int c = cmin; c <= cmax; ++c)
- {
- A(r, c) += mW[r] * mV[c];
- }
- }
- }
- void ReduceToUpperHessenberg()
- {
- for (int c = 0, cp1 = 1; c <= mSize - 3; ++c, ++cp1)
- {
- for (int r = cp1; r <= mSizeM1; ++r)
- {
- mX[r] = A(r, c);
- }
- House(cp1, mSizeM1);
- RowHouse(cp1, mSizeM1, c, mSizeM1);
- ColHouse(0, mSizeM1, cp1, mSizeM1);
- }
- }
- void FrancisQRStep(int rmin, int rmax)
- {
-
- int const i0 = rmax - 1, i1 = rmax;
- Real a00 = A(i0, i0);
- Real a01 = A(i0, i1);
- Real a10 = A(i1, i0);
- Real a11 = A(i1, i1);
- Real tr = a00 + a11;
- Real det = a00 * a11 - a01 * a10;
- int const j0 = rmin, j1 = j0 + 1, j2 = j1 + 1;
- Real b00 = A(j0, j0);
- Real b01 = A(j0, j1);
- Real b10 = A(j1, j0);
- Real b11 = A(j1, j1);
- Real b21 = A(j2, j1);
- mX[rmin] = b00 * (b00 - tr) + b01 * b10 + det;
- mX[rmin + 1] = b10 * (b00 + b11 - tr);
- mX[rmin + 2] = b10 * b21;
- House(rmin, rmin + 2);
- RowHouse(rmin, rmin + 2, rmin, rmax);
- ColHouse(rmin, std::min(rmax, rmin + 3), rmin, rmin + 2);
-
-
- for (int c = 0, cp1 = 1; c <= mSize - 3; ++c, ++cp1)
- {
- int kmax = std::min(cp1 + 2, mSizeM1);
- for (int r = cp1; r <= kmax; ++r)
- {
- mX[r] = A(r, c);
- }
- House(cp1, kmax);
- RowHouse(cp1, kmax, c, mSizeM1);
- ColHouse(0, mSizeM1, cp1, kmax);
- }
- }
- bool GetBlock(std::array<int, 2>& block)
- {
- for (int i = 0; i < mSizeM1; ++i)
- {
- Real a00 = A(i, i);
- Real a11 = A(i + 1, i + 1);
- Real a21 = A(i + 1, i);
- Real sum0 = a00 + a11;
- Real sum1 = sum0 + a21;
- mSubdiagonalFlag[i] = (sum1 != sum0 ? 1 : 0);
- }
- for (int i = 0; i < mSizeM1; ++i)
- {
- if (mSubdiagonalFlag[i] == 1)
- {
- block = { i, -1 };
- while (i < mSizeM1 && mSubdiagonalFlag[i] == 1)
- {
- block[1] = i++;
- }
- if (block[1] != block[0])
- {
- return true;
- }
- }
- }
- return false;
- }
-
- int32_t mSize, mSizeM1;
-
-
- uint32_t mMaxIterations;
-
- std::vector<Real> mMatrix;
-
- std::vector<Real> mX, mV, mScaledV, mW;
-
-
-
-
-
- std::vector<int> mFlagStorage;
- int* mSubdiagonalFlag;
- int mNumEigenvalues;
- std::vector<Real> mEigenvalues;
- };
- }
|