blob: 0a449cb99f1970ed74f13a55970c9736d1cd9a2d [file] [log] [blame]
Keir Mierle8ebb0732012-04-30 23:09:08 -07001// Ceres Solver - A fast non-linear least squares minimizer
2// Copyright 2010, 2011, 2012 Google Inc. All rights reserved.
3// http://code.google.com/p/ceres-solver/
4//
5// Redistribution and use in source and binary forms, with or without
6// modification, are permitted provided that the following conditions are met:
7//
8// * Redistributions of source code must retain the above copyright notice,
9// this list of conditions and the following disclaimer.
10// * Redistributions in binary form must reproduce the above copyright notice,
11// this list of conditions and the following disclaimer in the documentation
12// and/or other materials provided with the distribution.
13// * Neither the name of Google Inc. nor the names of its contributors may be
14// used to endorse or promote products derived from this software without
15// specific prior written permission.
16//
17// THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
18// AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
19// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
20// ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE
21// LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
22// CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
23// SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
24// INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
25// CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
26// ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
27// POSSIBILITY OF SUCH DAMAGE.
28//
29// Author: sameeragarwal@google.com (Sameer Agarwal)
30// keir@google.com (Keir Mierle)
31//
32// The Problem object is used to build and hold least squares problems.
33
34#ifndef CERES_PUBLIC_PROBLEM_H_
35#define CERES_PUBLIC_PROBLEM_H_
36
37#include <cstddef>
38#include <map>
39#include <set>
40#include <vector>
41
Keir Mierle8ebb0732012-04-30 23:09:08 -070042#include "ceres/internal/macros.h"
43#include "ceres/internal/port.h"
44#include "ceres/internal/scoped_ptr.h"
45#include "ceres/types.h"
Sameer Agarwal509f68c2013-02-20 01:39:03 -080046#include "glog/logging.h"
47
Keir Mierle8ebb0732012-04-30 23:09:08 -070048
49namespace ceres {
50
51class CostFunction;
52class LossFunction;
53class LocalParameterization;
Keir Mierle6196cba2012-06-18 11:03:40 -070054class Solver;
Sameer Agarwal509f68c2013-02-20 01:39:03 -080055struct CRSMatrix;
Keir Mierle8ebb0732012-04-30 23:09:08 -070056
57namespace internal {
58class Preprocessor;
59class ProblemImpl;
60class ParameterBlock;
61class ResidualBlock;
Keir Mierle8ebb0732012-04-30 23:09:08 -070062} // namespace internal
63
Keir Mierle04938ef2013-02-17 12:37:55 -080064// A ResidualBlockId is an opaque handle clients can use to remove residual
65// blocks from a Problem after adding them.
66typedef internal::ResidualBlock* ResidualBlockId;
Keir Mierle8ebb0732012-04-30 23:09:08 -070067
68// A class to represent non-linear least squares problems. Such
69// problems have a cost function that is a sum of error terms (known
70// as "residuals"), where each residual is a function of some subset
71// of the parameters. The cost function takes the form
72//
73// N 1
74// SUM --- loss( || r_i1, r_i2,..., r_ik ||^2 ),
75// i=1 2
76//
77// where
78//
79// r_ij is residual number i, component j; the residual is a
80// function of some subset of the parameters x1...xk. For
81// example, in a structure from motion problem a residual
82// might be the difference between a measured point in an
83// image and the reprojected position for the matching
84// camera, point pair. The residual would have two
85// components, error in x and error in y.
86//
87// loss(y) is the loss function; for example, squared error or
88// Huber L1 loss. If loss(y) = y, then the cost function is
89// non-robustified least squares.
90//
91// This class is specifically designed to address the important subset
92// of "sparse" least squares problems, where each component of the
93// residual depends only on a small number number of parameters, even
94// though the total number of residuals and parameters may be very
95// large. This property affords tremendous gains in scale, allowing
96// efficient solving of large problems that are otherwise
97// inaccessible.
98//
99// The canonical example of a sparse least squares problem is
100// "structure-from-motion" (SFM), where the parameters are points and
101// cameras, and residuals are reprojection errors. Typically a single
102// residual will depend only on 9 parameters (3 for the point, 6 for
103// the camera).
104//
105// To create a least squares problem, use the AddResidualBlock() and
106// AddParameterBlock() methods, documented below. Here is an example least
107// squares problem containing 3 parameter blocks of sizes 3, 4 and 5
108// respectively and two residual terms of size 2 and 6:
109//
110// double x1[] = { 1.0, 2.0, 3.0 };
111// double x2[] = { 1.0, 2.0, 3.0, 5.0 };
112// double x3[] = { 1.0, 2.0, 3.0, 6.0, 7.0 };
113//
114// Problem problem;
115//
116// problem.AddResidualBlock(new MyUnaryCostFunction(...), x1);
117// problem.AddResidualBlock(new MyBinaryCostFunction(...), x2, x3);
118//
119// Please see cost_function.h for details of the CostFunction object.
120class Problem {
121 public:
122 struct Options {
123 Options()
124 : cost_function_ownership(TAKE_OWNERSHIP),
125 loss_function_ownership(TAKE_OWNERSHIP),
Sameer Agarwal8e1f83c2013-02-15 08:35:40 -0800126 local_parameterization_ownership(TAKE_OWNERSHIP),
Keir Mierle04938ef2013-02-17 12:37:55 -0800127 enable_fast_parameter_block_removal(false),
Sameer Agarwal8e1f83c2013-02-15 08:35:40 -0800128 disable_all_safety_checks(false) {}
Keir Mierle8ebb0732012-04-30 23:09:08 -0700129
130 // These flags control whether the Problem object owns the cost
131 // functions, loss functions, and parameterizations passed into
132 // the Problem. If set to TAKE_OWNERSHIP, then the problem object
133 // will delete the corresponding cost or loss functions on
134 // destruction. The destructor is careful to delete the pointers
135 // only once, since sharing cost/loss/parameterizations is
136 // allowed.
137 Ownership cost_function_ownership;
138 Ownership loss_function_ownership;
139 Ownership local_parameterization_ownership;
Sameer Agarwal8e1f83c2013-02-15 08:35:40 -0800140
Keir Mierle04938ef2013-02-17 12:37:55 -0800141 // If true, trades memory for a faster RemoveParameterBlock() operation.
142 //
143 // RemoveParameterBlock() takes time proportional to the size of the entire
144 // Problem. If you only remove parameter blocks from the Problem
145 // occassionaly, this may be acceptable. However, if you are modifying the
146 // Problem frequently, and have memory to spare, then flip this switch to
147 // make RemoveParameterBlock() take time proportional to the number of
148 // residual blocks that depend on it. The increase in memory usage is an
149 // additonal hash set per parameter block containing all the residuals that
150 // depend on the parameter block.
151 bool enable_fast_parameter_block_removal;
152
153 // By default, Ceres performs a variety of safety checks when constructing
154 // the problem. There is a small but measurable performance penalty to
155 // these checks, typically around 5% of construction time. If you are sure
156 // your problem construction is correct, and 5% of the problem construction
Sameer Agarwal8e1f83c2013-02-15 08:35:40 -0800157 // time is truly an overhead you want to avoid, then you can set
158 // disable_all_safety_checks to true.
159 //
Keir Mierle04938ef2013-02-17 12:37:55 -0800160 // WARNING: Do not set this to true, unless you are absolutely sure of what
Sameer Agarwal8e1f83c2013-02-15 08:35:40 -0800161 // you are doing.
162 bool disable_all_safety_checks;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700163 };
164
165 // The default constructor is equivalent to the
166 // invocation Problem(Problem::Options()).
167 Problem();
168 explicit Problem(const Options& options);
169
170 ~Problem();
171
172 // Add a residual block to the overall cost function. The cost
173 // function carries with it information about the sizes of the
174 // parameter blocks it expects. The function checks that these match
175 // the sizes of the parameter blocks listed in parameter_blocks. The
176 // program aborts if a mismatch is detected. loss_function can be
177 // NULL, in which case the cost of the term is just the squared norm
178 // of the residuals.
179 //
180 // The user has the option of explicitly adding the parameter blocks
181 // using AddParameterBlock. This causes additional correctness
182 // checking; however, AddResidualBlock implicitly adds the parameter
183 // blocks if they are not present, so calling AddParameterBlock
184 // explicitly is not required.
185 //
186 // The Problem object by default takes ownership of the
187 // cost_function and loss_function pointers. These objects remain
188 // live for the life of the Problem object. If the user wishes to
189 // keep control over the destruction of these objects, then they can
190 // do this by setting the corresponding enums in the Options struct.
191 //
192 // Note: Even though the Problem takes ownership of cost_function
193 // and loss_function, it does not preclude the user from re-using
194 // them in another residual block. The destructor takes care to call
195 // delete on each cost_function or loss_function pointer only once,
196 // regardless of how many residual blocks refer to them.
197 //
198 // Example usage:
199 //
200 // double x1[] = {1.0, 2.0, 3.0};
201 // double x2[] = {1.0, 2.0, 5.0, 6.0};
202 // double x3[] = {3.0, 6.0, 2.0, 5.0, 1.0};
203 //
204 // Problem problem;
205 //
206 // problem.AddResidualBlock(new MyUnaryCostFunction(...), NULL, x1);
207 // problem.AddResidualBlock(new MyBinaryCostFunction(...), NULL, x2, x1);
208 //
209 ResidualBlockId AddResidualBlock(CostFunction* cost_function,
210 LossFunction* loss_function,
211 const vector<double*>& parameter_blocks);
212
213 // Convenience methods for adding residuals with a small number of
214 // parameters. This is the common case. Instead of specifying the
215 // parameter block arguments as a vector, list them as pointers.
216 ResidualBlockId AddResidualBlock(CostFunction* cost_function,
217 LossFunction* loss_function,
218 double* x0);
219 ResidualBlockId AddResidualBlock(CostFunction* cost_function,
220 LossFunction* loss_function,
221 double* x0, double* x1);
222 ResidualBlockId AddResidualBlock(CostFunction* cost_function,
223 LossFunction* loss_function,
224 double* x0, double* x1, double* x2);
225 ResidualBlockId AddResidualBlock(CostFunction* cost_function,
226 LossFunction* loss_function,
227 double* x0, double* x1, double* x2,
228 double* x3);
229 ResidualBlockId AddResidualBlock(CostFunction* cost_function,
230 LossFunction* loss_function,
231 double* x0, double* x1, double* x2,
232 double* x3, double* x4);
233 ResidualBlockId AddResidualBlock(CostFunction* cost_function,
234 LossFunction* loss_function,
235 double* x0, double* x1, double* x2,
236 double* x3, double* x4, double* x5);
Fisher12626e82012-10-21 14:12:04 -0400237 ResidualBlockId AddResidualBlock(CostFunction* cost_function,
238 LossFunction* loss_function,
239 double* x0, double* x1, double* x2,
240 double* x3, double* x4, double* x5,
241 double* x6);
242 ResidualBlockId AddResidualBlock(CostFunction* cost_function,
243 LossFunction* loss_function,
244 double* x0, double* x1, double* x2,
245 double* x3, double* x4, double* x5,
246 double* x6, double* x7);
247 ResidualBlockId AddResidualBlock(CostFunction* cost_function,
248 LossFunction* loss_function,
249 double* x0, double* x1, double* x2,
250 double* x3, double* x4, double* x5,
251 double* x6, double* x7, double* x8);
252 ResidualBlockId AddResidualBlock(CostFunction* cost_function,
253 LossFunction* loss_function,
254 double* x0, double* x1, double* x2,
255 double* x3, double* x4, double* x5,
256 double* x6, double* x7, double* x8,
257 double* x9);
Keir Mierle8ebb0732012-04-30 23:09:08 -0700258
259 // Add a parameter block with appropriate size to the problem.
260 // Repeated calls with the same arguments are ignored. Repeated
261 // calls with the same double pointer but a different size results
262 // in undefined behaviour.
263 void AddParameterBlock(double* values, int size);
264
265 // Add a parameter block with appropriate size and parameterization
266 // to the problem. Repeated calls with the same arguments are
267 // ignored. Repeated calls with the same double pointer but a
268 // different size results in undefined behaviour.
269 void AddParameterBlock(double* values,
270 int size,
271 LocalParameterization* local_parameterization);
272
Keir Mierle04938ef2013-02-17 12:37:55 -0800273 // Remove a parameter block from the problem. The parameterization of the
274 // parameter block, if it exists, will persist until the deletion of the
275 // problem (similar to cost/loss functions in residual block removal). Any
276 // residual blocks that depend on the parameter are also removed, as
277 // described above in RemoveResidualBlock().
278 //
279 // If Problem::Options::enable_fast_parameter_block_removal is true, then the
280 // removal is fast (almost constant time). Otherwise, removing a parameter
281 // block will incur a scan of the entire Problem object.
282 //
283 // WARNING: Removing a residual or parameter block will destroy the implicit
284 // ordering, rendering the jacobian or residuals returned from the solver
285 // uninterpretable. If you depend on the evaluated jacobian, do not use
286 // remove! This may change in a future release.
287 void RemoveParameterBlock(double* values);
288
289 // Remove a residual block from the problem. Any parameters that the residual
290 // block depends on are not removed. The cost and loss functions for the
291 // residual block will not get deleted immediately; won't happen until the
292 // problem itself is deleted.
293 //
294 // WARNING: Removing a residual or parameter block will destroy the implicit
295 // ordering, rendering the jacobian or residuals returned from the solver
296 // uninterpretable. If you depend on the evaluated jacobian, do not use
297 // remove! This may change in a future release.
298 void RemoveResidualBlock(ResidualBlockId residual_block);
299
Keir Mierle8ebb0732012-04-30 23:09:08 -0700300 // Hold the indicated parameter block constant during optimization.
301 void SetParameterBlockConstant(double* values);
302
303 // Allow the indicated parameter to vary during optimization.
304 void SetParameterBlockVariable(double* values);
305
306 // Set the local parameterization for one of the parameter blocks.
307 // The local_parameterization is owned by the Problem by default. It
308 // is acceptable to set the same parameterization for multiple
309 // parameters; the destructor is careful to delete local
310 // parameterizations only once. The local parameterization can only
311 // be set once per parameter, and cannot be changed once set.
312 void SetParameterization(double* values,
313 LocalParameterization* local_parameterization);
314
315 // Number of parameter blocks in the problem. Always equals
316 // parameter_blocks().size() and parameter_block_sizes().size().
317 int NumParameterBlocks() const;
318
319 // The size of the parameter vector obtained by summing over the
320 // sizes of all the parameter blocks.
321 int NumParameters() const;
322
323 // Number of residual blocks in the problem. Always equals
324 // residual_blocks().size().
325 int NumResidualBlocks() const;
326
327 // The size of the residual vector obtained by summing over the
328 // sizes of all of the residual blocks.
329 int NumResiduals() const;
330
Sameer Agarwal509f68c2013-02-20 01:39:03 -0800331 // Options struct to control Problem::Evaluate.
332 struct EvaluateOptions {
333 EvaluateOptions()
Sameer Agarwal039ff072013-02-26 09:15:39 -0800334 : apply_loss_function(true),
335 num_threads(1) {
Sameer Agarwal509f68c2013-02-20 01:39:03 -0800336 }
337
338 // The set of parameter blocks for which evaluation should be
339 // performed. This vector determines the order that parameter
340 // blocks occur in the gradient vector and in the columns of the
341 // jacobian matrix. If parameter_blocks is empty, then it is
342 // assumed to be equal to vector containing ALL the parameter
343 // blocks. Generally speaking the parameter blocks will occur in
344 // the order in which they were added to the problem. But, this
345 // may change if the user removes any parameter blocks from the
346 // problem.
347 //
348 // NOTE: This vector should contain the same pointers as the ones
Sameer Agarwal931c3092013-02-25 09:46:21 -0800349 // used to add parameter blocks to the Problem. These parameter
Sameer Agarwal509f68c2013-02-20 01:39:03 -0800350 // block should NOT point to new memory locations. Bad things will
351 // happen otherwise.
352 vector<double*> parameter_blocks;
353
354 // The set of residual blocks to evaluate. This vector determines
355 // the order in which the residuals occur, and how the rows of the
356 // jacobian are ordered. If residual_blocks is empty, then it is
357 // assumed to be equal to the vector containing all the residual
358 // blocks. If this vector is empty, then it is assumed to be equal
359 // to a vector containing ALL the residual blocks. Generally
360 // speaking the residual blocks will occur in the order in which
361 // they were added to the problem. But, this may change if the
362 // user removes any residual blocks from the problem.
363 vector<ResidualBlockId> residual_blocks;
Sameer Agarwal039ff072013-02-26 09:15:39 -0800364
365 // Even though the residual blocks in the problem may contain loss
366 // functions, setting apply_loss_function to false will turn off
367 // the application of the loss function to the output of the cost
368 // function. This is of use for example if the user wishes to
369 // analyse the solution quality by studying the distribution of
370 // residuals before and after the solve.
371 bool apply_loss_function;
372
Sameer Agarwal509f68c2013-02-20 01:39:03 -0800373 int num_threads;
374 };
375
376 // Evaluate Problem. Any of the output pointers can be NULL. Which
377 // residual blocks and parameter blocks are used is controlled by
378 // the EvaluateOptions struct above.
379 //
380 // Note 1: The evaluation will use the values stored in the memory
381 // locations pointed to by the parameter block pointers used at the
382 // time of the construction of the problem. i.e.,
383 //
384 // Problem problem;
385 // double x = 1;
Sameer Agarwal5e7ce8a2013-03-06 11:38:41 -0800386 // problem.AddResidualBlock(new MyCostFunction, NULL, &x);
Sameer Agarwal509f68c2013-02-20 01:39:03 -0800387 //
388 // double cost = 0.0;
389 // problem.Evaluate(Problem::EvaluateOptions(), &cost, NULL, NULL, NULL);
390 //
391 // The cost is evaluated at x = 1. If you wish to evaluate the
392 // problem at x = 2, then
393 //
394 // x = 2;
395 // problem.Evaluate(Problem::EvaluateOptions(), &cost, NULL, NULL, NULL);
396 //
397 // is the way to do so.
398 //
399 // Note 2: If no local parameterizations are used, then the size of
400 // the gradient vector (and the number of columns in the jacobian)
401 // is the sum of the sizes of all the parameter blocks. If a
402 // parameter block has a local parameterization, then it contributes
Sameer Agarwal931c3092013-02-25 09:46:21 -0800403 // "LocalSize" entries to the gradient vector (and the number of
Sameer Agarwal509f68c2013-02-20 01:39:03 -0800404 // columns in the jacobian).
405 bool Evaluate(const EvaluateOptions& options,
406 double* cost,
407 vector<double>* residuals,
408 vector<double>* gradient,
409 CRSMatrix* jacobian);
410
Keir Mierle8ebb0732012-04-30 23:09:08 -0700411 private:
Keir Mierle6196cba2012-06-18 11:03:40 -0700412 friend class Solver;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700413 internal::scoped_ptr<internal::ProblemImpl> problem_impl_;
Sameer Agarwal237d6592012-05-30 20:34:49 -0700414 CERES_DISALLOW_COPY_AND_ASSIGN(Problem);
Keir Mierle8ebb0732012-04-30 23:09:08 -0700415};
416
417} // namespace ceres
418
419#endif // CERES_PUBLIC_PROBLEM_H_