sampled_function.cc 3.3 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687
  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: sameeragarwal@google.com (Sameer Agarwal)
  30. //
  31. // A simple example of optimizing a sampled function by using cubic
  32. // interpolation.
  33. #include "ceres/ceres.h"
  34. #include "ceres/cubic_interpolation.h"
  35. #include "glog/logging.h"
  36. using Interpolator = ceres::CubicInterpolator<ceres::Grid1D<double>>;
  37. // A simple cost functor that interfaces an interpolated table of
  38. // values with automatic differentiation.
  39. struct InterpolatedCostFunctor {
  40. explicit InterpolatedCostFunctor(const Interpolator& interpolator)
  41. : interpolator(interpolator) {}
  42. template <typename T>
  43. bool operator()(const T* x, T* residuals) const {
  44. interpolator.Evaluate(*x, residuals);
  45. return true;
  46. }
  47. static ceres::CostFunction* Create(const Interpolator& interpolator) {
  48. return new ceres::AutoDiffCostFunction<InterpolatedCostFunctor, 1, 1>(
  49. new InterpolatedCostFunctor(interpolator));
  50. }
  51. private:
  52. const Interpolator& interpolator;
  53. };
  54. int main(int argc, char** argv) {
  55. google::InitGoogleLogging(argv[0]);
  56. // Evaluate the function f(x) = (x - 4.5)^2;
  57. const int kNumSamples = 10;
  58. double values[kNumSamples];
  59. for (int i = 0; i < kNumSamples; ++i) {
  60. values[i] = (i - 4.5) * (i - 4.5);
  61. }
  62. ceres::Grid1D<double> array(values, 0, kNumSamples);
  63. Interpolator interpolator(array);
  64. double x = 1.0;
  65. ceres::Problem problem;
  66. ceres::CostFunction* cost_function =
  67. InterpolatedCostFunctor::Create(interpolator);
  68. problem.AddResidualBlock(cost_function, nullptr, &x);
  69. ceres::Solver::Options options;
  70. options.minimizer_progress_to_stdout = true;
  71. ceres::Solver::Summary summary;
  72. ceres::Solve(options, &problem, &summary);
  73. std::cout << summary.BriefReport() << "\n";
  74. std::cout << "Expected x: 4.5. Actual x : " << x << std::endl;
  75. return 0;
  76. }