• Home
  • Line#
  • Scopes#
  • Navigate#
  • Raw
  • Download
1 // 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 <map>
36 #include <string>
37 #include <vector>
38 
39 #include "ceres/execution_summary.h"
40 #include "ceres/internal/port.h"
41 #include "ceres/types.h"
42 
43 namespace ceres {
44 
45 struct CRSMatrix;
46 
47 namespace internal {
48 
49 class Program;
50 class SparseMatrix;
51 
52 // The Evaluator interface offers a way to interact with a least squares cost
53 // function that is useful for an optimizer that wants to minimize the least
54 // squares objective. This insulates the optimizer from issues like Jacobian
55 // storage, parameterization, etc.
56 class Evaluator {
57  public:
58   virtual ~Evaluator();
59 
60   struct Options {
OptionsOptions61     Options()
62         : num_threads(1),
63           num_eliminate_blocks(-1),
64           linear_solver_type(DENSE_QR) {}
65 
66     int num_threads;
67     int num_eliminate_blocks;
68     LinearSolverType linear_solver_type;
69   };
70 
71   static Evaluator* Create(const Options& options,
72                            Program* program,
73                            string* error);
74 
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 
100   // 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 
119   // Options struct to control Evaluator::Evaluate;
120   struct EvaluateOptions {
EvaluateOptionsEvaluateOptions121     EvaluateOptions()
122         : apply_loss_function(true) {
123     }
124 
125     // If false, the loss function correction is not applied to the
126     // residual blocks.
127     bool apply_loss_function;
128   };
129 
130   // Evaluate the cost function for the given state. Returns the cost,
131   // residuals, and jacobian in the corresponding arguments. Both residuals and
132   // jacobian are optional; to avoid computing them, pass NULL.
133   //
134   // If non-NULL, the Jacobian must have a suitable sparsity pattern; only the
135   // values array of the jacobian is modified.
136   //
137   // state is an array of size NumParameters(), cost is a pointer to a single
138   // double, and residuals is an array of doubles of size NumResiduals().
139   virtual bool Evaluate(const EvaluateOptions& evaluate_options,
140                         const double* state,
141                         double* cost,
142                         double* residuals,
143                         double* gradient,
144                         SparseMatrix* jacobian) = 0;
145 
146   // Variant of Evaluator::Evaluate where the user wishes to use the
147   // default EvaluateOptions struct. This is mostly here as a
148   // convenience method.
Evaluate(const double * state,double * cost,double * residuals,double * gradient,SparseMatrix * jacobian)149   bool Evaluate(const double* state,
150                 double* cost,
151                 double* residuals,
152                 double* gradient,
153                 SparseMatrix* jacobian) {
154     return Evaluate(EvaluateOptions(),
155                     state,
156                     cost,
157                     residuals,
158                     gradient,
159                     jacobian);
160   }
161 
162   // Make a change delta (of size NumEffectiveParameters()) to state (of size
163   // NumParameters()) and store the result in state_plus_delta.
164   //
165   // In the case that there are no parameterizations used, this is equivalent to
166   //
167   //   state_plus_delta[i] = state[i] + delta[i] ;
168   //
169   // however, the mapping is more complicated in the case of parameterizations
170   // like quaternions. This is the same as the "Plus()" operation in
171   // local_parameterization.h, but operating over the entire state vector for a
172   // problem.
173   virtual bool Plus(const double* state,
174                     const double* delta,
175                     double* state_plus_delta) const = 0;
176 
177   // The number of parameters in the optimization problem.
178   virtual int NumParameters() const = 0;
179 
180   // This is the effective number of parameters that the optimizer may adjust.
181   // This applies when there are parameterizations on some of the parameters.
182   virtual int NumEffectiveParameters()  const = 0;
183 
184   // The number of residuals in the optimization problem.
185   virtual int NumResiduals() const = 0;
186 
187   // The following two methods return copies instead of references so
188   // that the base class implementation does not have to worry about
189   // life time issues. Further, these calls are not expected to be
190   // frequent or performance sensitive.
CallStatistics()191   virtual map<string, int> CallStatistics() const {
192     return map<string, int>();
193   }
194 
TimeStatistics()195   virtual map<string, double> TimeStatistics() const {
196     return map<string, double>();
197   }
198 };
199 
200 }  // namespace internal
201 }  // namespace ceres
202 
203 #endif  // CERES_INTERNAL_EVALUATOR_H_
204