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