blob: 3842e4ddca5ffde785c11453bb6406f018161bfd [file] [log] [blame]
Austin Schuh70cc9552019-01-21 19:46:48 -08001.. _chapter-solving_faqs:
2
3.. default-domain:: cpp
4
5.. cpp:namespace:: ceres
6
7=======
8Solving
9=======
10
11#. How do I evaluate the Jacobian for a solved problem?
12
13 Using :func:`Problem::Evaluate`.
14
15#. How do I choose the right linear solver?
16
17 When using the ``TRUST_REGION`` minimizer, the choice of linear
18 solver is an important decision. It affects solution quality and
19 runtime. Here is a simple way to reason about it.
20
21 1. For small (a few hundred parameters) or dense problems use
22 ``DENSE_QR``.
23
24 2. For general sparse problems (i.e., the Jacobian matrix has a
25 substantial number of zeros) use
Austin Schuh3de38b02024-06-25 18:25:10 -070026 ``SPARSE_NORMAL_CHOLESKY``.
Austin Schuh70cc9552019-01-21 19:46:48 -080027
28 3. For bundle adjustment problems with up to a hundred or so
29 cameras, use ``DENSE_SCHUR``.
30
31 4. For larger bundle adjustment problems with sparse Schur
Austin Schuh3de38b02024-06-25 18:25:10 -070032 Complement/Reduced camera matrices use ``SPARSE_SCHUR``.
Austin Schuh70cc9552019-01-21 19:46:48 -080033
34 If you do not have access to these libraries for whatever
35 reason, ``ITERATIVE_SCHUR`` with ``SCHUR_JACOBI`` is an
36 excellent alternative.
37
38 5. For large bundle adjustment problems (a few thousand cameras or
39 more) use the ``ITERATIVE_SCHUR`` solver. There are a number of
40 preconditioner choices here. ``SCHUR_JACOBI`` offers an
41 excellent balance of speed and accuracy. This is also the
42 recommended option if you are solving medium sized problems for
43 which ``DENSE_SCHUR`` is too slow but ``SuiteSparse`` is not
44 available.
45
46 .. NOTE::
47
48 If you are solving small to medium sized problems, consider
49 setting ``Solver::Options::use_explicit_schur_complement`` to
50 ``true``, it can result in a substantial performance boost.
51
52 If you are not satisfied with ``SCHUR_JACOBI``'s performance try
53 ``CLUSTER_JACOBI`` and ``CLUSTER_TRIDIAGONAL`` in that
54 order. They require that you have ``SuiteSparse``
55 installed. Both of these preconditioners use a clustering
56 algorithm. Use ``SINGLE_LINKAGE`` before ``CANONICAL_VIEWS``.
57
58#. Use :func:`Solver::Summary::FullReport` to diagnose performance problems.
59
60 When diagnosing Ceres performance issues - runtime and convergence,
61 the first place to start is by looking at the output of
62 ``Solver::Summary::FullReport``. Here is an example
63
64 .. code-block:: bash
65
66 ./bin/bundle_adjuster --input ../data/problem-16-22106-pre.txt
67
68 iter cost cost_change |gradient| |step| tr_ratio tr_radius ls_iter iter_time total_time
69 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
70 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
71 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
72 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
73 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
74 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
75
76 Ceres Solver v1.12.0 Solve Report
77 ----------------------------------
78 Original Reduced
79 Parameter blocks 22122 22122
80 Parameters 66462 66462
81 Residual blocks 83718 83718
82 Residual 167436 167436
83
84 Minimizer TRUST_REGION
85
86 Sparse linear algebra library SUITE_SPARSE
87 Trust region strategy LEVENBERG_MARQUARDT
88
89 Given Used
90 Linear solver SPARSE_SCHUR SPARSE_SCHUR
91 Threads 1 1
92 Linear solver threads 1 1
93 Linear solver ordering AUTOMATIC 22106, 16
94
95 Cost:
96 Initial 4.185660e+06
97 Final 1.803391e+04
98 Change 4.167626e+06
99
100 Minimizer iterations 5
101 Successful steps 5
102 Unsuccessful steps 0
103
104 Time (in seconds):
105 Preprocessor 0.283
106
107 Residual evaluation 0.061
108 Jacobian evaluation 0.361
109 Linear solver 0.382
110 Minimizer 0.895
111
112 Postprocessor 0.002
113 Total 1.220
114
115 Termination: NO_CONVERGENCE (Maximum number of iterations reached.)
116
117 Let us focus on run-time performance. The relevant lines to look at
118 are
119
120
121 .. code-block:: bash
122
123 Time (in seconds):
124 Preprocessor 0.283
125
126 Residual evaluation 0.061
127 Jacobian evaluation 0.361
128 Linear solver 0.382
129 Minimizer 0.895
130
131 Postprocessor 0.002
132 Total 1.220
133
134
135 Which tell us that of the total 1.2 seconds, about .3 seconds was
136 spent in the linear solver and the rest was mostly spent in
137 preprocessing and jacobian evaluation.
138
139 The preprocessing seems particularly expensive. Looking back at the
140 report, we observe
141
142 .. code-block:: bash
143
144 Linear solver ordering AUTOMATIC 22106, 16
145
146 Which indicates that we are using automatic ordering for the
147 ``SPARSE_SCHUR`` solver. This can be expensive at times. A straight
148 forward way to deal with this is to give the ordering manually. For
149 ``bundle_adjuster`` this can be done by passing the flag
150 ``-ordering=user``. Doing so and looking at the timing block of the
151 full report gives us
152
153 .. code-block:: bash
154
155 Time (in seconds):
156 Preprocessor 0.051
157
158 Residual evaluation 0.053
159 Jacobian evaluation 0.344
160 Linear solver 0.372
161 Minimizer 0.854
162
163 Postprocessor 0.002
164 Total 0.935
165
166
167
168 The preprocessor time has gone down by more than 5.5x!.