1 .. highlight:: c++ 2 3 .. default-domain:: cpp 4 5 .. _chapter-tutorial: 6 7 ======== 8 Tutorial 9 ======== 10 Ceres solves robustified non-linear least squares problems of the form 11 12 .. math:: \frac{1}{2}\sum_{i=1} \rho_i\left(\left\|f_i\left(x_{i_1}, ... ,x_{i_k}\right)\right\|^2\right). 13 :label: ceresproblem 14 15 The expression 16 :math:`\rho_i\left(\left\|f_i\left(x_{i_1},...,x_{i_k}\right)\right\|^2\right)` 17 is known as a ``ResidualBlock``, where :math:`f_i(\cdot)` is a 18 :class:`CostFunction` that depends on the parameter blocks 19 :math:`\left[x_{i_1},... , x_{i_k}\right]`. In most optimization 20 problems small groups of scalars occur together. For example the three 21 components of a translation vector and the four components of the 22 quaternion that define the pose of a camera. We refer to such a group 23 of small scalars as a ``ParameterBlock``. Of course a 24 ``ParameterBlock`` can just be a single parameter. 25 26 :math:`\rho_i` is a :class:`LossFunction`. A :class:`LossFunction` is 27 a scalar function that is used to reduce the influence of outliers on 28 the solution of non-linear least squares problems. As a special case, 29 when :math:`\rho_i(x) = x`, i.e., the identity function, we get the 30 more familiar `non-linear least squares problem 31 <http://en.wikipedia.org/wiki/Non-linear_least_squares>`_. 32 33 .. math:: \frac{1}{2}\sum_{i=1} \left\|f_i\left(x_{i_1}, ... ,x_{i_k}\right)\right\|^2. 34 :label: ceresproblem2 35 36 In this chapter we will learn how to solve :eq:`ceresproblem` using 37 Ceres Solver. Full working code for all the examples described in this 38 chapter and more can be found in the `examples 39 <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/>`_ 40 directory. 41 42 .. _section-hello-world: 43 44 Hello World! 45 ============ 46 47 To get started, consider the problem of finding the minimum of the 48 function 49 50 .. math:: \frac{1}{2}(10 -x)^2. 51 52 This is a trivial problem, whose minimum is located at :math:`x = 10`, 53 but it is a good place to start to illustrate the basics of solving a 54 problem with Ceres [#f1]_. 55 56 The first step is to write a functor that will evaluate this the 57 function :math:`f(x) = 10 - x`: 58 59 .. code-block:: c++ 60 61 struct CostFunctor { 62 template <typename T> 63 bool operator()(const T* const x, T* residual) const { 64 residual[0] = T(10.0) - x[0]; 65 return true; 66 } 67 }; 68 69 The important thing to note here is that ``operator()`` is a templated 70 method, which assumes that all its inputs and outputs are of some type 71 ``T``. The reason for using templates here is because Ceres will call 72 ``CostFunctor::operator<T>()``, with ``T=double`` when just the 73 residual is needed, and with a special type ``T=Jet`` when the 74 Jacobians are needed. In :ref:`section-derivatives` we discuss the 75 various ways of supplying derivatives to Ceres in more detail. 76 77 Once we have a way of computing the residual function, it is now time 78 to construct a non-linear least squares problem using it and have 79 Ceres solve it. 80 81 .. code-block:: c++ 82 83 int main(int argc, char** argv) { 84 google::InitGoogleLogging(argv[0]); 85 86 // The variable to solve for with its initial value. 87 double initial_x = 5.0; 88 double x = initial_x; 89 90 // Build the problem. 91 Problem problem; 92 93 // Set up the only cost function (also known as residual). This uses 94 // auto-differentiation to obtain the derivative (jacobian). 95 CostFunction* cost_function = 96 new AutoDiffCostFunction<CostFunctor, 1, 1>(new CostFunctor); 97 problem.AddResidualBlock(cost_function, NULL, &x); 98 99 // Run the solver! 100 Solver::Options options; 101 options.linear_solver_type = ceres::DENSE_QR; 102 options.minimizer_progress_to_stdout = true; 103 Solver::Summary summary; 104 Solve(options, &problem, &summary); 105 106 std::cout << summary.BriefReport() << "\n"; 107 std::cout << "x : " << initial_x 108 << " -> " << x << "\n"; 109 return 0; 110 } 111 112 :class:`AutoDiffCostFunction` takes a ``CostFunctor`` as input, 113 automatically differentiates it and gives it a :class:`CostFunction` 114 interface. 115 116 Compiling and running `examples/helloworld.cc 117 <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/helloworld.cc>`_ 118 gives us 119 120 .. code-block:: bash 121 122 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 123 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 124 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 125 Ceres Solver Report: Iterations: 2, Initial cost: 1.250000e+01, Final cost: 1.388518e-16, Termination: PARAMETER_TOLERANCE. 126 x : 5 -> 10 127 128 Starting from a :math:`x=5`, the solver in two iterations goes to 10 129 [#f2]_. The careful reader will note that this is a linear problem and 130 one linear solve should be enough to get the optimal value. The 131 default configuration of the solver is aimed at non-linear problems, 132 and for reasons of simplicity we did not change it in this example. It 133 is indeed possible to obtain the solution to this problem using Ceres 134 in one iteration. Also note that the solver did get very close to the 135 optimal function value of 0 in the very first iteration. We will 136 discuss these issues in greater detail when we talk about convergence 137 and parameter settings for Ceres. 138 139 .. rubric:: Footnotes 140 141 .. [#f1] `examples/helloworld.cc 142 <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/helloworld.cc>`_ 143 144 .. [#f2] Actually the solver ran for three iterations, and it was 145 by looking at the value returned by the linear solver in the third 146 iteration, it observed that the update to the parameter block was too 147 small and declared convergence. Ceres only prints out the display at 148 the end of an iteration, and terminates as soon as it detects 149 convergence, which is why you only see two iterations here and not 150 three. 151 152 .. _section-derivatives: 153 154 155 Derivatives 156 =========== 157 158 Ceres Solver like most optimization packages, depends on being able to 159 evaluate the value and the derivatives of each term in the objective 160 function at arbitrary parameter values. Doing so correctly and 161 efficiently is essential to getting good results. Ceres Solver 162 provides a number of ways of doing so. You have already seen one of 163 them in action -- 164 Automatic Differentiation in `examples/helloworld.cc 165 <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/helloworld.cc>`_ 166 167 We now consider the other two possibilities. Analytic and numeric 168 derivatives. 169 170 171 Numeric Derivatives 172 ------------------- 173 174 In some cases, its not possible to define a templated cost functor, 175 for example when the evaluation of the residual involves a call to a 176 library function that you do not have control over. In such a 177 situation, numerical differentiation can be used. The user defines a 178 functor which computes the residual value and construct a 179 :class:`NumericDiffCostFunction` using it. e.g., for :math:`f(x) = 10 - x` 180 the corresponding functor would be 181 182 .. code-block:: c++ 183 184 struct NumericDiffCostFunctor { 185 bool operator()(const double* const x, double* residual) const { 186 residual[0] = 10.0 - x[0]; 187 return true; 188 } 189 }; 190 191 Which is added to the :class:`Problem` as: 192 193 .. code-block:: c++ 194 195 CostFunction* cost_function = 196 new NumericDiffCostFunction<NumericDiffCostFunctor, ceres::CENTRAL, 1, 1, 1>( 197 new NumericDiffCostFunctor) 198 problem.AddResidualBlock(cost_function, NULL, &x); 199 200 Notice the parallel from when we were using automatic differentiation 201 202 .. code-block:: c++ 203 204 CostFunction* cost_function = 205 new AutoDiffCostFunction<CostFunctor, 1, 1>(new CostFunctor); 206 problem.AddResidualBlock(cost_function, NULL, &x); 207 208 The construction looks almost identical to the one used for automatic 209 differentiation, except for an extra template parameter that indicates 210 the kind of finite differencing scheme to be used for computing the 211 numerical derivatives [#f3]_. For more details see the documentation 212 for :class:`NumericDiffCostFunction`. 213 214 **Generally speaking we recommend automatic differentiation instead of 215 numeric differentiation. The use of C++ templates makes automatic 216 differentiation efficient, whereas numeric differentiation is 217 expensive, prone to numeric errors, and leads to slower convergence.** 218 219 220 Analytic Derivatives 221 -------------------- 222 223 In some cases, using automatic differentiation is not possible. For 224 example, it may be the case that it is more efficient to compute the 225 derivatives in closed form instead of relying on the chain rule used 226 by the automatic differentiation code. 227 228 In such cases, it is possible to supply your own residual and jacobian 229 computation code. To do this, define a subclass of 230 :class:`CostFunction` or :class:`SizedCostFunction` if you know the 231 sizes of the parameters and residuals at compile time. Here for 232 example is ``SimpleCostFunction`` that implements :math:`f(x) = 10 - 233 x`. 234 235 .. code-block:: c++ 236 237 class QuadraticCostFunction : public ceres::SizedCostFunction<1, 1> { 238 public: 239 virtual ~QuadraticCostFunction() {} 240 virtual bool Evaluate(double const* const* parameters, 241 double* residuals, 242 double** jacobians) const { 243 const double x = parameters[0][0]; 244 residuals[0] = 10 - x; 245 246 // Compute the Jacobian if asked for. 247 if (jacobians != NULL && jacobians[0] != NULL) { 248 jacobians[0][0] = -1; 249 } 250 return true; 251 } 252 }; 253 254 255 ``SimpleCostFunction::Evaluate`` is provided with an input array of 256 ``parameters``, an output array ``residuals`` for residuals and an 257 output array ``jacobians`` for Jacobians. The ``jacobians`` array is 258 optional, ``Evaluate`` is expected to check when it is non-null, and 259 if it is the case then fill it with the values of the derivative of 260 the residual function. In this case since the residual function is 261 linear, the Jacobian is constant [#f4]_ . 262 263 As can be seen from the above code fragments, implementing 264 :class:`CostFunction` objects is a bit tedious. We recommend that 265 unless you have a good reason to manage the jacobian computation 266 yourself, you use :class:`AutoDiffCostFunction` or 267 :class:`NumericDiffCostFunction` to construct your residual blocks. 268 269 More About Derivatives 270 ---------------------- 271 272 Computing derivatives is by far the most complicated part of using 273 Ceres, and depending on the circumstance the user may need more 274 sophisticated ways of computing derivatives. This section just 275 scratches the surface of how derivatives can be supplied to 276 Ceres. Once you are comfortable with using 277 :class:`NumericDiffCostFunction` and :class:`AutoDiffCostFunction` we 278 recommend taking a look at :class:`DynamicAutoDiffCostFunction`, 279 :class:`CostFunctionToFunctor`, :class:`NumericDiffFunctor` and 280 :class:`ConditionedCostFunction` for more advanced ways of 281 constructing and computing cost functions. 282 283 .. rubric:: Footnotes 284 285 .. [#f3] `examples/helloworld_numeric_diff.cc 286 <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/helloworld_numeric_diff.cc>`_. 287 288 .. [#f4] `examples/helloworld_analytic_diff.cc 289 <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/helloworld_analytic_diff.cc>`_. 290 291 292 .. _section-powell: 293 294 Powell's Function 295 ================= 296 297 Consider now a slightly more complicated example -- the minimization 298 of Powell's function. Let :math:`x = \left[x_1, x_2, x_3, x_4 \right]` 299 and 300 301 .. math:: 302 303 \begin{align} 304 f_1(x) &= x_1 + 10x_2 \\ 305 f_2(x) &= \sqrt{5} (x_3 - x_4)\\ 306 f_3(x) &= (x_2 - 2x_3)^2\\ 307 f_4(x) &= \sqrt{10} (x_1 - x_4)^2\\ 308 F(x) &= \left[f_1(x),\ f_2(x),\ f_3(x),\ f_4(x) \right] 309 \end{align} 310 311 312 :math:`F(x)` is a function of four parameters, has four residuals 313 and we wish to find :math:`x` such that :math:`\frac{1}{2}\|F(x)\|^2` 314 is minimized. 315 316 Again, the first step is to define functors that evaluate of the terms 317 in the objective functor. Here is the code for evaluating 318 :math:`f_4(x_1, x_4)`: 319 320 .. code-block:: c++ 321 322 struct F4 { 323 template <typename T> 324 bool operator()(const T* const x1, const T* const x4, T* residual) const { 325 residual[0] = T(sqrt(10.0)) * (x1[0] - x4[0]) * (x1[0] - x4[0]); 326 return true; 327 } 328 }; 329 330 331 Similarly, we can define classes ``F1``, ``F2`` and ``F4`` to evaluate 332 :math:`f_1(x_1, x_2)`, :math:`f_2(x_3, x_4)` and :math:`f_3(x_2, x_3)` 333 respectively. Using these, the problem can be constructed as follows: 334 335 336 .. code-block:: c++ 337 338 double x1 = 3.0; double x2 = -1.0; double x3 = 0.0; double x4 = 1.0; 339 340 Problem problem; 341 342 // Add residual terms to the problem using the using the autodiff 343 // wrapper to get the derivatives automatically. 344 problem.AddResidualBlock( 345 new AutoDiffCostFunction<F1, 1, 1, 1>(new F1), NULL, &x1, &x2); 346 problem.AddResidualBlock( 347 new AutoDiffCostFunction<F2, 1, 1, 1>(new F2), NULL, &x3, &x4); 348 problem.AddResidualBlock( 349 new AutoDiffCostFunction<F3, 1, 1, 1>(new F3), NULL, &x2, &x3) 350 problem.AddResidualBlock( 351 new AutoDiffCostFunction<F4, 1, 1, 1>(new F4), NULL, &x1, &x4); 352 353 354 Note that each ``ResidualBlock`` only depends on the two parameters 355 that the corresponding residual object depends on and not on all four 356 parameters. Compiling and running `examples/powell.cc 357 <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/powell.cc>`_ 358 gives us: 359 360 .. code-block:: bash 361 362 Initial x1 = 3, x2 = -1, x3 = 0, x4 = 1 363 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 364 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 365 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 366 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 367 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 368 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 369 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 370 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 371 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 372 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 373 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 374 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 375 Ceres Solver Report: Iterations: 12, Initial cost: 1.075000e+02, Final cost: 4.584044e-12, Termination: GRADIENT_TOLERANCE. 376 Final x1 = 0.00116741, x2 = -0.000116741, x3 = 0.000190535, x4 = 0.000190535 377 378 It is easy to see that the optimal solution to this problem is at 379 :math:`x_1=0, x_2=0, x_3=0, x_4=0` with an objective function value of 380 :math:`0`. In 10 iterations, Ceres finds a solution with an objective 381 function value of :math:`4\times 10^{-12}`. 382 383 384 .. rubric:: Footnotes 385 386 .. [#f5] `examples/powell.cc 387 <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/powell.cc>`_. 388 389 390 .. _section-fitting: 391 392 Curve Fitting 393 ============= 394 395 The examples we have seen until now are simple optimization problems 396 with no data. The original purpose of least squares and non-linear 397 least squares analysis was fitting curves to data. It is only 398 appropriate that we now consider an example of such a problem 399 [#f6]_. It contains data generated by sampling the curve :math:`y = 400 e^{0.3x + 0.1}` and adding Gaussian noise with standard deviation 401 :math:`\sigma = 0.2`. Let us fit some data to the curve 402 403 .. math:: y = e^{mx + c}. 404 405 We begin by defining a templated object to evaluate the 406 residual. There will be a residual for each observation. 407 408 .. code-block:: c++ 409 410 struct ExponentialResidual { 411 ExponentialResidual(double x, double y) 412 : x_(x), y_(y) {} 413 414 template <typename T> 415 bool operator()(const T* const m, const T* const c, T* residual) const { 416 residual[0] = T(y_) - exp(m[0] * T(x_) + c[0]); 417 return true; 418 } 419 420 private: 421 // Observations for a sample. 422 const double x_; 423 const double y_; 424 }; 425 426 Assuming the observations are in a :math:`2n` sized array called 427 ``data`` the problem construction is a simple matter of creating a 428 :class:`CostFunction` for every observation. 429 430 431 .. code-block:: c++ 432 433 double m = 0.0; 434 double c = 0.0; 435 436 Problem problem; 437 for (int i = 0; i < kNumObservations; ++i) { 438 CostFunction* cost_function = 439 new AutoDiffCostFunction<ExponentialResidual, 1, 1, 1>( 440 new ExponentialResidual(data[2 * i], data[2 * i + 1])); 441 problem.AddResidualBlock(cost_function, NULL, &m, &c); 442 } 443 444 Compiling and running `examples/curve_fitting.cc 445 <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/curve_fitting.cc>`_ 446 gives us: 447 448 .. code-block:: bash 449 450 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 451 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 452 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 453 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 454 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 455 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 456 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 457 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 458 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 459 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 460 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 461 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 462 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 463 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 464 Ceres Solver Report: Iterations: 13, Initial cost: 1.211734e+02, Final cost: 1.056751e+00, Termination: FUNCTION_TOLERANCE. 465 Initial m: 0 c: 0 466 Final m: 0.291861 c: 0.131439 467 468 469 Starting from parameter values :math:`m = 0, c=0` with an initial 470 objective function value of :math:`121.173` Ceres finds a solution 471 :math:`m= 0.291861, c = 0.131439` with an objective function value of 472 :math:`1.05675`. These values are a a bit different than the 473 parameters of the original model :math:`m=0.3, c= 0.1`, but this is 474 expected. When reconstructing a curve from noisy data, we expect to 475 see such deviations. Indeed, if you were to evaluate the objective 476 function for :math:`m=0.3, c=0.1`, the fit is worse with an objective 477 function value of :math:`1.082425`. The figure below illustrates the fit. 478 479 .. figure:: least_squares_fit.png 480 :figwidth: 500px 481 :height: 400px 482 :align: center 483 484 Least squares curve fitting. 485 486 487 .. rubric:: Footnotes 488 489 .. [#f6] `examples/curve_fitting.cc 490 <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/curve_fitting.cc>`_ 491 492 493 Robust Curve Fitting 494 ===================== 495 496 Now suppose the data we are given has some outliers, i.e., we have 497 some points that do not obey the noise model. If we were to use the 498 code above to fit such data, we would get a fit that looks as 499 below. Notice how the fitted curve deviates from the ground truth. 500 501 .. figure:: non_robust_least_squares_fit.png 502 :figwidth: 500px 503 :height: 400px 504 :align: center 505 506 To deal with outliers, a standard technique is to use a 507 :class:`LossFunction`. Loss functions, reduce the influence of 508 residual blocks with high residuals, usually the ones corresponding to 509 outliers. To associate a loss function in a residual block, we change 510 511 .. code-block:: c++ 512 513 problem.AddResidualBlock(cost_function, NULL , &m, &c); 514 515 to 516 517 .. code-block:: c++ 518 519 problem.AddResidualBlock(cost_function, new CauchyLoss(0.5) , &m, &c); 520 521 :class:`CauchyLoss` is one of the loss functions that ships with Ceres 522 Solver. The argument :math:`0.5` specifies the scale of the loss 523 function. As a result, we get the fit below [#f7]_. Notice how the 524 fitted curve moves back closer to the ground truth curve. 525 526 .. figure:: robust_least_squares_fit.png 527 :figwidth: 500px 528 :height: 400px 529 :align: center 530 531 Using :class:`LossFunction` to reduce the effect of outliers on a 532 least squares fit. 533 534 535 .. rubric:: Footnotes 536 537 .. [#f7] `examples/robust_curve_fitting.cc 538 <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/robust_curve_fitting.cc>`_ 539 540 541 Bundle Adjustment 542 ================= 543 544 One of the main reasons for writing Ceres was our need to solve large 545 scale bundle adjustment problems [HartleyZisserman]_, [Triggs]_. 546 547 Given a set of measured image feature locations and correspondences, 548 the goal of bundle adjustment is to find 3D point positions and camera 549 parameters that minimize the reprojection error. This optimization 550 problem is usually formulated as a non-linear least squares problem, 551 where the error is the squared :math:`L_2` norm of the difference between 552 the observed feature location and the projection of the corresponding 553 3D point on the image plane of the camera. Ceres has extensive support 554 for solving bundle adjustment problems. 555 556 Let us solve a problem from the `BAL 557 <http://grail.cs.washington.edu/projects/bal/>`_ dataset [#f8]_. 558 559 The first step as usual is to define a templated functor that computes 560 the reprojection error/residual. The structure of the functor is 561 similar to the ``ExponentialResidual``, in that there is an 562 instance of this object responsible for each image observation. 563 564 Each residual in a BAL problem depends on a three dimensional point 565 and a nine parameter camera. The nine parameters defining the camera 566 are: three for rotation as a Rodriques' axis-angle vector, three 567 for translation, one for focal length and two for radial distortion. 568 The details of this camera model can be found the `Bundler homepage 569 <http://phototour.cs.washington.edu/bundler/>`_ and the `BAL homepage 570 <http://grail.cs.washington.edu/projects/bal/>`_. 571 572 .. code-block:: c++ 573 574 struct SnavelyReprojectionError { 575 SnavelyReprojectionError(double observed_x, double observed_y) 576 : observed_x(observed_x), observed_y(observed_y) {} 577 578 template <typename T> 579 bool operator()(const T* const camera, 580 const T* const point, 581 T* residuals) const { 582 // camera[0,1,2] are the angle-axis rotation. 583 T p[3]; 584 ceres::AngleAxisRotatePoint(camera, point, p); 585 // camera[3,4,5] are the translation. 586 p[0] += camera[3]; p[1] += camera[4]; p[2] += camera[5]; 587 588 // Compute the center of distortion. The sign change comes from 589 // the camera model that Noah Snavely's Bundler assumes, whereby 590 // the camera coordinate system has a negative z axis. 591 T xp = - p[0] / p[2]; 592 T yp = - p[1] / p[2]; 593 594 // Apply second and fourth order radial distortion. 595 const T& l1 = camera[7]; 596 const T& l2 = camera[8]; 597 T r2 = xp*xp + yp*yp; 598 T distortion = T(1.0) + r2 * (l1 + l2 * r2); 599 600 // Compute final projected point position. 601 const T& focal = camera[6]; 602 T predicted_x = focal * distortion * xp; 603 T predicted_y = focal * distortion * yp; 604 605 // The error is the difference between the predicted and observed position. 606 residuals[0] = predicted_x - T(observed_x); 607 residuals[1] = predicted_y - T(observed_y); 608 return true; 609 } 610 611 // Factory to hide the construction of the CostFunction object from 612 // the client code. 613 static ceres::CostFunction* Create(const double observed_x, 614 const double observed_y) { 615 return (new ceres::AutoDiffCostFunction<SnavelyReprojectionError, 2, 9, 3>( 616 new SnavelyReprojectionError(observed_x, observed_y))); 617 } 618 619 double observed_x; 620 double observed_y; 621 }; 622 623 624 Note that unlike the examples before, this is a non-trivial function 625 and computing its analytic Jacobian is a bit of a pain. Automatic 626 differentiation makes life much simpler. The function 627 :func:`AngleAxisRotatePoint` and other functions for manipulating 628 rotations can be found in ``include/ceres/rotation.h``. 629 630 Given this functor, the bundle adjustment problem can be constructed 631 as follows: 632 633 .. code-block:: c++ 634 635 ceres::Problem problem; 636 for (int i = 0; i < bal_problem.num_observations(); ++i) { 637 ceres::CostFunction* cost_function = 638 new ceres::AutoDiffCostFunction<SnavelyReprojectionError, 2, 9, 3>( 639 new SnavelyReprojectionError( 640 bal_problem.observations()[2 * i + 0], 641 bal_problem.observations()[2 * i + 1])); 642 problem.AddResidualBlock(cost_function, 643 NULL /* squared loss */, 644 bal_problem.mutable_camera_for_observation(i), 645 bal_problem.mutable_point_for_observation(i)); 646 } 647 648 649 Notice that the problem construction for bundle adjustment is very 650 similar to the curve fitting example -- one term is added to the 651 objective function per observation. 652 653 Since this large sparse problem (well large for ``DENSE_QR`` anyways), 654 one way to solve this problem is to set 655 :member:`Solver::Options::linear_solver_type` to 656 ``SPARSE_NORMAL_CHOLESKY`` and call :member:`Solve`. And while this is 657 a reasonable thing to do, bundle adjustment problems have a special 658 sparsity structure that can be exploited to solve them much more 659 efficiently. Ceres provides three specialized solvers (collectively 660 known as Schur-based solvers) for this task. The example code uses the 661 simplest of them ``DENSE_SCHUR``. 662 663 .. code-block:: c++ 664 665 ceres::Solver::Options options; 666 options.linear_solver_type = ceres::DENSE_SCHUR; 667 options.minimizer_progress_to_stdout = true; 668 ceres::Solver::Summary summary; 669 ceres::Solve(options, &problem, &summary); 670 std::cout << summary.FullReport() << "\n"; 671 672 For a more sophisticated bundle adjustment example which demonstrates 673 the use of Ceres' more advanced features including its various linear 674 solvers, robust loss functions and local parameterizations see 675 `examples/bundle_adjuster.cc 676 <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/bundle_adjuster.cc>`_ 677 678 679 .. rubric:: Footnotes 680 681 .. [#f8] `examples/simple_bundle_adjuster.cc 682 <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/simple_bundle_adjuster.cc>`_ 683 684 685 Other Examples 686 ============== 687 688 Besides the examples in this chapter, the `example 689 <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/>`_ 690 directory contains a number of other examples: 691 692 #. `bundle_adjuster.cc 693 <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/bundle_adjuster.cc>`_ 694 shows how to use the various features of Ceres to solve bundle 695 adjustment problems. 696 697 #. `circle_fit.cc 698 <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/circle_fit.cc>`_ 699 shows how to fit data to a circle. 700 701 #. `denoising.cc 702 <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/denoising.cc>`_ 703 implements image denoising using the `Fields of Experts 704 <http://www.gris.informatik.tu-darmstadt.de/~sroth/research/foe/index.html>`_ 705 model. 706 707 #. `nist.cc 708 <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/nist.cc>`_ 709 implements and attempts to solves the `NIST 710 <http://www.itl.nist.gov/div898/strd/nls/nls_main.shtm>`_ 711 non-linear regression problems. 712 713 #. `libmv_bundle_adjuster.cc 714 <https://ceres-solver.googlesource.com/ceres-solver/+/master/examples/libmv_bundle_adjuster.cc>`_ 715 is the bundle adjustment algorithm used by `Blender <www.blender.org>`_/libmv. 716 717 718