blob: a8c2d8ef5fb17fc7266d27c22d1c389b1225d9df [file] [log] [blame]
Sameer Agarwal085cd4a2013-02-06 14:31:07 -08001.. highlight:: c++
2
3.. default-domain:: cpp
4
Sameer Agarwal3d87b722013-02-02 00:49:31 -08005.. _chapter-tutorial:
6
7========
8Tutorial
9========
10
Sameer Agarwalb1668062014-04-29 08:12:19 -070011Ceres solves robustified non-linear bounds constrained least squares
12problems of the form
13
14.. math:: :label: ceresproblem
15
16 \min_{\mathbf{x}} &\quad \frac{1}{2}\sum_{i} \rho_i\left(\left\|f_i\left(x_{i_1}, ... ,x_{i_k}\right)\right\|^2\right) \\
17 \text{s.t.} &\quad l_j \le x_j \le u_j
Sameer Agarwal3d87b722013-02-02 00:49:31 -080018
Sameer Agarwalc4cd29d2014-04-16 23:40:12 -070019Problems of this form comes up in a broad range of areas across
20science and engineering - from `fitting curves`_ in statistics, to
21constructing `3D models from photographs`_ in computer vision.
22
23.. _fitting curves: http://en.wikipedia.org/wiki/Nonlinear_regression
Sameer Agarwalb1668062014-04-29 08:12:19 -070024.. _3D models from photographs: http://en.wikipedia.org/wiki/Bundle_adjustment
Sameer Agarwalc4cd29d2014-04-16 23:40:12 -070025
26In this chapter we will learn how to solve :eq:`ceresproblem` using
27Ceres Solver. Full working code for all the examples described in this
28chapter and more can be found in the `examples
29<https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/>`_
30directory.
31
Sameer Agarwal085cd4a2013-02-06 14:31:07 -080032The expression
33:math:`\rho_i\left(\left\|f_i\left(x_{i_1},...,x_{i_k}\right)\right\|^2\right)`
34is known as a ``ResidualBlock``, where :math:`f_i(\cdot)` is a
35:class:`CostFunction` that depends on the parameter blocks
36:math:`\left[x_{i_1},... , x_{i_k}\right]`. In most optimization
37problems small groups of scalars occur together. For example the three
38components of a translation vector and the four components of the
39quaternion that define the pose of a camera. We refer to such a group
40of small scalars as a ``ParameterBlock``. Of course a
Sameer Agarwalb1668062014-04-29 08:12:19 -070041``ParameterBlock`` can just be a single parameter. :math:`l_j` and
42:math:`u_j` are bounds on the parameter block :math:`x_j`.
Sameer Agarwal3d87b722013-02-02 00:49:31 -080043
Sameer Agarwal085cd4a2013-02-06 14:31:07 -080044:math:`\rho_i` is a :class:`LossFunction`. A :class:`LossFunction` is
45a scalar function that is used to reduce the influence of outliers on
Sameer Agarwalb1668062014-04-29 08:12:19 -070046the solution of non-linear least squares problems.
47
48As a special case, when :math:`\rho_i(x) = x`, i.e., the identity
49function, and :math:`l_j = -\infty` and :math:`u_j = \infty` we get
50the more familiar `non-linear least squares problem
Sameer Agarwalefb47f32013-02-14 19:44:11 -080051<http://en.wikipedia.org/wiki/Non-linear_least_squares>`_.
Sameer Agarwal085cd4a2013-02-06 14:31:07 -080052
Sameer Agarwal7747bb02013-11-20 11:29:22 -080053.. math:: \frac{1}{2}\sum_{i} \left\|f_i\left(x_{i_1}, ... ,x_{i_k}\right)\right\|^2.
Sameer Agarwal085cd4a2013-02-06 14:31:07 -080054 :label: ceresproblem2
55
Sameer Agarwal085cd4a2013-02-06 14:31:07 -080056.. _section-hello-world:
57
Sameer Agarwal3d87b722013-02-02 00:49:31 -080058Hello World!
59============
60
Sameer Agarwal085cd4a2013-02-06 14:31:07 -080061To get started, consider the problem of finding the minimum of the
62function
Sameer Agarwal3d87b722013-02-02 00:49:31 -080063
64.. math:: \frac{1}{2}(10 -x)^2.
65
66This is a trivial problem, whose minimum is located at :math:`x = 10`,
67but it is a good place to start to illustrate the basics of solving a
68problem with Ceres [#f1]_.
69
Sameer Agarwal085cd4a2013-02-06 14:31:07 -080070The first step is to write a functor that will evaluate this the
71function :math:`f(x) = 10 - x`:
Sameer Agarwal3d87b722013-02-02 00:49:31 -080072
73.. code-block:: c++
74
Sameer Agarwal085cd4a2013-02-06 14:31:07 -080075 struct CostFunctor {
76 template <typename T>
77 bool operator()(const T* const x, T* residual) const {
78 residual[0] = T(10.0) - x[0];
79 return true;
80 }
81 };
Sameer Agarwal3d87b722013-02-02 00:49:31 -080082
Sameer Agarwal085cd4a2013-02-06 14:31:07 -080083The important thing to note here is that ``operator()`` is a templated
84method, which assumes that all its inputs and outputs are of some type
Sameer Agarwalb1668062014-04-29 08:12:19 -070085``T``. The use of templating here allows Ceres to call
86``CostFunctor::operator<T>()``, with ``T=double`` when just the value
87of the residual is needed, and with a special type ``T=Jet`` when the
88Jacobians are needed. In :ref:`section-derivatives` we will discuss the
Sameer Agarwal085cd4a2013-02-06 14:31:07 -080089various ways of supplying derivatives to Ceres in more detail.
90
91Once we have a way of computing the residual function, it is now time
92to construct a non-linear least squares problem using it and have
93Ceres solve it.
94
95.. code-block:: c++
96
97 int main(int argc, char** argv) {
98 google::InitGoogleLogging(argv[0]);
99
100 // The variable to solve for with its initial value.
101 double initial_x = 5.0;
102 double x = initial_x;
103
104 // Build the problem.
105 Problem problem;
106
107 // Set up the only cost function (also known as residual). This uses
108 // auto-differentiation to obtain the derivative (jacobian).
109 CostFunction* cost_function =
110 new AutoDiffCostFunction<CostFunctor, 1, 1>(new CostFunctor);
111 problem.AddResidualBlock(cost_function, NULL, &x);
112
113 // Run the solver!
114 Solver::Options options;
115 options.linear_solver_type = ceres::DENSE_QR;
116 options.minimizer_progress_to_stdout = true;
117 Solver::Summary summary;
118 Solve(options, &problem, &summary);
119
120 std::cout << summary.BriefReport() << "\n";
121 std::cout << "x : " << initial_x
122 << " -> " << x << "\n";
123 return 0;
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800124 }
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800125
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800126:class:`AutoDiffCostFunction` takes a ``CostFunctor`` as input,
127automatically differentiates it and gives it a :class:`CostFunction`
128interface.
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800129
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800130Compiling and running `examples/helloworld.cc
131<https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/helloworld.cc>`_
132gives us
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800133
134.. code-block:: bash
135
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800136 0: f: 1.250000e+01 d: 0.00e+00 g: 5.00e+00 h: 0.00e+00 rho: 0.00e+00 mu: 1.00e+04 li: 0 it: 6.91e-06 tt: 1.91e-03
137 1: f: 1.249750e-07 d: 1.25e+01 g: 5.00e-04 h: 5.00e+00 rho: 1.00e+00 mu: 3.00e+04 li: 1 it: 2.81e-05 tt: 1.99e-03
138 2: f: 1.388518e-16 d: 1.25e-07 g: 1.67e-08 h: 5.00e-04 rho: 1.00e+00 mu: 9.00e+04 li: 1 it: 1.00e-05 tt: 2.01e-03
Sameer Agarwal9bd23f22014-02-24 22:30:21 -0800139 Ceres Solver Report: Iterations: 2, Initial cost: 1.250000e+01, Final cost: 1.388518e-16, Termination: CONVERGENCE.
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800140 x : 5 -> 10
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800141
142Starting from a :math:`x=5`, the solver in two iterations goes to 10
143[#f2]_. The careful reader will note that this is a linear problem and
144one linear solve should be enough to get the optimal value. The
145default configuration of the solver is aimed at non-linear problems,
146and for reasons of simplicity we did not change it in this example. It
147is indeed possible to obtain the solution to this problem using Ceres
148in one iteration. Also note that the solver did get very close to the
149optimal function value of 0 in the very first iteration. We will
150discuss these issues in greater detail when we talk about convergence
151and parameter settings for Ceres.
152
153.. rubric:: Footnotes
154
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800155.. [#f1] `examples/helloworld.cc
156 <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/helloworld.cc>`_
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800157
158.. [#f2] Actually the solver ran for three iterations, and it was
159 by looking at the value returned by the linear solver in the third
160 iteration, it observed that the update to the parameter block was too
161 small and declared convergence. Ceres only prints out the display at
162 the end of an iteration, and terminates as soon as it detects
163 convergence, which is why you only see two iterations here and not
164 three.
165
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800166.. _section-derivatives:
167
168
169Derivatives
170===========
171
172Ceres Solver like most optimization packages, depends on being able to
173evaluate the value and the derivatives of each term in the objective
174function at arbitrary parameter values. Doing so correctly and
175efficiently is essential to getting good results. Ceres Solver
176provides a number of ways of doing so. You have already seen one of
177them in action --
178Automatic Differentiation in `examples/helloworld.cc
179<https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/helloworld.cc>`_
180
181We now consider the other two possibilities. Analytic and numeric
182derivatives.
183
184
185Numeric Derivatives
186-------------------
187
188In some cases, its not possible to define a templated cost functor,
189for example when the evaluation of the residual involves a call to a
190library function that you do not have control over. In such a
191situation, numerical differentiation can be used. The user defines a
192functor which computes the residual value and construct a
193:class:`NumericDiffCostFunction` using it. e.g., for :math:`f(x) = 10 - x`
194the corresponding functor would be
195
196.. code-block:: c++
197
198 struct NumericDiffCostFunctor {
199 bool operator()(const double* const x, double* residual) const {
200 residual[0] = 10.0 - x[0];
201 return true;
202 }
203 };
204
205Which is added to the :class:`Problem` as:
206
207.. code-block:: c++
208
209 CostFunction* cost_function =
Sameer Agarwal395b4e92013-05-23 11:01:16 -0700210 new NumericDiffCostFunction<NumericDiffCostFunctor, ceres::CENTRAL, 1, 1, 1>(
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800211 new NumericDiffCostFunctor)
212 problem.AddResidualBlock(cost_function, NULL, &x);
213
214Notice the parallel from when we were using automatic differentiation
215
216.. code-block:: c++
217
218 CostFunction* cost_function =
219 new AutoDiffCostFunction<CostFunctor, 1, 1>(new CostFunctor);
220 problem.AddResidualBlock(cost_function, NULL, &x);
221
Sameer Agarwaleb04dc12013-05-28 11:01:20 -0700222The construction looks almost identical to the one used for automatic
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800223differentiation, except for an extra template parameter that indicates
224the kind of finite differencing scheme to be used for computing the
225numerical derivatives [#f3]_. For more details see the documentation
226for :class:`NumericDiffCostFunction`.
227
228**Generally speaking we recommend automatic differentiation instead of
229numeric differentiation. The use of C++ templates makes automatic
230differentiation efficient, whereas numeric differentiation is
231expensive, prone to numeric errors, and leads to slower convergence.**
232
233
234Analytic Derivatives
235--------------------
236
237In some cases, using automatic differentiation is not possible. For
Sameer Agarwalebbb9842013-05-26 12:40:12 -0700238example, it may be the case that it is more efficient to compute the
239derivatives in closed form instead of relying on the chain rule used
240by the automatic differentiation code.
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800241
242In such cases, it is possible to supply your own residual and jacobian
243computation code. To do this, define a subclass of
244:class:`CostFunction` or :class:`SizedCostFunction` if you know the
245sizes of the parameters and residuals at compile time. Here for
246example is ``SimpleCostFunction`` that implements :math:`f(x) = 10 -
247x`.
248
249.. code-block:: c++
250
251 class QuadraticCostFunction : public ceres::SizedCostFunction<1, 1> {
252 public:
253 virtual ~QuadraticCostFunction() {}
254 virtual bool Evaluate(double const* const* parameters,
255 double* residuals,
256 double** jacobians) const {
257 const double x = parameters[0][0];
258 residuals[0] = 10 - x;
259
260 // Compute the Jacobian if asked for.
Sameer Agarwalf0b071b2013-05-31 13:22:51 -0700261 if (jacobians != NULL && jacobians[0] != NULL) {
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800262 jacobians[0][0] = -1;
263 }
264 return true;
265 }
266 };
267
268
269``SimpleCostFunction::Evaluate`` is provided with an input array of
270``parameters``, an output array ``residuals`` for residuals and an
271output array ``jacobians`` for Jacobians. The ``jacobians`` array is
272optional, ``Evaluate`` is expected to check when it is non-null, and
273if it is the case then fill it with the values of the derivative of
274the residual function. In this case since the residual function is
275linear, the Jacobian is constant [#f4]_ .
276
277As can be seen from the above code fragments, implementing
278:class:`CostFunction` objects is a bit tedious. We recommend that
279unless you have a good reason to manage the jacobian computation
280yourself, you use :class:`AutoDiffCostFunction` or
281:class:`NumericDiffCostFunction` to construct your residual blocks.
282
Sameer Agarwalebbb9842013-05-26 12:40:12 -0700283More About Derivatives
284----------------------
285
286Computing derivatives is by far the most complicated part of using
287Ceres, and depending on the circumstance the user may need more
288sophisticated ways of computing derivatives. This section just
289scratches the surface of how derivatives can be supplied to
290Ceres. Once you are comfortable with using
291:class:`NumericDiffCostFunction` and :class:`AutoDiffCostFunction` we
292recommend taking a look at :class:`DynamicAutoDiffCostFunction`,
293:class:`CostFunctionToFunctor`, :class:`NumericDiffFunctor` and
294:class:`ConditionedCostFunction` for more advanced ways of
295constructing and computing cost functions.
296
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800297.. rubric:: Footnotes
298
299.. [#f3] `examples/helloworld_numeric_diff.cc
300 <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/helloworld_numeric_diff.cc>`_.
301
302.. [#f4] `examples/helloworld_analytic_diff.cc
303 <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/helloworld_analytic_diff.cc>`_.
304
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800305
306.. _section-powell:
307
308Powell's Function
309=================
310
311Consider now a slightly more complicated example -- the minimization
312of Powell's function. Let :math:`x = \left[x_1, x_2, x_3, x_4 \right]`
313and
314
315.. math::
316
317 \begin{align}
318 f_1(x) &= x_1 + 10x_2 \\
319 f_2(x) &= \sqrt{5} (x_3 - x_4)\\
320 f_3(x) &= (x_2 - 2x_3)^2\\
321 f_4(x) &= \sqrt{10} (x_1 - x_4)^2\\
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800322 F(x) &= \left[f_1(x),\ f_2(x),\ f_3(x),\ f_4(x) \right]
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800323 \end{align}
324
325
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800326:math:`F(x)` is a function of four parameters, has four residuals
327and we wish to find :math:`x` such that :math:`\frac{1}{2}\|F(x)\|^2`
328is minimized.
329
330Again, the first step is to define functors that evaluate of the terms
331in the objective functor. Here is the code for evaluating
332:math:`f_4(x_1, x_4)`:
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800333
334.. code-block:: c++
335
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800336 struct F4 {
337 template <typename T>
338 bool operator()(const T* const x1, const T* const x4, T* residual) const {
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800339 residual[0] = T(sqrt(10.0)) * (x1[0] - x4[0]) * (x1[0] - x4[0]);
340 return true;
341 }
342 };
343
344
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800345Similarly, we can define classes ``F1``, ``F2`` and ``F4`` to evaluate
346:math:`f_1(x_1, x_2)`, :math:`f_2(x_3, x_4)` and :math:`f_3(x_2, x_3)`
347respectively. Using these, the problem can be constructed as follows:
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800348
349
350.. code-block:: c++
351
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800352 double x1 = 3.0; double x2 = -1.0; double x3 = 0.0; double x4 = 1.0;
353
354 Problem problem;
355
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800356 // Add residual terms to the problem using the using the autodiff
357 // wrapper to get the derivatives automatically.
358 problem.AddResidualBlock(
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800359 new AutoDiffCostFunction<F1, 1, 1, 1>(new F1), NULL, &x1, &x2);
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800360 problem.AddResidualBlock(
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800361 new AutoDiffCostFunction<F2, 1, 1, 1>(new F2), NULL, &x3, &x4);
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800362 problem.AddResidualBlock(
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800363 new AutoDiffCostFunction<F3, 1, 1, 1>(new F3), NULL, &x2, &x3)
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800364 problem.AddResidualBlock(
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800365 new AutoDiffCostFunction<F4, 1, 1, 1>(new F4), NULL, &x1, &x4);
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800366
367
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800368Note that each ``ResidualBlock`` only depends on the two parameters
369that the corresponding residual object depends on and not on all four
Sameer Agarwalebbb9842013-05-26 12:40:12 -0700370parameters. Compiling and running `examples/powell.cc
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800371<https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/powell.cc>`_
372gives us:
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800373
374.. code-block:: bash
375
376 Initial x1 = 3, x2 = -1, x3 = 0, x4 = 1
377 0: f: 1.075000e+02 d: 0.00e+00 g: 1.55e+02 h: 0.00e+00 rho: 0.00e+00 mu: 1.00e+04 li: 0 it: 0.00e+00 tt: 0.00e+00
378 1: f: 5.036190e+00 d: 1.02e+02 g: 2.00e+01 h: 2.16e+00 rho: 9.53e-01 mu: 3.00e+04 li: 1 it: 0.00e+00 tt: 0.00e+00
379 2: f: 3.148168e-01 d: 4.72e+00 g: 2.50e+00 h: 6.23e-01 rho: 9.37e-01 mu: 9.00e+04 li: 1 it: 0.00e+00 tt: 0.00e+00
380 3: f: 1.967760e-02 d: 2.95e-01 g: 3.13e-01 h: 3.08e-01 rho: 9.37e-01 mu: 2.70e+05 li: 1 it: 0.00e+00 tt: 0.00e+00
381 4: f: 1.229900e-03 d: 1.84e-02 g: 3.91e-02 h: 1.54e-01 rho: 9.37e-01 mu: 8.10e+05 li: 1 it: 0.00e+00 tt: 0.00e+00
382 5: f: 7.687123e-05 d: 1.15e-03 g: 4.89e-03 h: 7.69e-02 rho: 9.37e-01 mu: 2.43e+06 li: 1 it: 0.00e+00 tt: 0.00e+00
383 6: f: 4.804625e-06 d: 7.21e-05 g: 6.11e-04 h: 3.85e-02 rho: 9.37e-01 mu: 7.29e+06 li: 1 it: 0.00e+00 tt: 0.00e+00
384 7: f: 3.003028e-07 d: 4.50e-06 g: 7.64e-05 h: 1.92e-02 rho: 9.37e-01 mu: 2.19e+07 li: 1 it: 0.00e+00 tt: 0.00e+00
385 8: f: 1.877006e-08 d: 2.82e-07 g: 9.54e-06 h: 9.62e-03 rho: 9.37e-01 mu: 6.56e+07 li: 1 it: 0.00e+00 tt: 0.00e+00
386 9: f: 1.173223e-09 d: 1.76e-08 g: 1.19e-06 h: 4.81e-03 rho: 9.37e-01 mu: 1.97e+08 li: 1 it: 0.00e+00 tt: 0.00e+00
387 10: f: 7.333425e-11 d: 1.10e-09 g: 1.49e-07 h: 2.40e-03 rho: 9.37e-01 mu: 5.90e+08 li: 1 it: 0.00e+00 tt: 0.00e+00
388 11: f: 4.584044e-12 d: 6.88e-11 g: 1.86e-08 h: 1.20e-03 rho: 9.37e-01 mu: 1.77e+09 li: 1 it: 0.00e+00 tt: 0.00e+00
Sameer Agarwal9bd23f22014-02-24 22:30:21 -0800389 Ceres Solver Report: Iterations: 12, Initial cost: 1.075000e+02, Final cost: 4.584044e-12, Termination: CONVERGENCE.
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800390 Final x1 = 0.00116741, x2 = -0.000116741, x3 = 0.000190535, x4 = 0.000190535
391
392It is easy to see that the optimal solution to this problem is at
393:math:`x_1=0, x_2=0, x_3=0, x_4=0` with an objective function value of
394:math:`0`. In 10 iterations, Ceres finds a solution with an objective
395function value of :math:`4\times 10^{-12}`.
396
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800397
398.. rubric:: Footnotes
399
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800400.. [#f5] `examples/powell.cc
401 <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/powell.cc>`_.
402
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800403
404.. _section-fitting:
405
Sameer Agarwal08c891f2013-02-04 20:18:58 -0800406Curve Fitting
407=============
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800408
409The examples we have seen until now are simple optimization problems
410with no data. The original purpose of least squares and non-linear
411least squares analysis was fitting curves to data. It is only
412appropriate that we now consider an example of such a problem
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800413[#f6]_. It contains data generated by sampling the curve :math:`y =
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800414e^{0.3x + 0.1}` and adding Gaussian noise with standard deviation
Sameer Agarwal08c891f2013-02-04 20:18:58 -0800415:math:`\sigma = 0.2`. Let us fit some data to the curve
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800416
417.. math:: y = e^{mx + c}.
418
419We begin by defining a templated object to evaluate the
420residual. There will be a residual for each observation.
421
422.. code-block:: c++
423
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800424 struct ExponentialResidual {
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800425 ExponentialResidual(double x, double y)
426 : x_(x), y_(y) {}
427
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800428 template <typename T>
429 bool operator()(const T* const m, const T* const c, T* residual) const {
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800430 residual[0] = T(y_) - exp(m[0] * T(x_) + c[0]);
431 return true;
432 }
433
434 private:
435 // Observations for a sample.
436 const double x_;
437 const double y_;
438 };
439
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800440Assuming the observations are in a :math:`2n` sized array called
441``data`` the problem construction is a simple matter of creating a
442:class:`CostFunction` for every observation.
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800443
444
Sameer Agarwal08c891f2013-02-04 20:18:58 -0800445.. code-block:: c++
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800446
447 double m = 0.0;
448 double c = 0.0;
449
450 Problem problem;
451 for (int i = 0; i < kNumObservations; ++i) {
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800452 CostFunction* cost_function =
453 new AutoDiffCostFunction<ExponentialResidual, 1, 1, 1>(
454 new ExponentialResidual(data[2 * i], data[2 * i + 1]));
455 problem.AddResidualBlock(cost_function, NULL, &m, &c);
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800456 }
457
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800458Compiling and running `examples/curve_fitting.cc
459<https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/curve_fitting.cc>`_
460gives us:
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800461
462.. code-block:: bash
463
464 0: f: 1.211734e+02 d: 0.00e+00 g: 3.61e+02 h: 0.00e+00 rho: 0.00e+00 mu: 1.00e+04 li: 0 it: 0.00e+00 tt: 0.00e+00
465 1: f: 1.211734e+02 d:-2.21e+03 g: 3.61e+02 h: 7.52e-01 rho:-1.87e+01 mu: 5.00e+03 li: 1 it: 0.00e+00 tt: 0.00e+00
466 2: f: 1.211734e+02 d:-2.21e+03 g: 3.61e+02 h: 7.51e-01 rho:-1.86e+01 mu: 1.25e+03 li: 1 it: 0.00e+00 tt: 0.00e+00
467 3: f: 1.211734e+02 d:-2.19e+03 g: 3.61e+02 h: 7.48e-01 rho:-1.85e+01 mu: 1.56e+02 li: 1 it: 0.00e+00 tt: 0.00e+00
468 4: f: 1.211734e+02 d:-2.02e+03 g: 3.61e+02 h: 7.22e-01 rho:-1.70e+01 mu: 9.77e+00 li: 1 it: 0.00e+00 tt: 0.00e+00
469 5: f: 1.211734e+02 d:-7.34e+02 g: 3.61e+02 h: 5.78e-01 rho:-6.32e+00 mu: 3.05e-01 li: 1 it: 0.00e+00 tt: 0.00e+00
470 6: f: 3.306595e+01 d: 8.81e+01 g: 4.10e+02 h: 3.18e-01 rho: 1.37e+00 mu: 9.16e-01 li: 1 it: 0.00e+00 tt: 0.00e+00
471 7: f: 6.426770e+00 d: 2.66e+01 g: 1.81e+02 h: 1.29e-01 rho: 1.10e+00 mu: 2.75e+00 li: 1 it: 0.00e+00 tt: 0.00e+00
472 8: f: 3.344546e+00 d: 3.08e+00 g: 5.51e+01 h: 3.05e-02 rho: 1.03e+00 mu: 8.24e+00 li: 1 it: 0.00e+00 tt: 0.00e+00
473 9: f: 1.987485e+00 d: 1.36e+00 g: 2.33e+01 h: 8.87e-02 rho: 9.94e-01 mu: 2.47e+01 li: 1 it: 0.00e+00 tt: 0.00e+00
474 10: f: 1.211585e+00 d: 7.76e-01 g: 8.22e+00 h: 1.05e-01 rho: 9.89e-01 mu: 7.42e+01 li: 1 it: 0.00e+00 tt: 0.00e+00
475 11: f: 1.063265e+00 d: 1.48e-01 g: 1.44e+00 h: 6.06e-02 rho: 9.97e-01 mu: 2.22e+02 li: 1 it: 0.00e+00 tt: 0.00e+00
476 12: f: 1.056795e+00 d: 6.47e-03 g: 1.18e-01 h: 1.47e-02 rho: 1.00e+00 mu: 6.67e+02 li: 1 it: 0.00e+00 tt: 0.00e+00
477 13: f: 1.056751e+00 d: 4.39e-05 g: 3.79e-03 h: 1.28e-03 rho: 1.00e+00 mu: 2.00e+03 li: 1 it: 0.00e+00 tt: 0.00e+00
Sameer Agarwal9bd23f22014-02-24 22:30:21 -0800478 Ceres Solver Report: Iterations: 13, Initial cost: 1.211734e+02, Final cost: 1.056751e+00, Termination: CONVERGENCE.
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800479 Initial m: 0 c: 0
480 Final m: 0.291861 c: 0.131439
481
482
483Starting from parameter values :math:`m = 0, c=0` with an initial
484objective function value of :math:`121.173` Ceres finds a solution
485:math:`m= 0.291861, c = 0.131439` with an objective function value of
486:math:`1.05675`. These values are a a bit different than the
487parameters of the original model :math:`m=0.3, c= 0.1`, but this is
488expected. When reconstructing a curve from noisy data, we expect to
489see such deviations. Indeed, if you were to evaluate the objective
490function for :math:`m=0.3, c=0.1`, the fit is worse with an objective
491function value of :math:`1.082425`. The figure below illustrates the fit.
492
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800493.. figure:: least_squares_fit.png
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800494 :figwidth: 500px
495 :height: 400px
496 :align: center
497
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800498 Least squares curve fitting.
499
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800500
501.. rubric:: Footnotes
502
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800503.. [#f6] `examples/curve_fitting.cc
504 <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/curve_fitting.cc>`_
505
506
507Robust Curve Fitting
508=====================
509
510Now suppose the data we are given has some outliers, i.e., we have
511some points that do not obey the noise model. If we were to use the
512code above to fit such data, we would get a fit that looks as
513below. Notice how the fitted curve deviates from the ground truth.
514
515.. figure:: non_robust_least_squares_fit.png
516 :figwidth: 500px
517 :height: 400px
518 :align: center
519
520To deal with outliers, a standard technique is to use a
521:class:`LossFunction`. Loss functions, reduce the influence of
522residual blocks with high residuals, usually the ones corresponding to
523outliers. To associate a loss function in a residual block, we change
524
525.. code-block:: c++
526
527 problem.AddResidualBlock(cost_function, NULL , &m, &c);
528
529to
530
531.. code-block:: c++
532
533 problem.AddResidualBlock(cost_function, new CauchyLoss(0.5) , &m, &c);
534
535:class:`CauchyLoss` is one of the loss functions that ships with Ceres
536Solver. The argument :math:`0.5` specifies the scale of the loss
537function. As a result, we get the fit below [#f7]_. Notice how the
538fitted curve moves back closer to the ground truth curve.
539
540.. figure:: robust_least_squares_fit.png
541 :figwidth: 500px
542 :height: 400px
543 :align: center
544
545 Using :class:`LossFunction` to reduce the effect of outliers on a
546 least squares fit.
547
548
549.. rubric:: Footnotes
550
551.. [#f7] `examples/robust_curve_fitting.cc
552 <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/robust_curve_fitting.cc>`_
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800553
554
555Bundle Adjustment
556=================
557
558One of the main reasons for writing Ceres was our need to solve large
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800559scale bundle adjustment problems [HartleyZisserman]_, [Triggs]_.
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800560
561Given a set of measured image feature locations and correspondences,
562the goal of bundle adjustment is to find 3D point positions and camera
563parameters that minimize the reprojection error. This optimization
564problem is usually formulated as a non-linear least squares problem,
565where the error is the squared :math:`L_2` norm of the difference between
566the observed feature location and the projection of the corresponding
5673D point on the image plane of the camera. Ceres has extensive support
568for solving bundle adjustment problems.
569
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800570Let us solve a problem from the `BAL
571<http://grail.cs.washington.edu/projects/bal/>`_ dataset [#f8]_.
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800572
573The first step as usual is to define a templated functor that computes
574the reprojection error/residual. The structure of the functor is
575similar to the ``ExponentialResidual``, in that there is an
576instance of this object responsible for each image observation.
577
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800578Each residual in a BAL problem depends on a three dimensional point
579and a nine parameter camera. The nine parameters defining the camera
Pablo Specialedbc398d2013-05-23 15:16:24 -0700580are: three for rotation as a Rodriques' axis-angle vector, three
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800581for translation, one for focal length and two for radial distortion.
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800582The details of this camera model can be found the `Bundler homepage
583<http://phototour.cs.washington.edu/bundler/>`_ and the `BAL homepage
584<http://grail.cs.washington.edu/projects/bal/>`_.
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800585
586.. code-block:: c++
587
588 struct SnavelyReprojectionError {
589 SnavelyReprojectionError(double observed_x, double observed_y)
590 : observed_x(observed_x), observed_y(observed_y) {}
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800591
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800592 template <typename T>
593 bool operator()(const T* const camera,
594 const T* const point,
595 T* residuals) const {
596 // camera[0,1,2] are the angle-axis rotation.
597 T p[3];
598 ceres::AngleAxisRotatePoint(camera, point, p);
599 // camera[3,4,5] are the translation.
600 p[0] += camera[3]; p[1] += camera[4]; p[2] += camera[5];
601
602 // Compute the center of distortion. The sign change comes from
603 // the camera model that Noah Snavely's Bundler assumes, whereby
604 // the camera coordinate system has a negative z axis.
605 T xp = - p[0] / p[2];
606 T yp = - p[1] / p[2];
607
608 // Apply second and fourth order radial distortion.
609 const T& l1 = camera[7];
610 const T& l2 = camera[8];
611 T r2 = xp*xp + yp*yp;
612 T distortion = T(1.0) + r2 * (l1 + l2 * r2);
613
614 // Compute final projected point position.
615 const T& focal = camera[6];
616 T predicted_x = focal * distortion * xp;
617 T predicted_y = focal * distortion * yp;
618
619 // The error is the difference between the predicted and observed position.
620 residuals[0] = predicted_x - T(observed_x);
621 residuals[1] = predicted_y - T(observed_y);
622 return true;
623 }
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800624
625 // Factory to hide the construction of the CostFunction object from
626 // the client code.
627 static ceres::CostFunction* Create(const double observed_x,
628 const double observed_y) {
629 return (new ceres::AutoDiffCostFunction<SnavelyReprojectionError, 2, 9, 3>(
630 new SnavelyReprojectionError(observed_x, observed_y)));
631 }
632
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800633 double observed_x;
634 double observed_y;
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800635 };
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800636
637
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800638Note that unlike the examples before, this is a non-trivial function
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800639and computing its analytic Jacobian is a bit of a pain. Automatic
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800640differentiation makes life much simpler. The function
641:func:`AngleAxisRotatePoint` and other functions for manipulating
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800642rotations can be found in ``include/ceres/rotation.h``.
643
644Given this functor, the bundle adjustment problem can be constructed
645as follows:
646
647.. code-block:: c++
648
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800649 ceres::Problem problem;
650 for (int i = 0; i < bal_problem.num_observations(); ++i) {
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800651 ceres::CostFunction* cost_function =
Sameer Agarwalb1668062014-04-29 08:12:19 -0700652 SnavelyReprojectionError::Create(
653 bal_problem.observations()[2 * i + 0],
654 bal_problem.observations()[2 * i + 1]);
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800655 problem.AddResidualBlock(cost_function,
656 NULL /* squared loss */,
657 bal_problem.mutable_camera_for_observation(i),
658 bal_problem.mutable_point_for_observation(i));
659 }
660
661
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800662Notice that the problem construction for bundle adjustment is very
663similar to the curve fitting example -- one term is added to the
664objective function per observation.
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800665
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800666Since this large sparse problem (well large for ``DENSE_QR`` anyways),
667one way to solve this problem is to set
668:member:`Solver::Options::linear_solver_type` to
669``SPARSE_NORMAL_CHOLESKY`` and call :member:`Solve`. And while this is
670a reasonable thing to do, bundle adjustment problems have a special
671sparsity structure that can be exploited to solve them much more
672efficiently. Ceres provides three specialized solvers (collectively
673known as Schur-based solvers) for this task. The example code uses the
674simplest of them ``DENSE_SCHUR``.
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800675
676.. code-block:: c++
677
678 ceres::Solver::Options options;
679 options.linear_solver_type = ceres::DENSE_SCHUR;
680 options.minimizer_progress_to_stdout = true;
681 ceres::Solver::Summary summary;
682 ceres::Solve(options, &problem, &summary);
683 std::cout << summary.FullReport() << "\n";
684
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800685For a more sophisticated bundle adjustment example which demonstrates
686the use of Ceres' more advanced features including its various linear
687solvers, robust loss functions and local parameterizations see
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800688`examples/bundle_adjuster.cc
689<https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/bundle_adjuster.cc>`_
690
Sameer Agarwal3d87b722013-02-02 00:49:31 -0800691
692.. rubric:: Footnotes
693
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800694.. [#f8] `examples/simple_bundle_adjuster.cc
695 <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/simple_bundle_adjuster.cc>`_
Sameer Agarwald91b6712013-02-06 01:08:40 -0800696
697
698Other Examples
699==============
700
701Besides the examples in this chapter, the `example
702<https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/>`_
703directory contains a number of other examples:
704
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800705#. `bundle_adjuster.cc
706 <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/bundle_adjuster.cc>`_
707 shows how to use the various features of Ceres to solve bundle
708 adjustment problems.
709
Sameer Agarwald91b6712013-02-06 01:08:40 -0800710#. `circle_fit.cc
711 <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/circle_fit.cc>`_
712 shows how to fit data to a circle.
713
Sameer Agarwald91b6712013-02-06 01:08:40 -0800714#. `denoising.cc
715 <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/denoising.cc>`_
716 implements image denoising using the `Fields of Experts
717 <http://www.gris.informatik.tu-darmstadt.de/~sroth/research/foe/index.html>`_
718 model.
719
Sameer Agarwal085cd4a2013-02-06 14:31:07 -0800720#. `nist.cc
721 <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/nist.cc>`_
722 implements and attempts to solves the `NIST
723 <http://www.itl.nist.gov/div898/strd/nls/nls_main.shtm>`_
724 non-linear regression problems.
Sameer Agarwald91b6712013-02-06 01:08:40 -0800725
Sameer Agarwal44376392013-06-03 09:20:49 -0700726#. `libmv_bundle_adjuster.cc
727 <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/libmv_bundle_adjuster.cc>`_
728 is the bundle adjustment algorithm used by `Blender <www.blender.org>`_/libmv.