GaussNewtonMinimizer.h 9.2 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220
  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/CholeskyDecomposition.h>
  9. #include <functional>
  10. // Let F(p) = (F_{0}(p), F_{1}(p), ..., F_{n-1}(p)) be a vector-valued
  11. // function of the parameters p = (p_{0}, p_{1}, ..., p_{m-1}). The
  12. // nonlinear least-squares problem is to minimize the real-valued error
  13. // function E(p) = |F(p)|^2, which is the squared length of F(p).
  14. //
  15. // Let J = dF/dp = [dF_{r}/dp_{c}] denote the Jacobian matrix, which is the
  16. // matrix of first-order partial derivatives of F. The matrix has n rows and
  17. // m columns, and the indexing (r,c) refers to row r and column c. A
  18. // first-order approximation is F(p + d) = F(p) + J(p)d, where d is an m-by-1
  19. // vector with small length. Consequently, an approximation to E is E(p + d)
  20. // = |F(p + d)|^2 = |F(p) + J(p)d|^2. The goal is to choose d to minimize
  21. // |F(p) + J(p)d|^2 and, hopefully, with E(p + d) < E(p). Choosing an initial
  22. // p_{0}, the hope is that the algorithm generates a sequence p_{i} for which
  23. // E(p_{i+1}) < E(p_{i}) and, in the limit, E(p_{j}) approaches the global
  24. // minimum of E. The algorithm is referred to as Gauss-Newton iteration. If
  25. // E does not decrease for a step of the algorithm, one can modify the
  26. // algorithm to the Levenberg-Marquardt iteration. See
  27. // GteLevenbergMarquardtMinimizer.h for a description and an implementation.
  28. //
  29. // For a single Gauss-Newton iteration, we need to choose d to minimize
  30. // |F(p) + J(p)d|^2 where p is fixed. This is a linear least squares problem
  31. // which can be formulated using the normal equations
  32. // (J^T(p)*J(p))*d = -J^T(p)*F(p). The matrix J^T*J is positive semidefinite.
  33. // If it is invertible, then d = -(J^T(p)*J(p))^{-1}*F(p). If it is not
  34. // invertible, some other algorithm must be used to choose d; one option is
  35. // to use gradient descent for the step. A Cholesky decomposition can be
  36. // used to solve the linear system.
  37. //
  38. // Although an implementation can allow the caller to pass an array of
  39. // functions F_{i}(p) and an array of derivatives dF_{r}/dp_{c}, some
  40. // applications might involve a very large n that precludes storing all
  41. // the computed Jacobian matrix entries because of excessive memory
  42. // requirements. In such an application, it is better to compute instead
  43. // the entries of the m-by-m matrix J^T*J and the m-by-1 vector J^T*F.
  44. // Typically, m is small, so the memory requirements are not excessive. Also,
  45. // there might be additional structure to F for which the caller can take
  46. // advantage; for example, 3-tuples of components of F(p) might correspond to
  47. // vectors that can be manipulated using an already existing mathematics
  48. // library. The implementation here supports both approaches.
  49. namespace WwiseGTE
  50. {
  51. template <typename Real>
  52. class GaussNewtonMinimizer
  53. {
  54. public:
  55. // Convenient types for the domain vectors, the range vectors, the
  56. // function F and the Jacobian J.
  57. typedef GVector<Real> DVector; // numPDimensions
  58. typedef GVector<Real> RVector; // numFDimensions
  59. typedef GMatrix<Real> JMatrix; // numFDimensions-by-numPDimensions
  60. typedef GMatrix<Real> JTJMatrix; // numPDimensions-by-numPDimensions
  61. typedef GVector<Real> JTFVector; // numPDimensions
  62. typedef std::function<void(DVector const&, RVector&)> FFunction;
  63. typedef std::function<void(DVector const&, JMatrix&)> JFunction;
  64. typedef std::function<void(DVector const&, JTJMatrix&, JTFVector&)> JPlusFunction;
  65. // Create the minimizer that computes F(p) and J(p) directly.
  66. GaussNewtonMinimizer(int numPDimensions, int numFDimensions,
  67. FFunction const& inFFunction, JFunction const& inJFunction)
  68. :
  69. mNumPDimensions(numPDimensions),
  70. mNumFDimensions(numFDimensions),
  71. mFFunction(inFFunction),
  72. mJFunction(inJFunction),
  73. mF(mNumFDimensions),
  74. mJ(mNumFDimensions, mNumPDimensions),
  75. mJTJ(mNumPDimensions, mNumPDimensions),
  76. mNegJTF(mNumPDimensions),
  77. mDecomposer(mNumPDimensions),
  78. mUseJFunction(true)
  79. {
  80. LogAssert(mNumPDimensions > 0 && mNumFDimensions > 0, "Invalid dimensions.");
  81. }
  82. // Create the minimizer that computes J^T(p)*J(p) and -J(p)*F(p).
  83. GaussNewtonMinimizer(int numPDimensions, int numFDimensions,
  84. FFunction const& inFFunction, JPlusFunction const& inJPlusFunction)
  85. :
  86. mNumPDimensions(numPDimensions),
  87. mNumFDimensions(numFDimensions),
  88. mFFunction(inFFunction),
  89. mJPlusFunction(inJPlusFunction),
  90. mF(mNumFDimensions),
  91. mJ(mNumFDimensions, mNumPDimensions),
  92. mJTJ(mNumPDimensions, mNumPDimensions),
  93. mNegJTF(mNumPDimensions),
  94. mDecomposer(mNumPDimensions),
  95. mUseJFunction(false)
  96. {
  97. LogAssert(mNumPDimensions > 0 && mNumFDimensions > 0, "Invalid dimensions.");
  98. }
  99. // Disallow copy, assignment and move semantics.
  100. GaussNewtonMinimizer(GaussNewtonMinimizer const&) = delete;
  101. GaussNewtonMinimizer& operator=(GaussNewtonMinimizer const&) = delete;
  102. GaussNewtonMinimizer(GaussNewtonMinimizer&&) = delete;
  103. GaussNewtonMinimizer& operator=(GaussNewtonMinimizer&&) = delete;
  104. inline int GetNumPDimensions() const
  105. {
  106. return mNumPDimensions;
  107. }
  108. inline int GetNumFDimensions() const
  109. {
  110. return mNumFDimensions;
  111. }
  112. struct Result
  113. {
  114. DVector minLocation;
  115. Real minError;
  116. Real minErrorDifference;
  117. Real minUpdateLength;
  118. size_t numIterations;
  119. bool converged;
  120. };
  121. Result operator()(DVector const& p0, size_t maxIterations,
  122. Real updateLengthTolerance, Real errorDifferenceTolerance)
  123. {
  124. Result result;
  125. result.minLocation = p0;
  126. result.minError = std::numeric_limits<Real>::max();
  127. result.minErrorDifference = std::numeric_limits<Real>::max();
  128. result.minUpdateLength = (Real)0;
  129. result.numIterations = 0;
  130. result.converged = false;
  131. // As a simple precaution, ensure the tolerances are nonnegative.
  132. updateLengthTolerance = std::max(updateLengthTolerance, (Real)0);
  133. errorDifferenceTolerance = std::max(errorDifferenceTolerance, (Real)0);
  134. // Compute the initial error.
  135. mFFunction(p0, mF);
  136. result.minError = Dot(mF, mF);
  137. // Do the Gauss-Newton iterations.
  138. auto pCurrent = p0;
  139. for (result.numIterations = 1; result.numIterations <= maxIterations; ++result.numIterations)
  140. {
  141. ComputeLinearSystemInputs(pCurrent);
  142. if (!mDecomposer.Factor(mJTJ))
  143. {
  144. // TODO: The matrix mJTJ is positive semi-definite, so the
  145. // failure can occur when mJTJ has a zero eigenvalue in
  146. // which case mJTJ is not invertible. Generate an iterate
  147. // anyway, perhaps using gradient descent?
  148. return result;
  149. }
  150. mDecomposer.SolveLower(mJTJ, mNegJTF);
  151. mDecomposer.SolveUpper(mJTJ, mNegJTF);
  152. auto pNext = pCurrent + mNegJTF;
  153. mFFunction(pNext, mF);
  154. Real error = Dot(mF, mF);
  155. if (error < result.minError)
  156. {
  157. result.minErrorDifference = result.minError - error;
  158. result.minUpdateLength = Length(mNegJTF);
  159. result.minLocation = pNext;
  160. result.minError = error;
  161. if (result.minErrorDifference <= errorDifferenceTolerance
  162. || result.minUpdateLength <= updateLengthTolerance)
  163. {
  164. result.converged = true;
  165. return result;
  166. }
  167. }
  168. pCurrent = pNext;
  169. }
  170. return result;
  171. }
  172. private:
  173. void ComputeLinearSystemInputs(DVector const& pCurrent)
  174. {
  175. if (mUseJFunction)
  176. {
  177. mJFunction(pCurrent, mJ);
  178. mJTJ = MultiplyATB(mJ, mJ);
  179. mNegJTF = -(mF * mJ);
  180. }
  181. else
  182. {
  183. mJPlusFunction(pCurrent, mJTJ, mNegJTF);
  184. }
  185. }
  186. int mNumPDimensions, mNumFDimensions;
  187. FFunction mFFunction;
  188. JFunction mJFunction;
  189. JPlusFunction mJPlusFunction;
  190. // Storage for J^T(p)*J(p) and -J^T(p)*F(p) during the iterations.
  191. RVector mF;
  192. JMatrix mJ;
  193. JTJMatrix mJTJ;
  194. JTFVector mNegJTF;
  195. CholeskyDecomposition<Real> mDecomposer;
  196. bool mUseJFunction;
  197. };
  198. }