Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 1 | .. highlight:: c++ |
| 2 | |
| 3 | .. default-domain:: cpp |
| 4 | |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 5 | .. _chapter-tutorial: |
| 6 | |
| 7 | ======== |
| 8 | Tutorial |
| 9 | ======== |
| 10 | |
Sameer Agarwal | b166806 | 2014-04-29 08:12:19 -0700 | [diff] [blame] | 11 | Ceres solves robustified non-linear bounds constrained least squares |
| 12 | problems 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 Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 18 | |
Sameer Agarwal | c4cd29d | 2014-04-16 23:40:12 -0700 | [diff] [blame] | 19 | Problems of this form comes up in a broad range of areas across |
| 20 | science and engineering - from `fitting curves`_ in statistics, to |
| 21 | constructing `3D models from photographs`_ in computer vision. |
| 22 | |
| 23 | .. _fitting curves: http://en.wikipedia.org/wiki/Nonlinear_regression |
Sameer Agarwal | b166806 | 2014-04-29 08:12:19 -0700 | [diff] [blame] | 24 | .. _3D models from photographs: http://en.wikipedia.org/wiki/Bundle_adjustment |
Sameer Agarwal | c4cd29d | 2014-04-16 23:40:12 -0700 | [diff] [blame] | 25 | |
| 26 | In this chapter we will learn how to solve :eq:`ceresproblem` using |
| 27 | Ceres Solver. Full working code for all the examples described in this |
| 28 | chapter and more can be found in the `examples |
| 29 | <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/>`_ |
| 30 | directory. |
| 31 | |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 32 | The expression |
| 33 | :math:`\rho_i\left(\left\|f_i\left(x_{i_1},...,x_{i_k}\right)\right\|^2\right)` |
| 34 | is 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 |
| 37 | problems small groups of scalars occur together. For example the three |
| 38 | components of a translation vector and the four components of the |
| 39 | quaternion that define the pose of a camera. We refer to such a group |
| 40 | of small scalars as a ``ParameterBlock``. Of course a |
Sameer Agarwal | b166806 | 2014-04-29 08:12:19 -0700 | [diff] [blame] | 41 | ``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 Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 43 | |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 44 | :math:`\rho_i` is a :class:`LossFunction`. A :class:`LossFunction` is |
| 45 | a scalar function that is used to reduce the influence of outliers on |
Sameer Agarwal | b166806 | 2014-04-29 08:12:19 -0700 | [diff] [blame] | 46 | the solution of non-linear least squares problems. |
| 47 | |
| 48 | As a special case, when :math:`\rho_i(x) = x`, i.e., the identity |
| 49 | function, and :math:`l_j = -\infty` and :math:`u_j = \infty` we get |
| 50 | the more familiar `non-linear least squares problem |
Sameer Agarwal | efb47f3 | 2013-02-14 19:44:11 -0800 | [diff] [blame] | 51 | <http://en.wikipedia.org/wiki/Non-linear_least_squares>`_. |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 52 | |
Sameer Agarwal | 7747bb0 | 2013-11-20 11:29:22 -0800 | [diff] [blame] | 53 | .. math:: \frac{1}{2}\sum_{i} \left\|f_i\left(x_{i_1}, ... ,x_{i_k}\right)\right\|^2. |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 54 | :label: ceresproblem2 |
| 55 | |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 56 | .. _section-hello-world: |
| 57 | |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 58 | Hello World! |
| 59 | ============ |
| 60 | |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 61 | To get started, consider the problem of finding the minimum of the |
| 62 | function |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 63 | |
| 64 | .. math:: \frac{1}{2}(10 -x)^2. |
| 65 | |
| 66 | This is a trivial problem, whose minimum is located at :math:`x = 10`, |
| 67 | but it is a good place to start to illustrate the basics of solving a |
| 68 | problem with Ceres [#f1]_. |
| 69 | |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 70 | The first step is to write a functor that will evaluate this the |
| 71 | function :math:`f(x) = 10 - x`: |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 72 | |
| 73 | .. code-block:: c++ |
| 74 | |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 75 | 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 Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 82 | |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 83 | The important thing to note here is that ``operator()`` is a templated |
| 84 | method, which assumes that all its inputs and outputs are of some type |
Sameer Agarwal | b166806 | 2014-04-29 08:12:19 -0700 | [diff] [blame] | 85 | ``T``. The use of templating here allows Ceres to call |
| 86 | ``CostFunctor::operator<T>()``, with ``T=double`` when just the value |
| 87 | of the residual is needed, and with a special type ``T=Jet`` when the |
| 88 | Jacobians are needed. In :ref:`section-derivatives` we will discuss the |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 89 | various ways of supplying derivatives to Ceres in more detail. |
| 90 | |
| 91 | Once we have a way of computing the residual function, it is now time |
| 92 | to construct a non-linear least squares problem using it and have |
| 93 | Ceres 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 Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 124 | } |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 125 | |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 126 | :class:`AutoDiffCostFunction` takes a ``CostFunctor`` as input, |
| 127 | automatically differentiates it and gives it a :class:`CostFunction` |
| 128 | interface. |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 129 | |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 130 | Compiling and running `examples/helloworld.cc |
| 131 | <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/helloworld.cc>`_ |
| 132 | gives us |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 133 | |
| 134 | .. code-block:: bash |
| 135 | |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 136 | 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 Agarwal | 9bd23f2 | 2014-02-24 22:30:21 -0800 | [diff] [blame] | 139 | Ceres Solver Report: Iterations: 2, Initial cost: 1.250000e+01, Final cost: 1.388518e-16, Termination: CONVERGENCE. |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 140 | x : 5 -> 10 |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 141 | |
| 142 | Starting 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 |
| 144 | one linear solve should be enough to get the optimal value. The |
| 145 | default configuration of the solver is aimed at non-linear problems, |
| 146 | and for reasons of simplicity we did not change it in this example. It |
| 147 | is indeed possible to obtain the solution to this problem using Ceres |
| 148 | in one iteration. Also note that the solver did get very close to the |
| 149 | optimal function value of 0 in the very first iteration. We will |
| 150 | discuss these issues in greater detail when we talk about convergence |
| 151 | and parameter settings for Ceres. |
| 152 | |
| 153 | .. rubric:: Footnotes |
| 154 | |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 155 | .. [#f1] `examples/helloworld.cc |
| 156 | <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/helloworld.cc>`_ |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 157 | |
| 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 Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 166 | .. _section-derivatives: |
| 167 | |
| 168 | |
| 169 | Derivatives |
| 170 | =========== |
| 171 | |
| 172 | Ceres Solver like most optimization packages, depends on being able to |
| 173 | evaluate the value and the derivatives of each term in the objective |
| 174 | function at arbitrary parameter values. Doing so correctly and |
| 175 | efficiently is essential to getting good results. Ceres Solver |
| 176 | provides a number of ways of doing so. You have already seen one of |
| 177 | them in action -- |
| 178 | Automatic Differentiation in `examples/helloworld.cc |
| 179 | <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/helloworld.cc>`_ |
| 180 | |
| 181 | We now consider the other two possibilities. Analytic and numeric |
| 182 | derivatives. |
| 183 | |
| 184 | |
| 185 | Numeric Derivatives |
| 186 | ------------------- |
| 187 | |
| 188 | In some cases, its not possible to define a templated cost functor, |
| 189 | for example when the evaluation of the residual involves a call to a |
| 190 | library function that you do not have control over. In such a |
| 191 | situation, numerical differentiation can be used. The user defines a |
| 192 | functor which computes the residual value and construct a |
| 193 | :class:`NumericDiffCostFunction` using it. e.g., for :math:`f(x) = 10 - x` |
| 194 | the 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 | |
| 205 | Which is added to the :class:`Problem` as: |
| 206 | |
| 207 | .. code-block:: c++ |
| 208 | |
| 209 | CostFunction* cost_function = |
Sameer Agarwal | 395b4e9 | 2013-05-23 11:01:16 -0700 | [diff] [blame] | 210 | new NumericDiffCostFunction<NumericDiffCostFunctor, ceres::CENTRAL, 1, 1, 1>( |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 211 | new NumericDiffCostFunctor) |
| 212 | problem.AddResidualBlock(cost_function, NULL, &x); |
| 213 | |
| 214 | Notice 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 Agarwal | eb04dc1 | 2013-05-28 11:01:20 -0700 | [diff] [blame] | 222 | The construction looks almost identical to the one used for automatic |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 223 | differentiation, except for an extra template parameter that indicates |
| 224 | the kind of finite differencing scheme to be used for computing the |
| 225 | numerical derivatives [#f3]_. For more details see the documentation |
| 226 | for :class:`NumericDiffCostFunction`. |
| 227 | |
| 228 | **Generally speaking we recommend automatic differentiation instead of |
| 229 | numeric differentiation. The use of C++ templates makes automatic |
| 230 | differentiation efficient, whereas numeric differentiation is |
| 231 | expensive, prone to numeric errors, and leads to slower convergence.** |
| 232 | |
| 233 | |
| 234 | Analytic Derivatives |
| 235 | -------------------- |
| 236 | |
| 237 | In some cases, using automatic differentiation is not possible. For |
Sameer Agarwal | ebbb984 | 2013-05-26 12:40:12 -0700 | [diff] [blame] | 238 | example, it may be the case that it is more efficient to compute the |
| 239 | derivatives in closed form instead of relying on the chain rule used |
| 240 | by the automatic differentiation code. |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 241 | |
| 242 | In such cases, it is possible to supply your own residual and jacobian |
| 243 | computation code. To do this, define a subclass of |
| 244 | :class:`CostFunction` or :class:`SizedCostFunction` if you know the |
| 245 | sizes of the parameters and residuals at compile time. Here for |
| 246 | example is ``SimpleCostFunction`` that implements :math:`f(x) = 10 - |
| 247 | x`. |
| 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 Agarwal | f0b071b | 2013-05-31 13:22:51 -0700 | [diff] [blame] | 261 | if (jacobians != NULL && jacobians[0] != NULL) { |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 262 | 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 |
| 271 | output array ``jacobians`` for Jacobians. The ``jacobians`` array is |
| 272 | optional, ``Evaluate`` is expected to check when it is non-null, and |
| 273 | if it is the case then fill it with the values of the derivative of |
| 274 | the residual function. In this case since the residual function is |
| 275 | linear, the Jacobian is constant [#f4]_ . |
| 276 | |
| 277 | As can be seen from the above code fragments, implementing |
| 278 | :class:`CostFunction` objects is a bit tedious. We recommend that |
| 279 | unless you have a good reason to manage the jacobian computation |
| 280 | yourself, you use :class:`AutoDiffCostFunction` or |
| 281 | :class:`NumericDiffCostFunction` to construct your residual blocks. |
| 282 | |
Sameer Agarwal | ebbb984 | 2013-05-26 12:40:12 -0700 | [diff] [blame] | 283 | More About Derivatives |
| 284 | ---------------------- |
| 285 | |
| 286 | Computing derivatives is by far the most complicated part of using |
| 287 | Ceres, and depending on the circumstance the user may need more |
| 288 | sophisticated ways of computing derivatives. This section just |
| 289 | scratches the surface of how derivatives can be supplied to |
| 290 | Ceres. Once you are comfortable with using |
| 291 | :class:`NumericDiffCostFunction` and :class:`AutoDiffCostFunction` we |
| 292 | recommend taking a look at :class:`DynamicAutoDiffCostFunction`, |
| 293 | :class:`CostFunctionToFunctor`, :class:`NumericDiffFunctor` and |
| 294 | :class:`ConditionedCostFunction` for more advanced ways of |
| 295 | constructing and computing cost functions. |
| 296 | |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 297 | .. 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 Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 305 | |
| 306 | .. _section-powell: |
| 307 | |
| 308 | Powell's Function |
| 309 | ================= |
| 310 | |
| 311 | Consider now a slightly more complicated example -- the minimization |
| 312 | of Powell's function. Let :math:`x = \left[x_1, x_2, x_3, x_4 \right]` |
| 313 | and |
| 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 Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 322 | F(x) &= \left[f_1(x),\ f_2(x),\ f_3(x),\ f_4(x) \right] |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 323 | \end{align} |
| 324 | |
| 325 | |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 326 | :math:`F(x)` is a function of four parameters, has four residuals |
| 327 | and we wish to find :math:`x` such that :math:`\frac{1}{2}\|F(x)\|^2` |
| 328 | is minimized. |
| 329 | |
| 330 | Again, the first step is to define functors that evaluate of the terms |
| 331 | in the objective functor. Here is the code for evaluating |
| 332 | :math:`f_4(x_1, x_4)`: |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 333 | |
| 334 | .. code-block:: c++ |
| 335 | |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 336 | struct F4 { |
| 337 | template <typename T> |
| 338 | bool operator()(const T* const x1, const T* const x4, T* residual) const { |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 339 | residual[0] = T(sqrt(10.0)) * (x1[0] - x4[0]) * (x1[0] - x4[0]); |
| 340 | return true; |
| 341 | } |
| 342 | }; |
| 343 | |
| 344 | |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 345 | Similarly, 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)` |
| 347 | respectively. Using these, the problem can be constructed as follows: |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 348 | |
| 349 | |
| 350 | .. code-block:: c++ |
| 351 | |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 352 | double x1 = 3.0; double x2 = -1.0; double x3 = 0.0; double x4 = 1.0; |
| 353 | |
| 354 | Problem problem; |
| 355 | |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 356 | // Add residual terms to the problem using the using the autodiff |
| 357 | // wrapper to get the derivatives automatically. |
| 358 | problem.AddResidualBlock( |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 359 | new AutoDiffCostFunction<F1, 1, 1, 1>(new F1), NULL, &x1, &x2); |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 360 | problem.AddResidualBlock( |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 361 | new AutoDiffCostFunction<F2, 1, 1, 1>(new F2), NULL, &x3, &x4); |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 362 | problem.AddResidualBlock( |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 363 | new AutoDiffCostFunction<F3, 1, 1, 1>(new F3), NULL, &x2, &x3) |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 364 | problem.AddResidualBlock( |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 365 | new AutoDiffCostFunction<F4, 1, 1, 1>(new F4), NULL, &x1, &x4); |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 366 | |
| 367 | |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 368 | Note that each ``ResidualBlock`` only depends on the two parameters |
| 369 | that the corresponding residual object depends on and not on all four |
Sameer Agarwal | ebbb984 | 2013-05-26 12:40:12 -0700 | [diff] [blame] | 370 | parameters. Compiling and running `examples/powell.cc |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 371 | <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/powell.cc>`_ |
| 372 | gives us: |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 373 | |
| 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 Agarwal | 9bd23f2 | 2014-02-24 22:30:21 -0800 | [diff] [blame] | 389 | Ceres Solver Report: Iterations: 12, Initial cost: 1.075000e+02, Final cost: 4.584044e-12, Termination: CONVERGENCE. |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 390 | Final x1 = 0.00116741, x2 = -0.000116741, x3 = 0.000190535, x4 = 0.000190535 |
| 391 | |
| 392 | It 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 |
| 395 | function value of :math:`4\times 10^{-12}`. |
| 396 | |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 397 | |
| 398 | .. rubric:: Footnotes |
| 399 | |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 400 | .. [#f5] `examples/powell.cc |
| 401 | <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/powell.cc>`_. |
| 402 | |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 403 | |
| 404 | .. _section-fitting: |
| 405 | |
Sameer Agarwal | 08c891f | 2013-02-04 20:18:58 -0800 | [diff] [blame] | 406 | Curve Fitting |
| 407 | ============= |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 408 | |
| 409 | The examples we have seen until now are simple optimization problems |
| 410 | with no data. The original purpose of least squares and non-linear |
| 411 | least squares analysis was fitting curves to data. It is only |
| 412 | appropriate that we now consider an example of such a problem |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 413 | [#f6]_. It contains data generated by sampling the curve :math:`y = |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 414 | e^{0.3x + 0.1}` and adding Gaussian noise with standard deviation |
Sameer Agarwal | 08c891f | 2013-02-04 20:18:58 -0800 | [diff] [blame] | 415 | :math:`\sigma = 0.2`. Let us fit some data to the curve |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 416 | |
| 417 | .. math:: y = e^{mx + c}. |
| 418 | |
| 419 | We begin by defining a templated object to evaluate the |
| 420 | residual. There will be a residual for each observation. |
| 421 | |
| 422 | .. code-block:: c++ |
| 423 | |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 424 | struct ExponentialResidual { |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 425 | ExponentialResidual(double x, double y) |
| 426 | : x_(x), y_(y) {} |
| 427 | |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 428 | template <typename T> |
| 429 | bool operator()(const T* const m, const T* const c, T* residual) const { |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 430 | 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 Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 440 | Assuming 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 Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 443 | |
| 444 | |
Sameer Agarwal | 08c891f | 2013-02-04 20:18:58 -0800 | [diff] [blame] | 445 | .. code-block:: c++ |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 446 | |
| 447 | double m = 0.0; |
| 448 | double c = 0.0; |
| 449 | |
| 450 | Problem problem; |
| 451 | for (int i = 0; i < kNumObservations; ++i) { |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 452 | 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 Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 456 | } |
| 457 | |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 458 | Compiling and running `examples/curve_fitting.cc |
| 459 | <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/curve_fitting.cc>`_ |
| 460 | gives us: |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 461 | |
| 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 Agarwal | 9bd23f2 | 2014-02-24 22:30:21 -0800 | [diff] [blame] | 478 | Ceres Solver Report: Iterations: 13, Initial cost: 1.211734e+02, Final cost: 1.056751e+00, Termination: CONVERGENCE. |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 479 | Initial m: 0 c: 0 |
| 480 | Final m: 0.291861 c: 0.131439 |
| 481 | |
| 482 | |
| 483 | Starting from parameter values :math:`m = 0, c=0` with an initial |
| 484 | objective 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 |
| 487 | parameters of the original model :math:`m=0.3, c= 0.1`, but this is |
| 488 | expected. When reconstructing a curve from noisy data, we expect to |
| 489 | see such deviations. Indeed, if you were to evaluate the objective |
| 490 | function for :math:`m=0.3, c=0.1`, the fit is worse with an objective |
| 491 | function value of :math:`1.082425`. The figure below illustrates the fit. |
| 492 | |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 493 | .. figure:: least_squares_fit.png |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 494 | :figwidth: 500px |
| 495 | :height: 400px |
| 496 | :align: center |
| 497 | |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 498 | Least squares curve fitting. |
| 499 | |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 500 | |
| 501 | .. rubric:: Footnotes |
| 502 | |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 503 | .. [#f6] `examples/curve_fitting.cc |
| 504 | <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/curve_fitting.cc>`_ |
| 505 | |
| 506 | |
| 507 | Robust Curve Fitting |
| 508 | ===================== |
| 509 | |
| 510 | Now suppose the data we are given has some outliers, i.e., we have |
| 511 | some points that do not obey the noise model. If we were to use the |
| 512 | code above to fit such data, we would get a fit that looks as |
| 513 | below. 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 | |
| 520 | To deal with outliers, a standard technique is to use a |
| 521 | :class:`LossFunction`. Loss functions, reduce the influence of |
| 522 | residual blocks with high residuals, usually the ones corresponding to |
| 523 | outliers. 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 | |
| 529 | to |
| 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 |
| 536 | Solver. The argument :math:`0.5` specifies the scale of the loss |
| 537 | function. As a result, we get the fit below [#f7]_. Notice how the |
| 538 | fitted 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 Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 553 | |
| 554 | |
| 555 | Bundle Adjustment |
| 556 | ================= |
| 557 | |
| 558 | One of the main reasons for writing Ceres was our need to solve large |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 559 | scale bundle adjustment problems [HartleyZisserman]_, [Triggs]_. |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 560 | |
| 561 | Given a set of measured image feature locations and correspondences, |
| 562 | the goal of bundle adjustment is to find 3D point positions and camera |
| 563 | parameters that minimize the reprojection error. This optimization |
| 564 | problem is usually formulated as a non-linear least squares problem, |
| 565 | where the error is the squared :math:`L_2` norm of the difference between |
| 566 | the observed feature location and the projection of the corresponding |
| 567 | 3D point on the image plane of the camera. Ceres has extensive support |
| 568 | for solving bundle adjustment problems. |
| 569 | |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 570 | Let us solve a problem from the `BAL |
| 571 | <http://grail.cs.washington.edu/projects/bal/>`_ dataset [#f8]_. |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 572 | |
| 573 | The first step as usual is to define a templated functor that computes |
| 574 | the reprojection error/residual. The structure of the functor is |
| 575 | similar to the ``ExponentialResidual``, in that there is an |
| 576 | instance of this object responsible for each image observation. |
| 577 | |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 578 | Each residual in a BAL problem depends on a three dimensional point |
| 579 | and a nine parameter camera. The nine parameters defining the camera |
Pablo Speciale | dbc398d | 2013-05-23 15:16:24 -0700 | [diff] [blame] | 580 | are: three for rotation as a Rodriques' axis-angle vector, three |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 581 | for translation, one for focal length and two for radial distortion. |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 582 | The 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 Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 585 | |
| 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 Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 591 | |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 592 | 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 Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 624 | |
| 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 Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 633 | double observed_x; |
| 634 | double observed_y; |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 635 | }; |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 636 | |
| 637 | |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 638 | Note that unlike the examples before, this is a non-trivial function |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 639 | and computing its analytic Jacobian is a bit of a pain. Automatic |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 640 | differentiation makes life much simpler. The function |
| 641 | :func:`AngleAxisRotatePoint` and other functions for manipulating |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 642 | rotations can be found in ``include/ceres/rotation.h``. |
| 643 | |
| 644 | Given this functor, the bundle adjustment problem can be constructed |
| 645 | as follows: |
| 646 | |
| 647 | .. code-block:: c++ |
| 648 | |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 649 | ceres::Problem problem; |
| 650 | for (int i = 0; i < bal_problem.num_observations(); ++i) { |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 651 | ceres::CostFunction* cost_function = |
Sameer Agarwal | b166806 | 2014-04-29 08:12:19 -0700 | [diff] [blame] | 652 | SnavelyReprojectionError::Create( |
| 653 | bal_problem.observations()[2 * i + 0], |
| 654 | bal_problem.observations()[2 * i + 1]); |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 655 | 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 Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 662 | Notice that the problem construction for bundle adjustment is very |
| 663 | similar to the curve fitting example -- one term is added to the |
| 664 | objective function per observation. |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 665 | |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 666 | Since this large sparse problem (well large for ``DENSE_QR`` anyways), |
| 667 | one 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 |
| 670 | a reasonable thing to do, bundle adjustment problems have a special |
| 671 | sparsity structure that can be exploited to solve them much more |
| 672 | efficiently. Ceres provides three specialized solvers (collectively |
| 673 | known as Schur-based solvers) for this task. The example code uses the |
| 674 | simplest of them ``DENSE_SCHUR``. |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 675 | |
| 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 Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 685 | For a more sophisticated bundle adjustment example which demonstrates |
| 686 | the use of Ceres' more advanced features including its various linear |
| 687 | solvers, robust loss functions and local parameterizations see |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 688 | `examples/bundle_adjuster.cc |
| 689 | <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/bundle_adjuster.cc>`_ |
| 690 | |
Sameer Agarwal | 3d87b72 | 2013-02-02 00:49:31 -0800 | [diff] [blame] | 691 | |
| 692 | .. rubric:: Footnotes |
| 693 | |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 694 | .. [#f8] `examples/simple_bundle_adjuster.cc |
| 695 | <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/simple_bundle_adjuster.cc>`_ |
Sameer Agarwal | d91b671 | 2013-02-06 01:08:40 -0800 | [diff] [blame] | 696 | |
| 697 | |
| 698 | Other Examples |
| 699 | ============== |
| 700 | |
| 701 | Besides the examples in this chapter, the `example |
| 702 | <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/>`_ |
| 703 | directory contains a number of other examples: |
| 704 | |
Sameer Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 705 | #. `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 Agarwal | d91b671 | 2013-02-06 01:08:40 -0800 | [diff] [blame] | 710 | #. `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 Agarwal | d91b671 | 2013-02-06 01:08:40 -0800 | [diff] [blame] | 714 | #. `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 Agarwal | 085cd4a | 2013-02-06 14:31:07 -0800 | [diff] [blame] | 720 | #. `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 Agarwal | d91b671 | 2013-02-06 01:08:40 -0800 | [diff] [blame] | 725 | |
Sameer Agarwal | 4437639 | 2013-06-03 09:20:49 -0700 | [diff] [blame] | 726 | #. `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. |