123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277 |
- #pragma once
- #include <Mathematics/CholeskyDecomposition.h>
- #include <functional>
- namespace WwiseGTE
- {
- template <typename Real>
- class LevenbergMarquardtMinimizer
- {
- public:
-
-
- typedef GVector<Real> DVector;
- typedef GVector<Real> RVector;
- typedef GMatrix<Real> JMatrix;
- typedef GMatrix<Real> JTJMatrix;
- typedef GVector<Real> JTFVector;
- typedef std::function<void(DVector const&, RVector&)> FFunction;
- typedef std::function<void(DVector const&, JMatrix&)> JFunction;
- typedef std::function<void(DVector const&, JTJMatrix&, JTFVector&)> JPlusFunction;
-
- LevenbergMarquardtMinimizer(int numPDimensions, int numFDimensions,
- FFunction const& inFFunction, JFunction const& inJFunction)
- :
- mNumPDimensions(numPDimensions),
- mNumFDimensions(numFDimensions),
- mFFunction(inFFunction),
- mJFunction(inJFunction),
- mF(mNumFDimensions),
- mJ(mNumFDimensions, mNumPDimensions),
- mJTJ(mNumPDimensions, mNumPDimensions),
- mNegJTF(mNumPDimensions),
- mDecomposer(mNumPDimensions),
- mUseJFunction(true)
- {
- LogAssert(mNumPDimensions > 0 && mNumFDimensions > 0, "Invalid dimensions.");
- }
-
- LevenbergMarquardtMinimizer(int numPDimensions, int numFDimensions,
- FFunction const& inFFunction, JPlusFunction const& inJPlusFunction)
- :
- mNumPDimensions(numPDimensions),
- mNumFDimensions(numFDimensions),
- mFFunction(inFFunction),
- mJPlusFunction(inJPlusFunction),
- mF(mNumFDimensions),
- mJ(mNumFDimensions, mNumPDimensions),
- mJTJ(mNumPDimensions, mNumPDimensions),
- mNegJTF(mNumPDimensions),
- mDecomposer(mNumPDimensions),
- mUseJFunction(false)
- {
- LogAssert(mNumPDimensions > 0 && mNumFDimensions > 0, "Invalid dimensions.");
- }
-
- LevenbergMarquardtMinimizer(LevenbergMarquardtMinimizer const&) = delete;
- LevenbergMarquardtMinimizer& operator=(LevenbergMarquardtMinimizer const&) = delete;
- LevenbergMarquardtMinimizer(LevenbergMarquardtMinimizer&&) = delete;
- LevenbergMarquardtMinimizer& operator=(LevenbergMarquardtMinimizer&&) = delete;
- inline int GetNumPDimensions() const { return mNumPDimensions; }
- inline int GetNumFDimensions() const { return mNumFDimensions; }
-
-
-
-
-
- struct Result
- {
- DVector minLocation;
- Real minError;
- Real minErrorDifference;
- Real minUpdateLength;
- size_t numIterations;
- size_t numAdjustments;
- bool converged;
- };
- Result operator()(DVector const& p0, size_t maxIterations,
- Real updateLengthTolerance, Real errorDifferenceTolerance,
- Real lambdaFactor, Real lambdaAdjust, size_t maxAdjustments)
- {
- Result result;
- result.minLocation = p0;
- result.minError = std::numeric_limits<Real>::max();
- result.minErrorDifference = std::numeric_limits<Real>::max();
- result.minUpdateLength = (Real)0;
- result.numIterations = 0;
- result.numAdjustments = 0;
- result.converged = false;
-
-
- if (lambdaFactor <= (Real)0 || lambdaAdjust <= (Real)0)
- {
- maxAdjustments = 1;
- lambdaFactor = (Real)0;
- lambdaAdjust = (Real)1;
- }
-
- updateLengthTolerance = std::max(updateLengthTolerance, (Real)0);
- errorDifferenceTolerance = std::max(errorDifferenceTolerance, (Real)0);
-
- mFFunction(p0, mF);
- result.minError = Dot(mF, mF);
-
- auto pCurrent = p0;
- for (result.numIterations = 1; result.numIterations <= maxIterations; ++result.numIterations)
- {
- std::pair<bool, bool> status;
- DVector pNext;
- for (result.numAdjustments = 0; result.numAdjustments < maxAdjustments; ++result.numAdjustments)
- {
- status = DoIteration(pCurrent, lambdaFactor, updateLengthTolerance,
- errorDifferenceTolerance, pNext, result);
- if (status.first)
- {
-
-
-
-
- return result;
- }
- if (status.second)
- {
-
-
- break;
- }
- lambdaFactor *= lambdaAdjust;
- }
- if (result.numAdjustments < maxAdjustments)
- {
-
-
-
-
- lambdaFactor /= lambdaAdjust;
- }
- else
- {
-
-
-
-
-
- status = DoIteration(pCurrent, lambdaFactor, updateLengthTolerance,
- errorDifferenceTolerance, pNext, result);
- if (status.first)
- {
-
-
-
-
- return result;
- }
- }
- pCurrent = pNext;
- }
- return result;
- }
- private:
- void ComputeLinearSystemInputs(DVector const& pCurrent, Real lambda)
- {
- if (mUseJFunction)
- {
- mJFunction(pCurrent, mJ);
- mJTJ = MultiplyATB(mJ, mJ);
- mNegJTF = -(mF * mJ);
- }
- else
- {
- mJPlusFunction(pCurrent, mJTJ, mNegJTF);
- }
- Real diagonalSum(0);
- for (int i = 0; i < mNumPDimensions; ++i)
- {
- diagonalSum += mJTJ(i, i);
- }
- Real diagonalAdjust = lambda * diagonalSum / static_cast<Real>(mNumPDimensions);
- for (int i = 0; i < mNumPDimensions; ++i)
- {
- mJTJ(i, i) += diagonalAdjust;
- }
- }
-
-
-
-
-
-
- std::pair<bool, bool> DoIteration(DVector const& pCurrent, Real lambdaFactor,
- Real updateLengthTolerance, Real errorDifferenceTolerance, DVector& pNext,
- Result& result)
- {
- ComputeLinearSystemInputs(pCurrent, lambdaFactor);
- if (!mDecomposer.Factor(mJTJ))
- {
-
-
-
-
- return std::make_pair(true, false);
- }
- mDecomposer.SolveLower(mJTJ, mNegJTF);
- mDecomposer.SolveUpper(mJTJ, mNegJTF);
- pNext = pCurrent + mNegJTF;
- mFFunction(pNext, mF);
- Real error = Dot(mF, mF);
- if (error < result.minError)
- {
- result.minErrorDifference = result.minError - error;
- result.minUpdateLength = Length(mNegJTF);
- result.minLocation = pNext;
- result.minError = error;
- if (result.minErrorDifference <= errorDifferenceTolerance
- || result.minUpdateLength <= updateLengthTolerance)
- {
- result.converged = true;
- return std::make_pair(true, true);
- }
- else
- {
- return std::make_pair(false, true);
- }
- }
- else
- {
- return std::make_pair(false, false);
- }
- }
- int mNumPDimensions, mNumFDimensions;
- FFunction mFFunction;
- JFunction mJFunction;
- JPlusFunction mJPlusFunction;
-
- RVector mF;
- JMatrix mJ;
- JTJMatrix mJTJ;
- JTFVector mNegJTF;
- CholeskyDecomposition<Real> mDecomposer;
- bool mUseJFunction;
- };
- }
|