blob: 503e0fe950f01d9e749a390af16033c976fd85b0 [file] [log] [blame]
Austin Schuh70cc9552019-01-21 19:46:48 -08001// Ceres Solver - A fast non-linear least squares minimizer
2// Copyright 2015 Google Inc. All rights reserved.
3// http://ceres-solver.org/
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 <array>
38#include <cstddef>
39#include <map>
40#include <memory>
41#include <set>
42#include <vector>
43
44#include "ceres/context.h"
45#include "ceres/internal/disable_warnings.h"
46#include "ceres/internal/port.h"
47#include "ceres/types.h"
48#include "glog/logging.h"
49
50namespace ceres {
51
52class CostFunction;
53class LossFunction;
54class LocalParameterization;
55class Solver;
56struct CRSMatrix;
57
58namespace internal {
59class Preprocessor;
60class ProblemImpl;
61class ParameterBlock;
62class ResidualBlock;
63} // namespace internal
64
65// A ResidualBlockId is an opaque handle clients can use to remove residual
66// blocks from a Problem after adding them.
67typedef internal::ResidualBlock* ResidualBlockId;
68
69// A class to represent non-linear least squares problems. Such
70// problems have a cost function that is a sum of error terms (known
71// as "residuals"), where each residual is a function of some subset
72// of the parameters. The cost function takes the form
73//
74// N 1
75// SUM --- loss( || r_i1, r_i2,..., r_ik ||^2 ),
76// i=1 2
77//
78// where
79//
80// r_ij is residual number i, component j; the residual is a
81// function of some subset of the parameters x1...xk. For
82// example, in a structure from motion problem a residual
83// might be the difference between a measured point in an
84// image and the reprojected position for the matching
85// camera, point pair. The residual would have two
86// components, error in x and error in y.
87//
88// loss(y) is the loss function; for example, squared error or
89// Huber L1 loss. If loss(y) = y, then the cost function is
90// non-robustified least squares.
91//
92// This class is specifically designed to address the important subset
93// of "sparse" least squares problems, where each component of the
94// residual depends only on a small number number of parameters, even
95// though the total number of residuals and parameters may be very
96// large. This property affords tremendous gains in scale, allowing
97// efficient solving of large problems that are otherwise
98// inaccessible.
99//
100// The canonical example of a sparse least squares problem is
101// "structure-from-motion" (SFM), where the parameters are points and
102// cameras, and residuals are reprojection errors. Typically a single
103// residual will depend only on 9 parameters (3 for the point, 6 for
104// the camera).
105//
106// To create a least squares problem, use the AddResidualBlock() and
107// AddParameterBlock() methods, documented below. Here is an example least
108// squares problem containing 3 parameter blocks of sizes 3, 4 and 5
109// respectively and two residual terms of size 2 and 6:
110//
111// double x1[] = { 1.0, 2.0, 3.0 };
112// double x2[] = { 1.0, 2.0, 3.0, 5.0 };
113// double x3[] = { 1.0, 2.0, 3.0, 6.0, 7.0 };
114//
115// Problem problem;
116//
117// problem.AddResidualBlock(new MyUnaryCostFunction(...), NULL, x1);
118// problem.AddResidualBlock(new MyBinaryCostFunction(...), NULL, x2, x3);
119//
120// Please see cost_function.h for details of the CostFunction object.
121class CERES_EXPORT Problem {
122 public:
123 struct CERES_EXPORT Options {
124 // These flags control whether the Problem object owns the cost
125 // functions, loss functions, and parameterizations passed into
126 // the Problem. If set to TAKE_OWNERSHIP, then the problem object
127 // will delete the corresponding cost or loss functions on
128 // destruction. The destructor is careful to delete the pointers
129 // only once, since sharing cost/loss/parameterizations is
130 // allowed.
131 Ownership cost_function_ownership = TAKE_OWNERSHIP;
132 Ownership loss_function_ownership = TAKE_OWNERSHIP;
133 Ownership local_parameterization_ownership = TAKE_OWNERSHIP;
134
135 // If true, trades memory for faster RemoveResidualBlock() and
136 // RemoveParameterBlock() operations.
137 //
138 // By default, RemoveParameterBlock() and RemoveResidualBlock() take time
139 // proportional to the size of the entire problem. If you only ever remove
140 // parameters or residuals from the problem occasionally, this might be
141 // acceptable. However, if you have memory to spare, enable this option to
142 // make RemoveParameterBlock() take time proportional to the number of
143 // residual blocks that depend on it, and RemoveResidualBlock() take (on
144 // average) constant time.
145 //
146 // The increase in memory usage is twofold: an additional hash set per
147 // parameter block containing all the residuals that depend on the parameter
148 // block; and a hash set in the problem containing all residuals.
149 bool enable_fast_removal = false;
150
151 // By default, Ceres performs a variety of safety checks when constructing
152 // the problem. There is a small but measurable performance penalty to
153 // these checks, typically around 5% of construction time. If you are sure
154 // your problem construction is correct, and 5% of the problem construction
155 // time is truly an overhead you want to avoid, then you can set
156 // disable_all_safety_checks to true.
157 //
158 // WARNING: Do not set this to true, unless you are absolutely sure of what
159 // you are doing.
160 bool disable_all_safety_checks = false;
161
162 // A Ceres global context to use for solving this problem. This may help to
163 // reduce computation time as Ceres can reuse expensive objects to create.
164 // The context object can be NULL, in which case Ceres may create one.
165 //
166 // Ceres does NOT take ownership of the pointer.
167 Context* context = nullptr;
168 };
169
170 // The default constructor is equivalent to the
171 // invocation Problem(Problem::Options()).
172 Problem();
173 explicit Problem(const Options& options);
174 Problem(const Problem&) = delete;
175 void operator=(const Problem&) = delete;
176
177 ~Problem();
178
179 // Add a residual block to the overall cost function. The cost
180 // function carries with its information about the sizes of the
181 // parameter blocks it expects. The function checks that these match
182 // the sizes of the parameter blocks listed in parameter_blocks. The
183 // program aborts if a mismatch is detected. loss_function can be
184 // NULL, in which case the cost of the term is just the squared norm
185 // of the residuals.
186 //
187 // The user has the option of explicitly adding the parameter blocks
188 // using AddParameterBlock. This causes additional correctness
189 // checking; however, AddResidualBlock implicitly adds the parameter
190 // blocks if they are not present, so calling AddParameterBlock
191 // explicitly is not required.
192 //
193 // The Problem object by default takes ownership of the
194 // cost_function and loss_function pointers. These objects remain
195 // live for the life of the Problem object. If the user wishes to
196 // keep control over the destruction of these objects, then they can
197 // do this by setting the corresponding enums in the Options struct.
198 //
199 // Note: Even though the Problem takes ownership of cost_function
200 // and loss_function, it does not preclude the user from re-using
201 // them in another residual block. The destructor takes care to call
202 // delete on each cost_function or loss_function pointer only once,
203 // regardless of how many residual blocks refer to them.
204 //
205 // Example usage:
206 //
207 // double x1[] = {1.0, 2.0, 3.0};
208 // double x2[] = {1.0, 2.0, 5.0, 6.0};
209 // double x3[] = {3.0, 6.0, 2.0, 5.0, 1.0};
210 //
211 // Problem problem;
212 //
213 // problem.AddResidualBlock(new MyUnaryCostFunction(...), NULL, x1);
214 // problem.AddResidualBlock(new MyBinaryCostFunction(...), NULL, x2, x1);
215 //
216 // Add a residual block by listing the parameter block pointers
217 // directly instead of wapping them in a container.
218 template <typename... Ts>
219 ResidualBlockId AddResidualBlock(CostFunction* cost_function,
220 LossFunction* loss_function,
221 double* x0,
222 Ts*... xs) {
223 const std::array<double*, sizeof...(Ts) + 1> parameter_blocks{{x0, xs...}};
224 return AddResidualBlock(cost_function, loss_function,
225 parameter_blocks.data(),
226 static_cast<int>(parameter_blocks.size()));
227 }
228
229 // Add a residual block by providing a vector of parameter blocks.
230 ResidualBlockId AddResidualBlock(
231 CostFunction* cost_function,
232 LossFunction* loss_function,
233 const std::vector<double*>& parameter_blocks);
234
235 // Add a residual block by providing a pointer to the parameter block array
236 // and the number of parameter blocks.
237 ResidualBlockId AddResidualBlock(
238 CostFunction* cost_function,
239 LossFunction* loss_function,
240 double* const* const parameter_blocks,
241 int num_parameter_blocks);
242
243 // Add a parameter block with appropriate size to the problem.
244 // Repeated calls with the same arguments are ignored. Repeated
245 // calls with the same double pointer but a different size results
246 // in undefined behaviour.
247 void AddParameterBlock(double* values, int size);
248
249 // Add a parameter block with appropriate size and parameterization
250 // to the problem. Repeated calls with the same arguments are
251 // ignored. Repeated calls with the same double pointer but a
252 // different size results in undefined behaviour.
253 void AddParameterBlock(double* values,
254 int size,
255 LocalParameterization* local_parameterization);
256
257 // Remove a parameter block from the problem. The parameterization of the
258 // parameter block, if it exists, will persist until the deletion of the
259 // problem (similar to cost/loss functions in residual block removal). Any
260 // residual blocks that depend on the parameter are also removed, as
261 // described above in RemoveResidualBlock().
262 //
263 // If Problem::Options::enable_fast_removal is true, then the
264 // removal is fast (almost constant time). Otherwise, removing a parameter
265 // block will incur a scan of the entire Problem object.
266 //
267 // WARNING: Removing a residual or parameter block will destroy the implicit
268 // ordering, rendering the jacobian or residuals returned from the solver
269 // uninterpretable. If you depend on the evaluated jacobian, do not use
270 // remove! This may change in a future release.
271 void RemoveParameterBlock(double* values);
272
273 // Remove a residual block from the problem. Any parameters that the residual
274 // block depends on are not removed. The cost and loss functions for the
275 // residual block will not get deleted immediately; won't happen until the
276 // problem itself is deleted.
277 //
278 // WARNING: Removing a residual or parameter block will destroy the implicit
279 // ordering, rendering the jacobian or residuals returned from the solver
280 // uninterpretable. If you depend on the evaluated jacobian, do not use
281 // remove! This may change in a future release.
282 void RemoveResidualBlock(ResidualBlockId residual_block);
283
284 // Hold the indicated parameter block constant during optimization.
285 void SetParameterBlockConstant(double* values);
286
287 // Allow the indicated parameter block to vary during optimization.
288 void SetParameterBlockVariable(double* values);
289
290 // Returns true if a parameter block is set constant, and false otherwise.
291 bool IsParameterBlockConstant(double* values) const;
292
293 // Set the local parameterization for one of the parameter blocks.
294 // The local_parameterization is owned by the Problem by default. It
295 // is acceptable to set the same parameterization for multiple
296 // parameters; the destructor is careful to delete local
297 // parameterizations only once. The local parameterization can only
298 // be set once per parameter, and cannot be changed once set.
299 void SetParameterization(double* values,
300 LocalParameterization* local_parameterization);
301
302 // Get the local parameterization object associated with this
303 // parameter block. If there is no parameterization object
304 // associated then NULL is returned.
305 const LocalParameterization* GetParameterization(double* values) const;
306
307 // Set the lower/upper bound for the parameter at position "index".
308 void SetParameterLowerBound(double* values, int index, double lower_bound);
309 void SetParameterUpperBound(double* values, int index, double upper_bound);
310
311 // Get the lower/upper bound for the parameter at position
312 // "index". If the parameter is not bounded by the user, then its
313 // lower bound is -std::numeric_limits<double>::max() and upper
314 // bound is std::numeric_limits<double>::max().
315 double GetParameterLowerBound(double* values, int index) const;
316 double GetParameterUpperBound(double* values, int index) const;
317
318 // Number of parameter blocks in the problem. Always equals
319 // parameter_blocks().size() and parameter_block_sizes().size().
320 int NumParameterBlocks() const;
321
322 // The size of the parameter vector obtained by summing over the
323 // sizes of all the parameter blocks.
324 int NumParameters() const;
325
326 // Number of residual blocks in the problem. Always equals
327 // residual_blocks().size().
328 int NumResidualBlocks() const;
329
330 // The size of the residual vector obtained by summing over the
331 // sizes of all of the residual blocks.
332 int NumResiduals() const;
333
334 // The size of the parameter block.
335 int ParameterBlockSize(const double* values) const;
336
337 // The size of local parameterization for the parameter block. If
338 // there is no local parameterization associated with this parameter
339 // block, then ParameterBlockLocalSize = ParameterBlockSize.
340 int ParameterBlockLocalSize(const double* values) const;
341
342 // Is the given parameter block present in this problem or not?
343 bool HasParameterBlock(const double* values) const;
344
345 // Fills the passed parameter_blocks vector with pointers to the
346 // parameter blocks currently in the problem. After this call,
347 // parameter_block.size() == NumParameterBlocks.
348 void GetParameterBlocks(std::vector<double*>* parameter_blocks) const;
349
350 // Fills the passed residual_blocks vector with pointers to the
351 // residual blocks currently in the problem. After this call,
352 // residual_blocks.size() == NumResidualBlocks.
353 void GetResidualBlocks(std::vector<ResidualBlockId>* residual_blocks) const;
354
355 // Get all the parameter blocks that depend on the given residual block.
356 void GetParameterBlocksForResidualBlock(
357 const ResidualBlockId residual_block,
358 std::vector<double*>* parameter_blocks) const;
359
360 // Get the CostFunction for the given residual block.
361 const CostFunction* GetCostFunctionForResidualBlock(
362 const ResidualBlockId residual_block) const;
363
364 // Get the LossFunction for the given residual block. Returns NULL
365 // if no loss function is associated with this residual block.
366 const LossFunction* GetLossFunctionForResidualBlock(
367 const ResidualBlockId residual_block) const;
368
369 // Get all the residual blocks that depend on the given parameter block.
370 //
371 // If Problem::Options::enable_fast_removal is true, then
372 // getting the residual blocks is fast and depends only on the number of
373 // residual blocks. Otherwise, getting the residual blocks for a parameter
374 // block will incur a scan of the entire Problem object.
375 void GetResidualBlocksForParameterBlock(
376 const double* values,
377 std::vector<ResidualBlockId>* residual_blocks) const;
378
379 // Options struct to control Problem::Evaluate.
380 struct EvaluateOptions {
381 // The set of parameter blocks for which evaluation should be
382 // performed. This vector determines the order that parameter
383 // blocks occur in the gradient vector and in the columns of the
384 // jacobian matrix. If parameter_blocks is empty, then it is
385 // assumed to be equal to vector containing ALL the parameter
386 // blocks. Generally speaking the parameter blocks will occur in
387 // the order in which they were added to the problem. But, this
388 // may change if the user removes any parameter blocks from the
389 // problem.
390 //
391 // NOTE: This vector should contain the same pointers as the ones
392 // used to add parameter blocks to the Problem. These parameter
393 // block should NOT point to new memory locations. Bad things will
394 // happen otherwise.
395 std::vector<double*> parameter_blocks;
396
397 // The set of residual blocks to evaluate. This vector determines
398 // the order in which the residuals occur, and how the rows of the
399 // jacobian are ordered. If residual_blocks is empty, then it is
400 // assumed to be equal to the vector containing ALL the residual
401 // blocks. Generally speaking the residual blocks will occur in
402 // the order in which they were added to the problem. But, this
403 // may change if the user removes any residual blocks from the
404 // problem.
405 std::vector<ResidualBlockId> residual_blocks;
406
407 // Even though the residual blocks in the problem may contain loss
408 // functions, setting apply_loss_function to false will turn off
409 // the application of the loss function to the output of the cost
410 // function. This is of use for example if the user wishes to
411 // analyse the solution quality by studying the distribution of
412 // residuals before and after the solve.
413 bool apply_loss_function = true;
414
415 int num_threads = 1;
416 };
417
418 // Evaluate Problem. Any of the output pointers can be NULL. Which
419 // residual blocks and parameter blocks are used is controlled by
420 // the EvaluateOptions struct above.
421 //
422 // Note 1: The evaluation will use the values stored in the memory
423 // locations pointed to by the parameter block pointers used at the
424 // time of the construction of the problem. i.e.,
425 //
426 // Problem problem;
427 // double x = 1;
428 // problem.AddResidualBlock(new MyCostFunction, NULL, &x);
429 //
430 // double cost = 0.0;
431 // problem.Evaluate(Problem::EvaluateOptions(), &cost, NULL, NULL, NULL);
432 //
433 // The cost is evaluated at x = 1. If you wish to evaluate the
434 // problem at x = 2, then
435 //
436 // x = 2;
437 // problem.Evaluate(Problem::EvaluateOptions(), &cost, NULL, NULL, NULL);
438 //
439 // is the way to do so.
440 //
441 // Note 2: If no local parameterizations are used, then the size of
442 // the gradient vector (and the number of columns in the jacobian)
443 // is the sum of the sizes of all the parameter blocks. If a
444 // parameter block has a local parameterization, then it contributes
445 // "LocalSize" entries to the gradient vector (and the number of
446 // columns in the jacobian).
447 //
448 // Note 3: This function cannot be called while the problem is being
449 // solved, for example it cannot be called from an IterationCallback
450 // at the end of an iteration during a solve.
451 bool Evaluate(const EvaluateOptions& options,
452 double* cost,
453 std::vector<double>* residuals,
454 std::vector<double>* gradient,
455 CRSMatrix* jacobian);
456
457 private:
458 friend class Solver;
459 friend class Covariance;
460 std::unique_ptr<internal::ProblemImpl> problem_impl_;
461};
462
463} // namespace ceres
464
465#include "ceres/internal/reenable_warnings.h"
466
467#endif // CERES_PUBLIC_PROBLEM_H_