blob: 76cf1e0718e20d5d63e5170f96079493ce2080d8 [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#ifndef CERES_INTERNAL_EVALUATOR_H_
33#define CERES_INTERNAL_EVALUATOR_H_
34
35#include <string>
Sameer Agarwal4997cbc2012-07-02 12:44:34 -070036#include <vector>
Sameer Agarwalbdd87c02013-01-29 16:24:31 -080037
38#include "ceres/execution_summary.h"
Keir Mierle8ebb0732012-04-30 23:09:08 -070039#include "ceres/internal/port.h"
40#include "ceres/types.h"
41
42namespace ceres {
Sameer Agarwal4997cbc2012-07-02 12:44:34 -070043
Sameer Agarwald3ace022012-09-23 18:19:49 -070044struct CRSMatrix;
Sameer Agarwal4997cbc2012-07-02 12:44:34 -070045
Keir Mierle8ebb0732012-04-30 23:09:08 -070046namespace internal {
47
48class Program;
49class SparseMatrix;
50
51// The Evaluator interface offers a way to interact with a least squares cost
52// function that is useful for an optimizer that wants to minimize the least
53// squares objective. This insulates the optimizer from issues like Jacobian
54// storage, parameterization, etc.
55class Evaluator {
56 public:
57 virtual ~Evaluator();
58
59 struct Options {
60 Options()
61 : num_threads(1),
62 num_eliminate_blocks(-1),
63 linear_solver_type(DENSE_QR) {}
64
65 int num_threads;
66 int num_eliminate_blocks;
67 LinearSolverType linear_solver_type;
68 };
69
70 static Evaluator* Create(const Options& options,
71 Program* program,
72 string* error);
73
Sameer Agarwal4997cbc2012-07-02 12:44:34 -070074
75 // This is used for computing the cost, residual and Jacobian for
76 // returning to the user. For actually solving the optimization
77 // problem, the optimization algorithm uses the ProgramEvaluator
78 // objects directly.
79 //
80 // The residual, gradients and jacobian pointers can be NULL, in
81 // which case they will not be evaluated. cost cannot be NULL.
82 //
83 // The parallelism of the evaluator is controlled by num_threads; it
84 // should be at least 1.
85 //
86 // Note: That this function does not take a parameter vector as
87 // input. The parameter blocks are evaluated on the values contained
88 // in the arrays pointed to by their user_state pointers.
89 //
90 // Also worth noting is that this function mutates program by
91 // calling Program::SetParameterOffsetsAndIndex() on it so that an
92 // evaluator object can be constructed.
93 static bool Evaluate(Program* program,
94 int num_threads,
95 double* cost,
96 vector<double>* residuals,
97 vector<double>* gradient,
98 CRSMatrix* jacobian);
99
Keir Mierle8ebb0732012-04-30 23:09:08 -0700100 // Build and return a sparse matrix for storing and working with the Jacobian
101 // of the objective function. The jacobian has dimensions
102 // NumEffectiveParameters() by NumParameters(), and is typically extremely
103 // sparse. Since the sparsity pattern of the Jacobian remains constant over
104 // the lifetime of the optimization problem, this method is used to
105 // instantiate a SparseMatrix object with the appropriate sparsity structure
106 // (which can be an expensive operation) and then reused by the optimization
107 // algorithm and the various linear solvers.
108 //
109 // It is expected that the classes implementing this interface will be aware
110 // of their client's requirements for the kind of sparse matrix storage and
111 // layout that is needed for an efficient implementation. For example
112 // CompressedRowOptimizationProblem creates a compressed row representation of
113 // the jacobian for use with CHOLMOD, where as BlockOptimizationProblem
114 // creates a BlockSparseMatrix representation of the jacobian for use in the
115 // Schur complement based methods.
116 virtual SparseMatrix* CreateJacobian() const = 0;
117
118 // Evaluate the cost function for the given state. Returns the cost,
119 // residuals, and jacobian in the corresponding arguments. Both residuals and
120 // jacobian are optional; to avoid computing them, pass NULL.
121 //
122 // If non-NULL, the Jacobian must have a suitable sparsity pattern; only the
123 // values array of the jacobian is modified.
124 //
125 // state is an array of size NumParameters(), cost is a pointer to a single
126 // double, and residuals is an array of doubles of size NumResiduals().
127 virtual bool Evaluate(const double* state,
128 double* cost,
129 double* residuals,
Keir Mierlef44907f2012-07-06 13:52:32 -0700130 double* gradient,
Keir Mierle8ebb0732012-04-30 23:09:08 -0700131 SparseMatrix* jacobian) = 0;
132
133 // Make a change delta (of size NumEffectiveParameters()) to state (of size
134 // NumParameters()) and store the result in state_plus_delta.
135 //
136 // In the case that there are no parameterizations used, this is equivalent to
137 //
138 // state_plus_delta[i] = state[i] + delta[i] ;
139 //
140 // however, the mapping is more complicated in the case of parameterizations
141 // like quaternions. This is the same as the "Plus()" operation in
142 // local_parameterization.h, but operating over the entire state vector for a
143 // problem.
144 virtual bool Plus(const double* state,
145 const double* delta,
146 double* state_plus_delta) const = 0;
147
148 // The number of parameters in the optimization problem.
149 virtual int NumParameters() const = 0;
150
151 // This is the effective number of parameters that the optimizer may adjust.
152 // This applies when there are parameterizations on some of the parameters.
153 virtual int NumEffectiveParameters() const = 0;
154
155 // The number of residuals in the optimization problem.
156 virtual int NumResiduals() const = 0;
Sameer Agarwalbdd87c02013-01-29 16:24:31 -0800157
Sameer Agarwal42a84b82013-02-01 12:22:53 -0800158 // The following two methods return copies instead of references so
159 // that the base class implementation does not have to worry about
160 // life time issues. Further, these calls are not expected to be
161 // frequent or performance sensitive.
162 virtual map<string, int> CallStatistics() const {
163 return map<string, int>();
164 }
165
166 virtual map<string, double> TimeStatistics() const {
167 return map<string, double>();
168 }
Keir Mierle8ebb0732012-04-30 23:09:08 -0700169};
170
171} // namespace internal
172} // namespace ceres
173
174#endif // CERES_INTERNAL_EVALUATOR_H_