helloworld_analytic_diff.cc 3.7 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899
  1. // Ceres Solver - A fast non-linear least squares minimizer
  2. // Copyright 2023 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. // A CostFunction implementing analytically derivatives for the
  38. // function f(x) = 10 - x.
  39. class QuadraticCostFunction
  40. : public ceres::SizedCostFunction<1 /* number of residuals */,
  41. 1 /* size of first parameter */> {
  42. public:
  43. bool Evaluate(double const* const* parameters,
  44. double* residuals,
  45. double** jacobians) const override {
  46. double x = parameters[0][0];
  47. // f(x) = 10 - x.
  48. residuals[0] = 10 - x;
  49. // f'(x) = -1. Since there's only 1 parameter and that parameter
  50. // has 1 dimension, there is only 1 element to fill in the
  51. // jacobians.
  52. //
  53. // Since the Evaluate function can be called with the jacobians
  54. // pointer equal to nullptr, the Evaluate function must check to see
  55. // if jacobians need to be computed.
  56. //
  57. // For this simple problem it is overkill to check if jacobians[0]
  58. // is nullptr, but in general when writing more complex
  59. // CostFunctions, it is possible that Ceres may only demand the
  60. // derivatives w.r.t. a subset of the parameter blocks.
  61. if (jacobians != nullptr && jacobians[0] != nullptr) {
  62. jacobians[0][0] = -1;
  63. }
  64. return true;
  65. }
  66. };
  67. int main(int argc, char** argv) {
  68. google::InitGoogleLogging(argv[0]);
  69. // The variable to solve for with its initial value. It will be
  70. // mutated in place by the solver.
  71. double x = 0.5;
  72. const double initial_x = x;
  73. // Build the problem.
  74. ceres::Problem problem;
  75. // Set up the only cost function (also known as residual).
  76. ceres::CostFunction* cost_function = new QuadraticCostFunction;
  77. problem.AddResidualBlock(cost_function, nullptr, &x);
  78. // Run the solver!
  79. ceres::Solver::Options options;
  80. options.minimizer_progress_to_stdout = true;
  81. ceres::Solver::Summary summary;
  82. ceres::Solve(options, &problem, &summary);
  83. std::cout << summary.BriefReport() << "\n";
  84. std::cout << "x : " << initial_x << " -> " << x << "\n";
  85. return 0;
  86. }