Imported Upstream version ceres 1.13.0
[platform/upstream/ceres-solver.git] / internal / ceres / program.h
1 // 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: keir@google.com (Keir Mierle)
30
31 #ifndef CERES_INTERNAL_PROGRAM_H_
32 #define CERES_INTERNAL_PROGRAM_H_
33
34 #include <set>
35 #include <string>
36 #include <vector>
37 #include "ceres/internal/port.h"
38
39 namespace ceres {
40 namespace internal {
41
42 class ParameterBlock;
43 class ProblemImpl;
44 class ResidualBlock;
45 class TripletSparseMatrix;
46
47 // A nonlinear least squares optimization problem. This is different from the
48 // similarly-named "Problem" object, which offers a mutation interface for
49 // adding and modifying parameters and residuals. The Program contains the core
50 // part of the Problem, which is the parameters and the residuals, stored in a
51 // particular ordering. The ordering is critical, since it defines the mapping
52 // between (residual, parameter) pairs and a position in the jacobian of the
53 // objective function. Various parts of Ceres transform one Program into
54 // another; for example, the first stage of solving involves stripping all
55 // constant parameters and residuals. This is in contrast with Problem, which is
56 // not built for transformation.
57 class Program {
58  public:
59   Program();
60   explicit Program(const Program& program);
61
62   // The ordered parameter and residual blocks for the program.
63   const std::vector<ParameterBlock*>& parameter_blocks() const;
64   const std::vector<ResidualBlock*>& residual_blocks() const;
65   std::vector<ParameterBlock*>* mutable_parameter_blocks();
66   std::vector<ResidualBlock*>* mutable_residual_blocks();
67
68   // Serialize to/from the program and update states.
69   //
70   // NOTE: Setting the state of a parameter block can trigger the
71   // computation of the Jacobian of its local parameterization. If
72   // this computation fails for some reason, then this method returns
73   // false and the state of the parameter blocks cannot be trusted.
74   bool StateVectorToParameterBlocks(const double *state);
75   void ParameterBlocksToStateVector(double *state) const;
76
77   // Copy internal state to the user's parameters.
78   void CopyParameterBlockStateToUserState();
79
80   // Set the parameter block pointers to the user pointers. Since this
81   // runs parameter block set state internally, which may call local
82   // parameterizations, this can fail. False is returned on failure.
83   bool SetParameterBlockStatePtrsToUserStatePtrs();
84
85   // Update a state vector for the program given a delta.
86   bool Plus(const double* state,
87             const double* delta,
88             double* state_plus_delta) const;
89
90   // Set the parameter indices and offsets. This permits mapping backward
91   // from a ParameterBlock* to an index in the parameter_blocks() vector. For
92   // any parameter block p, after calling SetParameterOffsetsAndIndex(), it
93   // is true that
94   //
95   //   parameter_blocks()[p->index()] == p
96   //
97   // If a parameter appears in a residual but not in the parameter block, then
98   // it will have an index of -1.
99   //
100   // This also updates p->state_offset() and p->delta_offset(), which are the
101   // position of the parameter in the state and delta vector respectively.
102   void SetParameterOffsetsAndIndex();
103
104   // Check if the internal state of the program (the indexing and the
105   // offsets) are correct.
106   bool IsValid() const;
107
108   bool ParameterBlocksAreFinite(std::string* message) const;
109
110   // Returns true if the program has any non-constant parameter blocks
111   // which have non-trivial bounds constraints.
112   bool IsBoundsConstrained() const;
113
114   // Returns false, if the program has any constant parameter blocks
115   // which are not feasible, or any variable parameter blocks which
116   // have a lower bound greater than or equal to the upper bound.
117   bool IsFeasible(std::string* message) const;
118
119   // Loop over each residual block and ensure that no two parameter
120   // blocks in the same residual block are part of
121   // parameter_blocks as that would violate the assumption that it
122   // is an independent set in the Hessian matrix.
123   bool IsParameterBlockSetIndependent(
124       const std::set<double*>& independent_set) const;
125
126   // Create a TripletSparseMatrix which contains the zero-one
127   // structure corresponding to the block sparsity of the transpose of
128   // the Jacobian matrix.
129   //
130   // Caller owns the result.
131   TripletSparseMatrix* CreateJacobianBlockSparsityTranspose() const;
132
133   // Create a copy of this program and removes constant parameter
134   // blocks and residual blocks with no varying parameter blocks while
135   // preserving their relative order.
136   //
137   // removed_parameter_blocks on exit will contain the list of
138   // parameter blocks that were removed.
139   //
140   // fixed_cost will be equal to the sum of the costs of the residual
141   // blocks that were removed.
142   //
143   // If there was a problem, then the function will return a NULL
144   // pointer and error will contain a human readable description of
145   // the problem.
146   Program* CreateReducedProgram(std::vector<double*>* removed_parameter_blocks,
147                                 double* fixed_cost,
148                                 std::string* error) const;
149
150   // See problem.h for what these do.
151   int NumParameterBlocks() const;
152   int NumParameters() const;
153   int NumEffectiveParameters() const;
154   int NumResidualBlocks() const;
155   int NumResiduals() const;
156
157   int MaxScratchDoublesNeededForEvaluate() const;
158   int MaxDerivativesPerResidualBlock() const;
159   int MaxParametersPerResidualBlock() const;
160   int MaxResidualsPerResidualBlock() const;
161
162   // A human-readable dump of the parameter blocks for debugging.
163   // TODO(keir): If necessary, also dump the residual blocks.
164   std::string ToString() const;
165
166  private:
167   // Remove constant parameter blocks and residual blocks with no
168   // varying parameter blocks while preserving their relative order.
169   //
170   // removed_parameter_blocks on exit will contain the list of
171   // parameter blocks that were removed.
172   //
173   // fixed_cost will be equal to the sum of the costs of the residual
174   // blocks that were removed.
175   //
176   // If there was a problem, then the function will return false and
177   // error will contain a human readable description of the problem.
178   bool RemoveFixedBlocks(std::vector<double*>* removed_parameter_blocks,
179                          double* fixed_cost,
180                          std::string* message);
181
182   // The Program does not own the ParameterBlock or ResidualBlock objects.
183   std::vector<ParameterBlock*> parameter_blocks_;
184   std::vector<ResidualBlock*> residual_blocks_;
185
186   friend class ProblemImpl;
187 };
188
189 }  // namespace internal
190 }  // namespace ceres
191
192 #endif  // CERES_INTERNAL_PROGRAM_H_