1 // Ceres Solver - A fast non-linear least squares minimizer
2 // Copyright 2015 Google Inc. All rights reserved.
3 // http://ceres-solver.org/
5 // Redistribution and use in source and binary forms, with or without
6 // modification, are permitted provided that the following conditions are met:
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.
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.
29 // Author: sameeragarwal@google.com (Sameer Agarwal)
31 #ifndef CERES_PUBLIC_GRADIENT_PROBLEM_H_
32 #define CERES_PUBLIC_GRADIENT_PROBLEM_H_
34 #include "ceres/internal/macros.h"
35 #include "ceres/internal/port.h"
36 #include "ceres/internal/scoped_ptr.h"
37 #include "ceres/local_parameterization.h"
41 class FirstOrderFunction;
43 // Instances of GradientProblem represent general non-linear
44 // optimization problems that must be solved using just the value of
45 // the objective function and its gradient. Unlike the Problem class,
46 // which can only be used to model non-linear least squares problems,
47 // instances of GradientProblem not restricted in the form of the
48 // objective function.
50 // Structurally GradientProblem is a composition of a
51 // FirstOrderFunction and optionally a LocalParameterization.
53 // The FirstOrderFunction is responsible for evaluating the cost and
54 // gradient of the objective function.
56 // The LocalParameterization is responsible for going back and forth
57 // between the ambient space and the local tangent space. (See
58 // local_parameterization.h for more details). When a
59 // LocalParameterization is not provided, then the tangent space is
60 // assumed to coincide with the ambient Euclidean space that the
61 // gradient vector lives in.
65 // The following demonstrate the problem construction for Rosenbrock's function
67 // f(x,y) = (1-x)^2 + 100(y - x^2)^2;
69 // class Rosenbrock : public ceres::FirstOrderFunction {
71 // virtual ~Rosenbrock() {}
73 // virtual bool Evaluate(const double* parameters,
75 // double* gradient) const {
76 // const double x = parameters[0];
77 // const double y = parameters[1];
79 // cost[0] = (1.0 - x) * (1.0 - x) + 100.0 * (y - x * x) * (y - x * x);
80 // if (gradient != NULL) {
81 // gradient[0] = -2.0 * (1.0 - x) - 200.0 * (y - x * x) * 2.0 * x;
82 // gradient[1] = 200.0 * (y - x * x);
87 // virtual int NumParameters() const { return 2; };
90 // ceres::GradientProblem problem(new Rosenbrock());
91 class CERES_EXPORT GradientProblem {
93 // Takes ownership of the function.
94 explicit GradientProblem(FirstOrderFunction* function);
96 // Takes ownership of the function and the parameterization.
97 GradientProblem(FirstOrderFunction* function,
98 LocalParameterization* parameterization);
100 int NumParameters() const;
101 int NumLocalParameters() const;
103 // This call is not thread safe.
104 bool Evaluate(const double* parameters, double* cost, double* gradient) const;
105 bool Plus(const double* x, const double* delta, double* x_plus_delta) const;
108 internal::scoped_ptr<FirstOrderFunction> function_;
109 internal::scoped_ptr<LocalParameterization> parameterization_;
110 internal::scoped_array<double> scratch_;
113 // A FirstOrderFunction object implements the evaluation of a function
115 class CERES_EXPORT FirstOrderFunction {
117 virtual ~FirstOrderFunction() {}
118 // cost is never NULL. gradient may be null.
119 virtual bool Evaluate(const double* const parameters,
121 double* gradient) const = 0;
122 virtual int NumParameters() const = 0;
127 #endif // CERES_PUBLIC_GRADIENT_PROBLEM_H_