polynomial.cc 13 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389
  1. // Ceres Solver - A fast non-linear least squares minimizer
  2. // Copyright 2015 Google Inc. All rights reserved.
  3. // http://ceres-solver.org/
  4. //
  5. // Redistribution and use in source and binary forms, with or without
  6. // modification, are permitted provided that the following conditions are met:
  7. //
  8. // * Redistributions of source code must retain the above copyright notice,
  9. // this list of conditions and the following disclaimer.
  10. // * Redistributions in binary form must reproduce the above copyright notice,
  11. // this list of conditions and the following disclaimer in the documentation
  12. // and/or other materials provided with the distribution.
  13. // * Neither the name of Google Inc. nor the names of its contributors may be
  14. // used to endorse or promote products derived from this software without
  15. // specific prior written permission.
  16. //
  17. // THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
  18. // AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
  19. // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
  20. // ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE
  21. // LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
  22. // CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
  23. // SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
  24. // INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
  25. // CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
  26. // ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
  27. // POSSIBILITY OF SUCH DAMAGE.
  28. //
  29. // Author: moll.markus@arcor.de (Markus Moll)
  30. // sameeragarwal@google.com (Sameer Agarwal)
  31. #include "ceres/polynomial.h"
  32. #include <cmath>
  33. #include <cstddef>
  34. #include <vector>
  35. #include "Eigen/Dense"
  36. #include "ceres/function_sample.h"
  37. #include "ceres/internal/export.h"
  38. #include "glog/logging.h"
  39. namespace ceres::internal {
  40. namespace {
  41. // Balancing function as described by B. N. Parlett and C. Reinsch,
  42. // "Balancing a Matrix for Calculation of Eigenvalues and Eigenvectors".
  43. // In: Numerische Mathematik, Volume 13, Number 4 (1969), 293-304,
  44. // Springer Berlin / Heidelberg. DOI: 10.1007/BF02165404
  45. void BalanceCompanionMatrix(Matrix* companion_matrix_ptr) {
  46. CHECK(companion_matrix_ptr != nullptr);
  47. Matrix& companion_matrix = *companion_matrix_ptr;
  48. Matrix companion_matrix_offdiagonal = companion_matrix;
  49. companion_matrix_offdiagonal.diagonal().setZero();
  50. const int degree = companion_matrix.rows();
  51. // gamma <= 1 controls how much a change in the scaling has to
  52. // lower the 1-norm of the companion matrix to be accepted.
  53. //
  54. // gamma = 1 seems to lead to cycles (numerical issues?), so
  55. // we set it slightly lower.
  56. const double gamma = 0.9;
  57. // Greedily scale row/column pairs until there is no change.
  58. bool scaling_has_changed;
  59. do {
  60. scaling_has_changed = false;
  61. for (int i = 0; i < degree; ++i) {
  62. const double row_norm = companion_matrix_offdiagonal.row(i).lpNorm<1>();
  63. const double col_norm = companion_matrix_offdiagonal.col(i).lpNorm<1>();
  64. // Decompose row_norm/col_norm into mantissa * 2^exponent,
  65. // where 0.5 <= mantissa < 1. Discard mantissa (return value
  66. // of frexp), as only the exponent is needed.
  67. int exponent = 0;
  68. std::frexp(row_norm / col_norm, &exponent);
  69. exponent /= 2;
  70. if (exponent != 0) {
  71. const double scaled_col_norm = std::ldexp(col_norm, exponent);
  72. const double scaled_row_norm = std::ldexp(row_norm, -exponent);
  73. if (scaled_col_norm + scaled_row_norm < gamma * (col_norm + row_norm)) {
  74. // Accept the new scaling. (Multiplication by powers of 2 should not
  75. // introduce rounding errors (ignoring non-normalized numbers and
  76. // over- or underflow))
  77. scaling_has_changed = true;
  78. companion_matrix_offdiagonal.row(i) *= std::ldexp(1.0, -exponent);
  79. companion_matrix_offdiagonal.col(i) *= std::ldexp(1.0, exponent);
  80. }
  81. }
  82. }
  83. } while (scaling_has_changed);
  84. companion_matrix_offdiagonal.diagonal() = companion_matrix.diagonal();
  85. companion_matrix = companion_matrix_offdiagonal;
  86. VLOG(3) << "Balanced companion matrix is\n" << companion_matrix;
  87. }
  88. void BuildCompanionMatrix(const Vector& polynomial,
  89. Matrix* companion_matrix_ptr) {
  90. CHECK(companion_matrix_ptr != nullptr);
  91. Matrix& companion_matrix = *companion_matrix_ptr;
  92. const int degree = polynomial.size() - 1;
  93. companion_matrix.resize(degree, degree);
  94. companion_matrix.setZero();
  95. companion_matrix.diagonal(-1).setOnes();
  96. companion_matrix.col(degree - 1) = -polynomial.reverse().head(degree);
  97. }
  98. // Remove leading terms with zero coefficients.
  99. Vector RemoveLeadingZeros(const Vector& polynomial_in) {
  100. int i = 0;
  101. while (i < (polynomial_in.size() - 1) && polynomial_in(i) == 0.0) {
  102. ++i;
  103. }
  104. return polynomial_in.tail(polynomial_in.size() - i);
  105. }
  106. void FindLinearPolynomialRoots(const Vector& polynomial,
  107. Vector* real,
  108. Vector* imaginary) {
  109. CHECK_EQ(polynomial.size(), 2);
  110. if (real != nullptr) {
  111. real->resize(1);
  112. (*real)(0) = -polynomial(1) / polynomial(0);
  113. }
  114. if (imaginary != nullptr) {
  115. imaginary->setZero(1);
  116. }
  117. }
  118. void FindQuadraticPolynomialRoots(const Vector& polynomial,
  119. Vector* real,
  120. Vector* imaginary) {
  121. CHECK_EQ(polynomial.size(), 3);
  122. const double a = polynomial(0);
  123. const double b = polynomial(1);
  124. const double c = polynomial(2);
  125. const double D = b * b - 4 * a * c;
  126. const double sqrt_D = sqrt(fabs(D));
  127. if (real != nullptr) {
  128. real->setZero(2);
  129. }
  130. if (imaginary != nullptr) {
  131. imaginary->setZero(2);
  132. }
  133. // Real roots.
  134. if (D >= 0) {
  135. if (real != nullptr) {
  136. // Stable quadratic roots according to BKP Horn.
  137. // http://people.csail.mit.edu/bkph/articles/Quadratics.pdf
  138. if (b >= 0) {
  139. (*real)(0) = (-b - sqrt_D) / (2.0 * a);
  140. (*real)(1) = (2.0 * c) / (-b - sqrt_D);
  141. } else {
  142. (*real)(0) = (2.0 * c) / (-b + sqrt_D);
  143. (*real)(1) = (-b + sqrt_D) / (2.0 * a);
  144. }
  145. }
  146. return;
  147. }
  148. // Use the normal quadratic formula for the complex case.
  149. if (real != nullptr) {
  150. (*real)(0) = -b / (2.0 * a);
  151. (*real)(1) = -b / (2.0 * a);
  152. }
  153. if (imaginary != nullptr) {
  154. (*imaginary)(0) = sqrt_D / (2.0 * a);
  155. (*imaginary)(1) = -sqrt_D / (2.0 * a);
  156. }
  157. }
  158. } // namespace
  159. bool FindPolynomialRoots(const Vector& polynomial_in,
  160. Vector* real,
  161. Vector* imaginary) {
  162. if (polynomial_in.size() == 0) {
  163. LOG(ERROR) << "Invalid polynomial of size 0 passed to FindPolynomialRoots";
  164. return false;
  165. }
  166. Vector polynomial = RemoveLeadingZeros(polynomial_in);
  167. const int degree = polynomial.size() - 1;
  168. VLOG(3) << "Input polynomial: " << polynomial_in.transpose();
  169. if (polynomial.size() != polynomial_in.size()) {
  170. VLOG(3) << "Trimmed polynomial: " << polynomial.transpose();
  171. }
  172. // Is the polynomial constant?
  173. if (degree == 0) {
  174. LOG(WARNING) << "Trying to extract roots from a constant "
  175. << "polynomial in FindPolynomialRoots";
  176. // We return true with no roots, not false, as if the polynomial is constant
  177. // it is correct that there are no roots. It is not the case that they were
  178. // there, but that we have failed to extract them.
  179. return true;
  180. }
  181. // Linear
  182. if (degree == 1) {
  183. FindLinearPolynomialRoots(polynomial, real, imaginary);
  184. return true;
  185. }
  186. // Quadratic
  187. if (degree == 2) {
  188. FindQuadraticPolynomialRoots(polynomial, real, imaginary);
  189. return true;
  190. }
  191. // The degree is now known to be at least 3. For cubic or higher
  192. // roots we use the method of companion matrices.
  193. // Divide by leading term
  194. const double leading_term = polynomial(0);
  195. polynomial /= leading_term;
  196. // Build and balance the companion matrix to the polynomial.
  197. Matrix companion_matrix(degree, degree);
  198. BuildCompanionMatrix(polynomial, &companion_matrix);
  199. BalanceCompanionMatrix(&companion_matrix);
  200. // Find its (complex) eigenvalues.
  201. Eigen::EigenSolver<Matrix> solver(companion_matrix, false);
  202. if (solver.info() != Eigen::Success) {
  203. LOG(ERROR) << "Failed to extract eigenvalues from companion matrix.";
  204. return false;
  205. }
  206. // Output roots
  207. if (real != nullptr) {
  208. *real = solver.eigenvalues().real();
  209. } else {
  210. LOG(WARNING) << "nullptr pointer passed as real argument to "
  211. << "FindPolynomialRoots. Real parts of the roots will not "
  212. << "be returned.";
  213. }
  214. if (imaginary != nullptr) {
  215. *imaginary = solver.eigenvalues().imag();
  216. }
  217. return true;
  218. }
  219. Vector DifferentiatePolynomial(const Vector& polynomial) {
  220. const int degree = polynomial.rows() - 1;
  221. CHECK_GE(degree, 0);
  222. // Degree zero polynomials are constants, and their derivative does
  223. // not result in a smaller degree polynomial, just a degree zero
  224. // polynomial with value zero.
  225. if (degree == 0) {
  226. return Eigen::VectorXd::Zero(1);
  227. }
  228. Vector derivative(degree);
  229. for (int i = 0; i < degree; ++i) {
  230. derivative(i) = (degree - i) * polynomial(i);
  231. }
  232. return derivative;
  233. }
  234. void MinimizePolynomial(const Vector& polynomial,
  235. const double x_min,
  236. const double x_max,
  237. double* optimal_x,
  238. double* optimal_value) {
  239. // Find the minimum of the polynomial at the two ends.
  240. //
  241. // We start by inspecting the middle of the interval. Technically
  242. // this is not needed, but we do this to make this code as close to
  243. // the minFunc package as possible.
  244. *optimal_x = (x_min + x_max) / 2.0;
  245. *optimal_value = EvaluatePolynomial(polynomial, *optimal_x);
  246. const double x_min_value = EvaluatePolynomial(polynomial, x_min);
  247. if (x_min_value < *optimal_value) {
  248. *optimal_value = x_min_value;
  249. *optimal_x = x_min;
  250. }
  251. const double x_max_value = EvaluatePolynomial(polynomial, x_max);
  252. if (x_max_value < *optimal_value) {
  253. *optimal_value = x_max_value;
  254. *optimal_x = x_max;
  255. }
  256. // If the polynomial is linear or constant, we are done.
  257. if (polynomial.rows() <= 2) {
  258. return;
  259. }
  260. const Vector derivative = DifferentiatePolynomial(polynomial);
  261. Vector roots_real;
  262. if (!FindPolynomialRoots(derivative, &roots_real, nullptr)) {
  263. LOG(WARNING) << "Unable to find the critical points of "
  264. << "the interpolating polynomial.";
  265. return;
  266. }
  267. // This is a bit of an overkill, as some of the roots may actually
  268. // have a complex part, but its simpler to just check these values.
  269. for (int i = 0; i < roots_real.rows(); ++i) {
  270. const double root = roots_real(i);
  271. if ((root < x_min) || (root > x_max)) {
  272. continue;
  273. }
  274. const double value = EvaluatePolynomial(polynomial, root);
  275. if (value < *optimal_value) {
  276. *optimal_value = value;
  277. *optimal_x = root;
  278. }
  279. }
  280. }
  281. Vector FindInterpolatingPolynomial(const std::vector<FunctionSample>& samples) {
  282. const int num_samples = samples.size();
  283. int num_constraints = 0;
  284. for (int i = 0; i < num_samples; ++i) {
  285. if (samples[i].value_is_valid) {
  286. ++num_constraints;
  287. }
  288. if (samples[i].gradient_is_valid) {
  289. ++num_constraints;
  290. }
  291. }
  292. const int degree = num_constraints - 1;
  293. Matrix lhs = Matrix::Zero(num_constraints, num_constraints);
  294. Vector rhs = Vector::Zero(num_constraints);
  295. int row = 0;
  296. for (int i = 0; i < num_samples; ++i) {
  297. const FunctionSample& sample = samples[i];
  298. if (sample.value_is_valid) {
  299. for (int j = 0; j <= degree; ++j) {
  300. lhs(row, j) = pow(sample.x, degree - j);
  301. }
  302. rhs(row) = sample.value;
  303. ++row;
  304. }
  305. if (sample.gradient_is_valid) {
  306. for (int j = 0; j < degree; ++j) {
  307. lhs(row, j) = (degree - j) * pow(sample.x, degree - j - 1);
  308. }
  309. rhs(row) = sample.gradient;
  310. ++row;
  311. }
  312. }
  313. // TODO(sameeragarwal): This is a hack.
  314. // https://github.com/ceres-solver/ceres-solver/issues/248
  315. Eigen::FullPivLU<Matrix> lu(lhs);
  316. return lu.setThreshold(0.0).solve(rhs);
  317. }
  318. void MinimizeInterpolatingPolynomial(const std::vector<FunctionSample>& samples,
  319. double x_min,
  320. double x_max,
  321. double* optimal_x,
  322. double* optimal_value) {
  323. const Vector polynomial = FindInterpolatingPolynomial(samples);
  324. MinimizePolynomial(polynomial, x_min, x_max, optimal_x, optimal_value);
  325. for (const auto& sample : samples) {
  326. if ((sample.x < x_min) || (sample.x > x_max)) {
  327. continue;
  328. }
  329. const double value = EvaluatePolynomial(polynomial, sample.x);
  330. if (value < *optimal_value) {
  331. *optimal_x = sample.x;
  332. *optimal_value = value;
  333. }
  334. }
  335. }
  336. } // namespace ceres::internal