Imported Upstream version ceres 1.13.0
[platform/upstream/ceres-solver.git] / internal / ceres / problem_impl.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 // This is the implementation of the public Problem API. The pointer to
32 // implementation (PIMPL) idiom makes it possible for Ceres internal code to
33 // refer to the private data members without needing to exposing it to the
34 // world. An alternative to PIMPL is to have a factory which returns instances
35 // of a virtual base class; while that approach would work, it requires clients
36 // to always put a Problem object into a scoped pointer; this needlessly muddies
37 // client code for little benefit. Therefore, the PIMPL comprise was chosen.
38
39 #ifndef CERES_PUBLIC_PROBLEM_IMPL_H_
40 #define CERES_PUBLIC_PROBLEM_IMPL_H_
41
42 #include <map>
43 #include <vector>
44
45 #include "ceres/internal/macros.h"
46 #include "ceres/internal/port.h"
47 #include "ceres/internal/scoped_ptr.h"
48 #include "ceres/collections_port.h"
49 #include "ceres/problem.h"
50 #include "ceres/types.h"
51
52 namespace ceres {
53
54 class CostFunction;
55 class LossFunction;
56 class LocalParameterization;
57 struct CRSMatrix;
58
59 namespace internal {
60
61 class Program;
62 class ResidualBlock;
63
64 class ProblemImpl {
65  public:
66   typedef std::map<double*, ParameterBlock*> ParameterMap;
67   typedef HashSet<ResidualBlock*> ResidualBlockSet;
68
69   ProblemImpl();
70   explicit ProblemImpl(const Problem::Options& options);
71
72   ~ProblemImpl();
73
74   // See the public problem.h file for description of these methods.
75   ResidualBlockId AddResidualBlock(
76       CostFunction* cost_function,
77       LossFunction* loss_function,
78       const std::vector<double*>& parameter_blocks);
79   ResidualBlockId AddResidualBlock(CostFunction* cost_function,
80                                    LossFunction* loss_function,
81                                    double* x0);
82   ResidualBlockId AddResidualBlock(CostFunction* cost_function,
83                                    LossFunction* loss_function,
84                                    double* x0, double* x1);
85   ResidualBlockId AddResidualBlock(CostFunction* cost_function,
86                                    LossFunction* loss_function,
87                                    double* x0, double* x1, double* x2);
88   ResidualBlockId AddResidualBlock(CostFunction* cost_function,
89                                    LossFunction* loss_function,
90                                    double* x0, double* x1, double* x2,
91                                    double* x3);
92   ResidualBlockId AddResidualBlock(CostFunction* cost_function,
93                                    LossFunction* loss_function,
94                                    double* x0, double* x1, double* x2,
95                                    double* x3, double* x4);
96   ResidualBlockId AddResidualBlock(CostFunction* cost_function,
97                                    LossFunction* loss_function,
98                                    double* x0, double* x1, double* x2,
99                                    double* x3, double* x4, double* x5);
100   ResidualBlockId AddResidualBlock(CostFunction* cost_function,
101                                    LossFunction* loss_function,
102                                    double* x0, double* x1, double* x2,
103                                    double* x3, double* x4, double* x5,
104                                    double* x6);
105   ResidualBlockId AddResidualBlock(CostFunction* cost_function,
106                                    LossFunction* loss_function,
107                                    double* x0, double* x1, double* x2,
108                                    double* x3, double* x4, double* x5,
109                                    double* x6, double* x7);
110   ResidualBlockId AddResidualBlock(CostFunction* cost_function,
111                                    LossFunction* loss_function,
112                                    double* x0, double* x1, double* x2,
113                                    double* x3, double* x4, double* x5,
114                                    double* x6, double* x7, double* x8);
115   ResidualBlockId AddResidualBlock(CostFunction* cost_function,
116                                    LossFunction* loss_function,
117                                    double* x0, double* x1, double* x2,
118                                    double* x3, double* x4, double* x5,
119                                    double* x6, double* x7, double* x8,
120                                    double* x9);
121   void AddParameterBlock(double* values, int size);
122   void AddParameterBlock(double* values,
123                          int size,
124                          LocalParameterization* local_parameterization);
125
126   void RemoveResidualBlock(ResidualBlock* residual_block);
127   void RemoveParameterBlock(double* values);
128
129   void SetParameterBlockConstant(double* values);
130   void SetParameterBlockVariable(double* values);
131   bool IsParameterBlockConstant(double* values) const;
132
133   void SetParameterization(double* values,
134                            LocalParameterization* local_parameterization);
135   const LocalParameterization* GetParameterization(double* values) const;
136
137   void SetParameterLowerBound(double* values, int index, double lower_bound);
138   void SetParameterUpperBound(double* values, int index, double upper_bound);
139
140   bool Evaluate(const Problem::EvaluateOptions& options,
141                 double* cost,
142                 std::vector<double>* residuals,
143                 std::vector<double>* gradient,
144                 CRSMatrix* jacobian);
145
146   int NumParameterBlocks() const;
147   int NumParameters() const;
148   int NumResidualBlocks() const;
149   int NumResiduals() const;
150
151   int ParameterBlockSize(const double* parameter_block) const;
152   int ParameterBlockLocalSize(const double* parameter_block) const;
153
154   bool HasParameterBlock(const double* parameter_block) const;
155
156   void GetParameterBlocks(std::vector<double*>* parameter_blocks) const;
157   void GetResidualBlocks(std::vector<ResidualBlockId>* residual_blocks) const;
158
159   void GetParameterBlocksForResidualBlock(
160       const ResidualBlockId residual_block,
161       std::vector<double*>* parameter_blocks) const;
162
163   const CostFunction* GetCostFunctionForResidualBlock(
164       const ResidualBlockId residual_block) const;
165   const LossFunction* GetLossFunctionForResidualBlock(
166       const ResidualBlockId residual_block) const;
167
168   void GetResidualBlocksForParameterBlock(
169       const double* values,
170       std::vector<ResidualBlockId>* residual_blocks) const;
171
172   const Program& program() const { return *program_; }
173   Program* mutable_program() { return program_.get(); }
174
175   const ParameterMap& parameter_map() const { return parameter_block_map_; }
176   const ResidualBlockSet& residual_block_set() const {
177     CHECK(options_.enable_fast_removal)
178         << "Fast removal not enabled, residual_block_set is not maintained.";
179     return residual_block_set_;
180   }
181
182  private:
183   ParameterBlock* InternalAddParameterBlock(double* values, int size);
184   void InternalRemoveResidualBlock(ResidualBlock* residual_block);
185
186   bool InternalEvaluate(Program* program,
187                         double* cost,
188                         std::vector<double>* residuals,
189                         std::vector<double>* gradient,
190                         CRSMatrix* jacobian);
191
192   // Delete the arguments in question. These differ from the Remove* functions
193   // in that they do not clean up references to the block to delete; they
194   // merely delete them.
195   template<typename Block>
196   void DeleteBlockInVector(std::vector<Block*>* mutable_blocks,
197                            Block* block_to_remove);
198   void DeleteBlock(ResidualBlock* residual_block);
199   void DeleteBlock(ParameterBlock* parameter_block);
200
201   const Problem::Options options_;
202
203   // The mapping from user pointers to parameter blocks.
204   std::map<double*, ParameterBlock*> parameter_block_map_;
205
206   // Iff enable_fast_removal is enabled, contains the current residual blocks.
207   ResidualBlockSet residual_block_set_;
208
209   // The actual parameter and residual blocks.
210   internal::scoped_ptr<internal::Program> program_;
211
212   // When removing residual and parameter blocks, cost/loss functions and
213   // parameterizations have ambiguous ownership. Instead of scanning the entire
214   // problem to see if the cost/loss/parameterization is shared with other
215   // residual or parameter blocks, buffer them until destruction.
216   //
217   // TODO(keir): See if it makes sense to use sets instead.
218   std::vector<CostFunction*> cost_functions_to_delete_;
219   std::vector<LossFunction*> loss_functions_to_delete_;
220   std::vector<LocalParameterization*> local_parameterizations_to_delete_;
221
222   CERES_DISALLOW_COPY_AND_ASSIGN(ProblemImpl);
223 };
224
225 }  // namespace internal
226 }  // namespace ceres
227
228 #endif  // CERES_PUBLIC_PROBLEM_IMPL_H_