GoogleGit

blob: 9dbbc7686339f349a8fb660300ede8dc87e6df2b [file] [log] [blame]
  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: keir@google.com (Keir Mierle)
  30. //
  31. // A simple example of using the Ceres minimizer.
  32. //
  33. // Minimize 0.5 (10 - x)^2 using analytic jacobian matrix.
  34. #include <vector>
  35. #include "ceres/ceres.h"
  36. #include "glog/logging.h"
  37. using ceres::CostFunction;
  38. using ceres::SizedCostFunction;
  39. using ceres::Problem;
  40. using ceres::Solver;
  41. using ceres::Solve;
  42. // A CostFunction implementing analytically derivatives for the
  43. // function f(x) = 10 - x.
  44. class QuadraticCostFunction
  45. : public SizedCostFunction<1 /* number of residuals */,
  46. 1 /* size of first parameter */> {
  47. public:
  48. virtual ~QuadraticCostFunction() {}
  49. virtual bool Evaluate(double const* const* parameters,
  50. double* residuals,
  51. double** jacobians) const {
  52. double x = parameters[0][0];
  53. // f(x) = 10 - x.
  54. residuals[0] = 10 - x;
  55. // f'(x) = -1. Since there's only 1 parameter and that parameter
  56. // has 1 dimension, there is only 1 element to fill in the
  57. // jacobians.
  58. //
  59. // Since the Evaluate function can be called with the jacobians
  60. // pointer equal to NULL, the Evaluate function must check to see
  61. // if jacobians need to be computed.
  62. //
  63. // For this simple problem it is overkill to check if jacobians[0]
  64. // is NULL, but in general when writing more complex
  65. // CostFunctions, it is possible that Ceres may only demand the
  66. // derivatives w.r.t. a subset of the parameter blocks.
  67. if (jacobians != NULL && jacobians[0] != NULL) {
  68. jacobians[0][0] = -1;
  69. }
  70. return true;
  71. }
  72. };
  73. int main(int argc, char** argv) {
  74. google::InitGoogleLogging(argv[0]);
  75. // The variable to solve for with its initial value. It will be
  76. // mutated in place by the solver.
  77. double x = 0.5;
  78. const double initial_x = x;
  79. // Build the problem.
  80. Problem problem;
  81. // Set up the only cost function (also known as residual).
  82. CostFunction* cost_function = new QuadraticCostFunction;
  83. problem.AddResidualBlock(cost_function, NULL, &x);
  84. // Run the solver!
  85. Solver::Options options;
  86. options.minimizer_progress_to_stdout = true;
  87. Solver::Summary summary;
  88. Solve(options, &problem, &summary);
  89. std::cout << summary.BriefReport() << "\n";
  90. std::cout << "x : " << initial_x
  91. << " -> " << x << "\n";
  92. return 0;
  93. }