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
31#ifndef CERES_INTERNAL_MINIMIZER_H_
32#define CERES_INTERNAL_MINIMIZER_H_
33
34#include <string>
35#include <vector>
36#include "ceres/internal/port.h"
37#include "ceres/iteration_callback.h"
38#include "ceres/solver.h"
39
40namespace ceres {
41namespace internal {
42
43class Evaluator;
44class LinearSolver;
45class SparseMatrix;
46class TrustRegionStrategy;
47
48// Interface for non-linear least squares solvers.
49class Minimizer {
50 public:
51  // Options struct to control the behaviour of the Minimizer. Please
52  // see solver.h for detailed information about the meaning and
53  // default values of each of these parameters.
54  struct Options {
55    Options() {
56      Init(Solver::Options());
57    }
58
59    explicit Options(const Solver::Options& options) {
60      Init(options);
61    }
62
63    void Init(const Solver::Options& options) {
64      num_threads = options.num_threads;
65      max_num_iterations = options.max_num_iterations;
66      max_solver_time_in_seconds = options.max_solver_time_in_seconds;
67      max_step_solver_retries = 5;
68      gradient_tolerance = options.gradient_tolerance;
69      parameter_tolerance = options.parameter_tolerance;
70      function_tolerance = options.function_tolerance;
71      min_relative_decrease = options.min_relative_decrease;
72      eta = options.eta;
73      jacobi_scaling = options.jacobi_scaling;
74      use_nonmonotonic_steps = options.use_nonmonotonic_steps;
75      max_consecutive_nonmonotonic_steps =
76          options.max_consecutive_nonmonotonic_steps;
77      trust_region_problem_dump_directory =
78          options.trust_region_problem_dump_directory;
79      trust_region_minimizer_iterations_to_dump =
80          options.trust_region_minimizer_iterations_to_dump;
81      trust_region_problem_dump_format_type =
82          options.trust_region_problem_dump_format_type;
83      max_num_consecutive_invalid_steps =
84          options.max_num_consecutive_invalid_steps;
85      min_trust_region_radius = options.min_trust_region_radius;
86      line_search_direction_type = options.line_search_direction_type;
87      line_search_type = options.line_search_type;
88      nonlinear_conjugate_gradient_type =
89          options.nonlinear_conjugate_gradient_type;
90      max_lbfgs_rank = options.max_lbfgs_rank;
91      use_approximate_eigenvalue_bfgs_scaling =
92          options.use_approximate_eigenvalue_bfgs_scaling;
93      line_search_interpolation_type =
94          options.line_search_interpolation_type;
95      min_line_search_step_size = options.min_line_search_step_size;
96      line_search_sufficient_function_decrease =
97          options.line_search_sufficient_function_decrease;
98      max_line_search_step_contraction =
99          options.max_line_search_step_contraction;
100      min_line_search_step_contraction =
101          options.min_line_search_step_contraction;
102      max_num_line_search_step_size_iterations =
103          options.max_num_line_search_step_size_iterations;
104      max_num_line_search_direction_restarts =
105          options.max_num_line_search_direction_restarts;
106      line_search_sufficient_curvature_decrease =
107          options.line_search_sufficient_curvature_decrease;
108      max_line_search_step_expansion =
109          options.max_line_search_step_expansion;
110      evaluator = NULL;
111      trust_region_strategy = NULL;
112      jacobian = NULL;
113      callbacks = options.callbacks;
114      inner_iteration_minimizer = NULL;
115      inner_iteration_tolerance = options.inner_iteration_tolerance;
116    }
117
118    int max_num_iterations;
119    double max_solver_time_in_seconds;
120    int num_threads;
121
122    // Number of times the linear solver should be retried in case of
123    // numerical failure. The retries are done by exponentially scaling up
124    // mu at each retry. This leads to stronger and stronger
125    // regularization making the linear least squares problem better
126    // conditioned at each retry.
127    int max_step_solver_retries;
128    double gradient_tolerance;
129    double parameter_tolerance;
130    double function_tolerance;
131    double min_relative_decrease;
132    double eta;
133    bool jacobi_scaling;
134    bool use_nonmonotonic_steps;
135    int max_consecutive_nonmonotonic_steps;
136    vector<int> trust_region_minimizer_iterations_to_dump;
137    DumpFormatType trust_region_problem_dump_format_type;
138    string trust_region_problem_dump_directory;
139    int max_num_consecutive_invalid_steps;
140    double min_trust_region_radius;
141    LineSearchDirectionType line_search_direction_type;
142    LineSearchType line_search_type;
143    NonlinearConjugateGradientType nonlinear_conjugate_gradient_type;
144    int max_lbfgs_rank;
145    bool use_approximate_eigenvalue_bfgs_scaling;
146    LineSearchInterpolationType line_search_interpolation_type;
147    double min_line_search_step_size;
148    double line_search_sufficient_function_decrease;
149    double max_line_search_step_contraction;
150    double min_line_search_step_contraction;
151    int max_num_line_search_step_size_iterations;
152    int max_num_line_search_direction_restarts;
153    double line_search_sufficient_curvature_decrease;
154    double max_line_search_step_expansion;
155
156
157    // List of callbacks that are executed by the Minimizer at the end
158    // of each iteration.
159    //
160    // The Options struct does not own these pointers.
161    vector<IterationCallback*> callbacks;
162
163    // Object responsible for evaluating the cost, residuals and
164    // Jacobian matrix. The Options struct does not own this pointer.
165    Evaluator* evaluator;
166
167    // Object responsible for actually computing the trust region
168    // step, and sizing the trust region radius. The Options struct
169    // does not own this pointer.
170    TrustRegionStrategy* trust_region_strategy;
171
172    // Object holding the Jacobian matrix. It is assumed that the
173    // sparsity structure of the matrix has already been initialized
174    // and will remain constant for the life time of the
175    // optimization. The Options struct does not own this pointer.
176    SparseMatrix* jacobian;
177
178    Minimizer* inner_iteration_minimizer;
179    double inner_iteration_tolerance;
180  };
181
182  static bool RunCallbacks(const vector<IterationCallback*> callbacks,
183                           const IterationSummary& iteration_summary,
184                           Solver::Summary* summary);
185
186  virtual ~Minimizer();
187  // Note: The minimizer is expected to update the state of the
188  // parameters array every iteration. This is required for the
189  // StateUpdatingCallback to work.
190  virtual void Minimize(const Options& options,
191                        double* parameters,
192                        Solver::Summary* summary) = 0;
193};
194
195}  // namespace internal
196}  // namespace ceres
197
198#endif  // CERES_INTERNAL_MINIMIZER_H_
199