• Home
  • Line#
  • Scopes#
  • Navigate#
  • Raw
  • Download
1.. _chapter-tricks:
2
3===================
4FAQS, Tips & Tricks
5===================
6
7Answers to frequently asked questions, tricks of the trade and general
8wisdom.
9
10Building
11========
12
13#. Use `google-glog <http://code.google.com/p/google-glog>`_.
14
15   Ceres has extensive support for logging detailed information about
16   memory allocations and time consumed in various parts of the solve,
17   internal error conditions etc. This is done logging using the
18   `google-glog <http://code.google.com/p/google-glog>`_ library. We
19   use it extensively to observe and analyze Ceres's
20   performance. `google-glog <http://code.google.com/p/google-glog>`_
21   allows you to control its behaviour from the command line `flags
22   <http://google-glog.googlecode.com/svn/trunk/doc/glog.html>`_. Starting
23   with ``-logtostdterr`` you can add ``-v=N`` for increasing values
24   of ``N`` to get more and more verbose and detailed information
25   about Ceres internals.
26
27   In an attempt to reduce dependencies, it is tempting to use
28   `miniglog` - a minimal implementation of the ``glog`` interface
29   that ships with Ceres. This is a bad idea. ``miniglog`` was written
30   primarily for building and using Ceres on Android because the
31   current version of `google-glog
32   <http://code.google.com/p/google-glog>`_ does not build using the
33   NDK. It has worse performance than the full fledged glog library
34   and is much harder to control and use.
35
36
37Modeling
38========
39
40#. Use analytical/automatic derivatives.
41
42   This is the single most important piece of advice we can give to
43   you. It is tempting to take the easy way out and use numeric
44   differentiation. This is a bad idea. Numeric differentiation is
45   slow, ill-behaved, hard to get right, and results in poor
46   convergence behaviour.
47
48   Ceres allows the user to define templated functors which will
49   be automatically differentiated. For most situations this is enough
50   and we recommend using this facility. In some cases the derivatives
51   are simple enough or the performance considerations are such that
52   the overhead of automatic differentiation is too much. In such
53   cases, analytic derivatives are recommended.
54
55   The use of numerical derivatives should be a measure of last
56   resort, where it is simply not possible to write a templated
57   implementation of the cost function.
58
59   In many cases it is not possible to do analytic or automatic
60   differentiation of the entire cost function, but it is generally
61   the case that it is possible to decompose the cost function into
62   parts that need to be numerically differentiated and parts that can
63   be automatically or analytically differentiated.
64
65   To this end, Ceres has extensive support for mixing analytic,
66   automatic and numeric differentiation. See
67   :class:`NumericDiffFunctor` and :class:`CostFunctionToFunctor`.
68
69#. Putting `Inverse Function Theorem
70   <http://en.wikipedia.org/wiki/Inverse_function_theorem>`_ to use.
71
72   Every now and then we have to deal with functions which cannot be
73   evaluated analytically. Computing the Jacobian in such cases is
74   tricky. A particularly interesting case is where the inverse of the
75   function is easy to compute analytically. An example of such a
76   function is the Coordinate transformation between the `ECEF
77   <http://en.wikipedia.org/wiki/ECEF>`_ and the `WGS84
78   <http://en.wikipedia.org/wiki/World_Geodetic_System>`_ where the
79   conversion from WGS84 from ECEF is analytic, but the conversion
80   back to ECEF uses an iterative algorithm. So how do you compute the
81   derivative of the ECEF to WGS84 transformation?
82
83   One obvious approach would be to numerically
84   differentiate the conversion function. This is not a good idea. For
85   one, it will be slow, but it will also be numerically quite
86   bad.
87
88   Turns out you can use the `Inverse Function Theorem
89   <http://en.wikipedia.org/wiki/Inverse_function_theorem>`_ in this
90   case to compute the derivatives more or less analytically.
91
92   The key result here is. If :math:`x = f^{-1}(y)`, and :math:`Df(x)`
93   is the invertible Jacobian of :math:`f` at :math:`x`. Then the
94   Jacobian :math:`Df^{-1}(y) = [Df(x)]^{-1}`, i.e., the Jacobian of
95   the :math:`f^{-1}` is the inverse of the Jacobian of :math:`f`.
96
97   Algorithmically this means that given :math:`y`, compute :math:`x =
98   f^{-1}(y)` by whatever means you can. Evaluate the Jacobian of
99   :math:`f` at :math:`x`. If the Jacobian matrix is invertible, then
100   the inverse is the Jacobian of the inverse at :math:`y`.
101
102   One can put this into practice with the following code fragment.
103
104   .. code-block:: c++
105
106      Eigen::Vector3d ecef; // Fill some values
107      // Iterative computation.
108      Eigen::Vector3d lla = ECEFToLLA(ecef);
109      // Analytic derivatives
110      Eigen::Matrix3d lla_to_ecef_jacobian = LLAToECEFJacobian(lla);
111      bool invertible;
112      Eigen::Matrix3d ecef_to_lla_jacobian;
113      lla_to_ecef_jacobian.computeInverseWithCheck(ecef_to_lla_jacobian, invertible);
114
115#. When using Quaternions, use :class:`QuaternionParameterization`.
116
117   TBD
118
119#. How to choose a parameter block size?
120
121   TBD
122
123Solving
124=======
125
126#. Choosing a linear solver.
127
128   When using the ``TRUST_REGION`` minimizer, the choice of linear
129   solver is an important decision. It affects solution quality and
130   runtime. Here is a simple way to reason about it.
131
132   1. For small (a few hundred parameters) or dense problems use
133      ``DENSE_QR``.
134
135   2. For general sparse problems (i.e., the Jacobian matrix has a
136      substantial number of zeros) use
137      ``SPARSE_NORMAL_CHOLESKY``. This requires that you have
138      ``SuiteSparse`` or ``CXSparse`` installed.
139
140   3. For bundle adjustment problems with up to a hundred or so
141      cameras, use ``DENSE_SCHUR``.
142
143   4. For larger bundle adjustment problems with sparse Schur
144      Complement/Reduced camera matrices use ``SPARSE_SCHUR``. This
145      requires that you have ``SuiteSparse`` or ``CXSparse``
146      installed.
147
148   5. For large bundle adjustment problems (a few thousand cameras or
149      more) use the ``ITERATIVE_SCHUR`` solver. There are a number of
150      preconditioner choices here. ``SCHUR_JACOBI`` offers an
151      excellent balance of speed and accuracy. This is also the
152      recommended option if you are solving medium sized problems for
153      which ``DENSE_SCHUR`` is too slow but ``SuiteSparse`` is not
154      available.
155
156      If you are not satisfied with ``SCHUR_JACOBI``'s performance try
157      ``CLUSTER_JACOBI`` and ``CLUSTER_TRIDIAGONAL`` in that
158      order. They require that you have ``SuiteSparse``
159      installed. Both of these preconditioners use a clustering
160      algorithm. Use ``SINGLE_LINKAGE`` before ``CANONICAL_VIEWS``.
161
162#. Use `Solver::Summary::FullReport` to diagnose performance problems.
163
164   When diagnosing Ceres performance issues - runtime and convergence,
165   the first place to start is by looking at the output of
166   ``Solver::Summary::FullReport``. Here is an example
167
168   .. code-block:: bash
169
170     ./bin/bundle_adjuster --input ../data/problem-16-22106-pre.txt
171
172     iter      cost      cost_change  |gradient|   |step|    tr_ratio  tr_radius  ls_iter  iter_time  total_time
173        0  4.185660e+06    0.00e+00    2.16e+07   0.00e+00   0.00e+00  1.00e+04       0    7.50e-02    3.58e-01
174        1  1.980525e+05    3.99e+06    5.34e+06   2.40e+03   9.60e-01  3.00e+04       1    1.84e-01    5.42e-01
175        2  5.086543e+04    1.47e+05    2.11e+06   1.01e+03   8.22e-01  4.09e+04       1    1.53e-01    6.95e-01
176        3  1.859667e+04    3.23e+04    2.87e+05   2.64e+02   9.85e-01  1.23e+05       1    1.71e-01    8.66e-01
177        4  1.803857e+04    5.58e+02    2.69e+04   8.66e+01   9.93e-01  3.69e+05       1    1.61e-01    1.03e+00
178        5  1.803391e+04    4.66e+00    3.11e+02   1.02e+01   1.00e+00  1.11e+06       1    1.49e-01    1.18e+00
179
180     Ceres Solver v1.10.0 Solve Report
181     ----------------------------------
182                                          Original                  Reduced
183     Parameter blocks                        22122                    22122
184     Parameters                              66462                    66462
185     Residual blocks                         83718                    83718
186     Residual                               167436                   167436
187
188     Minimizer                        TRUST_REGION
189
190     Sparse linear algebra library    SUITE_SPARSE
191     Trust region strategy     LEVENBERG_MARQUARDT
192
193                                             Given                     Used
194     Linear solver                    SPARSE_SCHUR             SPARSE_SCHUR
195     Threads                                     1                        1
196     Linear solver threads                       1                        1
197     Linear solver ordering              AUTOMATIC                22106, 16
198
199     Cost:
200     Initial                          4.185660e+06
201     Final                            1.803391e+04
202     Change                           4.167626e+06
203
204     Minimizer iterations                        5
205     Successful steps                            5
206     Unsuccessful steps                          0
207
208     Time (in seconds):
209     Preprocessor                            0.283
210
211       Residual evaluation                   0.061
212       Jacobian evaluation                   0.361
213       Linear solver                         0.382
214     Minimizer                               0.895
215
216     Postprocessor                           0.002
217     Total                                   1.220
218
219     Termination:                   NO_CONVERGENCE (Maximum number of iterations reached.)
220
221  Let us focus on run-time performance. The relevant lines to look at
222  are
223
224
225   .. code-block:: bash
226
227     Time (in seconds):
228     Preprocessor                            0.283
229
230       Residual evaluation                   0.061
231       Jacobian evaluation                   0.361
232       Linear solver                         0.382
233     Minimizer                               0.895
234
235     Postprocessor                           0.002
236     Total                                   1.220
237
238
239  Which tell us that of the total 1.2 seconds, about .3 seconds was
240  spent in the linear solver and the rest was mostly spent in
241  preprocessing and jacobian evaluation.
242
243  The preprocessing seems particularly expensive. Looking back at the
244  report, we observe
245
246   .. code-block:: bash
247
248     Linear solver ordering              AUTOMATIC                22106, 16
249
250  Which indicates that we are using automatic ordering for the
251  ``SPARSE_SCHUR`` solver. This can be expensive at times. A straight
252  forward way to deal with this is to give the ordering manually. For
253  ``bundle_adjuster`` this can be done by passing the flag
254  ``-ordering=user``. Doing so and looking at the timing block of the
255  full report gives us
256
257   .. code-block:: bash
258
259     Time (in seconds):
260     Preprocessor                            0.051
261
262       Residual evaluation                   0.053
263       Jacobian evaluation                   0.344
264       Linear solver                         0.372
265     Minimizer                               0.854
266
267     Postprocessor                           0.002
268     Total                                   0.935
269
270
271
272  The preprocessor time has gone down by more than 5.5x!.
273
274Further Reading
275===============
276
277For a short but informative introduction to the subject we recommend
278the booklet by [Madsen]_ . For a general introduction to non-linear
279optimization we recommend [NocedalWright]_. [Bjorck]_ remains the
280seminal reference on least squares problems. [TrefethenBau]_ book is
281our favorite text on introductory numerical linear algebra. [Triggs]_
282provides a thorough coverage of the bundle adjustment problem.
283