1 // Ceres Solver - A fast non-linear least squares minimizer
2 // Copyright 2010, 2011, 2012 Google Inc. All rights reserved.
3 // http://code.google.com/p/ceres-solver/
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 // An example program that minimizes Powell's singular function.
32 //
33 // F = 1/2 (f1^2 + f2^2 + f3^2 + f4^2)
34 //
35 // f1 = x1 + 10*x2;
36 // f2 = sqrt(5) * (x3 - x4)
37 // f3 = (x2 - 2*x3)^2
38 // f4 = sqrt(10) * (x1 - x4)^2
39 //
40 // The starting values are x1 = 3, x2 = -1, x3 = 0, x4 = 1.
41 // The minimum is 0 at (x1, x2, x3, x4) = 0.
42 //
43 // From: Testing Unconstrained Optimization Software by Jorge J. More, Burton S.
44 // Garbow and Kenneth E. Hillstrom in ACM Transactions on Mathematical Software,
45 // Vol 7(1), March 1981.
46
47 #include <vector>
48 #include "ceres/ceres.h"
49 #include "gflags/gflags.h"
50 #include "glog/logging.h"
51
52 using ceres::AutoDiffCostFunction;
53 using ceres::CostFunction;
54 using ceres::Problem;
55 using ceres::Solver;
56 using ceres::Solve;
57
58 class F1 {
59 public:
operator ()(const T * const x1,const T * const x2,T * residual) const60 template <typename T> bool operator()(const T* const x1,
61 const T* const x2,
62 T* residual) const {
63 // f1 = x1 + 10 * x2;
64 residual[0] = x1[0] + T(10.0) * x2[0];
65 return true;
66 }
67 };
68
69 class F2 {
70 public:
operator ()(const T * const x3,const T * const x4,T * residual) const71 template <typename T> bool operator()(const T* const x3,
72 const T* const x4,
73 T* residual) const {
74 // f2 = sqrt(5) (x3 - x4)
75 residual[0] = T(sqrt(5.0)) * (x3[0] - x4[0]);
76 return true;
77 }
78 };
79
80 class F3 {
81 public:
operator ()(const T * const x2,const T * const x4,T * residual) const82 template <typename T> bool operator()(const T* const x2,
83 const T* const x4,
84 T* residual) const {
85 // f3 = (x2 - 2 x3)^2
86 residual[0] = (x2[0] - T(2.0) * x4[0]) * (x2[0] - T(2.0) * x4[0]);
87 return true;
88 }
89 };
90
91 class F4 {
92 public:
operator ()(const T * const x1,const T * const x4,T * residual) const93 template <typename T> bool operator()(const T* const x1,
94 const T* const x4,
95 T* residual) const {
96 // f4 = sqrt(10) (x1 - x4)^2
97 residual[0] = T(sqrt(10.0)) * (x1[0] - x4[0]) * (x1[0] - x4[0]);
98 return true;
99 }
100 };
101
main(int argc,char ** argv)102 int main(int argc, char** argv) {
103 google::ParseCommandLineFlags(&argc, &argv, true);
104 google::InitGoogleLogging(argv[0]);
105
106 double x1 = 3.0;
107 double x2 = -1.0;
108 double x3 = 0.0;
109 double x4 = 1.0;
110
111 Problem problem;
112 // Add residual terms to the problem using the using the autodiff
113 // wrapper to get the derivatives automatically. The parameters, x1 through
114 // x4, are modified in place.
115 problem.AddResidualBlock(new AutoDiffCostFunction<F1, 1, 1, 1>(new F1),
116 NULL,
117 &x1, &x2);
118 problem.AddResidualBlock(new AutoDiffCostFunction<F2, 1, 1, 1>(new F2),
119 NULL,
120 &x3, &x4);
121 problem.AddResidualBlock(new AutoDiffCostFunction<F3, 1, 1, 1>(new F3),
122 NULL,
123 &x2, &x3);
124 problem.AddResidualBlock(new AutoDiffCostFunction<F4, 1, 1, 1>(new F4),
125 NULL,
126 &x1, &x4);
127
128 // Run the solver!
129 Solver::Options options;
130 options.max_num_iterations = 30;
131 options.linear_solver_type = ceres::DENSE_QR;
132 options.minimizer_progress_to_stdout = true;
133
134 Solver::Summary summary;
135
136 std::cout << "Initial x1 = " << x1
137 << ", x2 = " << x2
138 << ", x3 = " << x3
139 << ", x4 = " << x4
140 << "\n";
141
142 Solve(options, &problem, &summary);
143
144 std::cout << summary.BriefReport() << "\n";
145 std::cout << "Final x1 = " << x1
146 << ", x2 = " << x2
147 << ", x3 = " << x3
148 << ", x4 = " << x4
149 << "\n";
150 return 0;
151 }
152