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