blob: 88943dfbcff536eab149554fa175591a3f2b1c02 [file] [log] [blame]
Austin Schuh70cc9552019-01-21 19:46:48 -08001// Ceres Solver - A fast non-linear least squares minimizer
Austin Schuh1d1e6ea2020-12-23 21:56:30 -08002// Copyright 2019 Google Inc. All rights reserved.
Austin Schuh70cc9552019-01-21 19:46:48 -08003// 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_PARAMETER_BLOCK_H_
32#define CERES_INTERNAL_PARAMETER_BLOCK_H_
33
34#include <algorithm>
35#include <cstdint>
36#include <cstdlib>
37#include <limits>
38#include <memory>
39#include <string>
40#include <unordered_set>
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080041
Austin Schuh70cc9552019-01-21 19:46:48 -080042#include "ceres/array_utils.h"
43#include "ceres/internal/eigen.h"
44#include "ceres/internal/port.h"
45#include "ceres/local_parameterization.h"
46#include "ceres/stringprintf.h"
47#include "glog/logging.h"
48
49namespace ceres {
50namespace internal {
51
52class ProblemImpl;
53class ResidualBlock;
54
55// The parameter block encodes the location of the user's original value, and
56// also the "current state" of the parameter. The evaluator uses whatever is in
57// the current state of the parameter when evaluating. This is inlined since the
58// methods are performance sensitive.
59//
60// The class is not thread-safe, unless only const methods are called. The
61// parameter block may also hold a pointer to a local parameterization; the
62// parameter block does not take ownership of this pointer, so the user is
63// responsible for the proper disposal of the local parameterization.
64class ParameterBlock {
65 public:
66 typedef std::unordered_set<ResidualBlock*> ResidualBlockSet;
67
68 // Create a parameter block with the user state, size, and index specified.
69 // The size is the size of the parameter block and the index is the position
70 // of the parameter block inside a Program (if any).
71 ParameterBlock(double* user_state, int size, int index)
72 : user_state_(user_state),
Austin Schuh70cc9552019-01-21 19:46:48 -080073 size_(size),
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080074 state_(user_state),
Austin Schuh70cc9552019-01-21 19:46:48 -080075 index_(index) {}
76
77 ParameterBlock(double* user_state,
78 int size,
79 int index,
80 LocalParameterization* local_parameterization)
81 : user_state_(user_state),
Austin Schuh70cc9552019-01-21 19:46:48 -080082 size_(size),
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080083 state_(user_state),
Austin Schuh70cc9552019-01-21 19:46:48 -080084 index_(index) {
85 if (local_parameterization != nullptr) {
86 SetParameterization(local_parameterization);
87 }
88 }
89
90 // The size of the parameter block.
91 int Size() const { return size_; }
92
93 // Manipulate the parameter state.
94 bool SetState(const double* x) {
95 CHECK(x != nullptr) << "Tried to set the state of constant parameter "
96 << "with user location " << user_state_;
97 CHECK(!IsConstant()) << "Tried to set the state of constant parameter "
98 << "with user location " << user_state_;
99
100 state_ = x;
101 return UpdateLocalParameterizationJacobian();
102 }
103
104 // Copy the current parameter state out to x. This is "GetState()" rather than
105 // simply "state()" since it is actively copying the data into the passed
106 // pointer.
107 void GetState(double* x) const {
108 if (x != state_) {
109 std::copy(state_, state_ + size_, x);
110 }
111 }
112
113 // Direct pointers to the current state.
114 const double* state() const { return state_; }
115 const double* user_state() const { return user_state_; }
116 double* mutable_user_state() { return user_state_; }
117 const LocalParameterization* local_parameterization() const {
118 return local_parameterization_;
119 }
120 LocalParameterization* mutable_local_parameterization() {
121 return local_parameterization_;
122 }
123
124 // Set this parameter block to vary or not.
125 void SetConstant() { is_set_constant_ = true; }
126 void SetVarying() { is_set_constant_ = false; }
Austin Schuh70cc9552019-01-21 19:46:48 -0800127 bool IsConstant() const { return (is_set_constant_ || LocalSize() == 0); }
128
129 double UpperBound(int index) const {
130 return (upper_bounds_ ? upper_bounds_[index]
131 : std::numeric_limits<double>::max());
132 }
133
134 double LowerBound(int index) const {
135 return (lower_bounds_ ? lower_bounds_[index]
136 : -std::numeric_limits<double>::max());
137 }
138
139 bool IsUpperBounded() const { return (upper_bounds_ == nullptr); }
140 bool IsLowerBounded() const { return (lower_bounds_ == nullptr); }
141
142 // This parameter block's index in an array.
143 int index() const { return index_; }
144 void set_index(int index) { index_ = index; }
145
146 // This parameter offset inside a larger state vector.
147 int state_offset() const { return state_offset_; }
148 void set_state_offset(int state_offset) { state_offset_ = state_offset; }
149
150 // This parameter offset inside a larger delta vector.
151 int delta_offset() const { return delta_offset_; }
152 void set_delta_offset(int delta_offset) { delta_offset_ = delta_offset; }
153
154 // Methods relating to the parameter block's parameterization.
155
156 // The local to global jacobian. Returns nullptr if there is no local
157 // parameterization for this parameter block. The returned matrix is row-major
158 // and has Size() rows and LocalSize() columns.
159 const double* LocalParameterizationJacobian() const {
160 return local_parameterization_jacobian_.get();
161 }
162
163 int LocalSize() const {
164 return (local_parameterization_ == nullptr)
165 ? size_
166 : local_parameterization_->LocalSize();
167 }
168
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800169 // Set the parameterization. The parameter block does not take
170 // ownership of the parameterization.
Austin Schuh70cc9552019-01-21 19:46:48 -0800171 void SetParameterization(LocalParameterization* new_parameterization) {
Austin Schuh70cc9552019-01-21 19:46:48 -0800172 // Nothing to do if the new parameterization is the same as the
173 // old parameterization.
174 if (new_parameterization == local_parameterization_) {
175 return;
176 }
177
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800178 if (new_parameterization == nullptr) {
179 local_parameterization_ = nullptr;
180 return;
181 }
Austin Schuh70cc9552019-01-21 19:46:48 -0800182
183 CHECK(new_parameterization->GlobalSize() == size_)
184 << "Invalid parameterization for parameter block. The parameter block "
185 << "has size " << size_ << " while the parameterization has a global "
186 << "size of " << new_parameterization->GlobalSize() << ". Did you "
187 << "accidentally use the wrong parameter block or parameterization?";
188
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800189 CHECK_GE(new_parameterization->LocalSize(), 0)
Austin Schuh70cc9552019-01-21 19:46:48 -0800190 << "Invalid parameterization. Parameterizations must have a "
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800191 << "non-negative dimensional tangent space.";
Austin Schuh70cc9552019-01-21 19:46:48 -0800192
193 local_parameterization_ = new_parameterization;
194 local_parameterization_jacobian_.reset(
195 new double[local_parameterization_->GlobalSize() *
196 local_parameterization_->LocalSize()]);
197 CHECK(UpdateLocalParameterizationJacobian())
198 << "Local parameterization Jacobian computation failed for x: "
199 << ConstVectorRef(state_, Size()).transpose();
200 }
201
202 void SetUpperBound(int index, double upper_bound) {
203 CHECK_LT(index, size_);
204
205 if (upper_bound >= std::numeric_limits<double>::max() && !upper_bounds_) {
206 return;
207 }
208
209 if (!upper_bounds_) {
210 upper_bounds_.reset(new double[size_]);
211 std::fill(upper_bounds_.get(),
212 upper_bounds_.get() + size_,
213 std::numeric_limits<double>::max());
214 }
215
216 upper_bounds_[index] = upper_bound;
217 }
218
219 void SetLowerBound(int index, double lower_bound) {
220 CHECK_LT(index, size_);
221
222 if (lower_bound <= -std::numeric_limits<double>::max() && !lower_bounds_) {
223 return;
224 }
225
226 if (!lower_bounds_) {
227 lower_bounds_.reset(new double[size_]);
228 std::fill(lower_bounds_.get(),
229 lower_bounds_.get() + size_,
230 -std::numeric_limits<double>::max());
231 }
232
233 lower_bounds_[index] = lower_bound;
234 }
235
236 // Generalization of the addition operation. This is the same as
237 // LocalParameterization::Plus() followed by projection onto the
238 // hyper cube implied by the bounds constraints.
239 bool Plus(const double* x, const double* delta, double* x_plus_delta) {
240 if (local_parameterization_ != nullptr) {
241 if (!local_parameterization_->Plus(x, delta, x_plus_delta)) {
242 return false;
243 }
244 } else {
245 VectorRef(x_plus_delta, size_) =
246 ConstVectorRef(x, size_) + ConstVectorRef(delta, size_);
247 }
248
249 // Project onto the box constraints.
250 if (lower_bounds_.get() != nullptr) {
251 for (int i = 0; i < size_; ++i) {
252 x_plus_delta[i] = std::max(x_plus_delta[i], lower_bounds_[i]);
253 }
254 }
255
256 if (upper_bounds_.get() != nullptr) {
257 for (int i = 0; i < size_; ++i) {
258 x_plus_delta[i] = std::min(x_plus_delta[i], upper_bounds_[i]);
259 }
260 }
261
262 return true;
263 }
264
265 std::string ToString() const {
266 return StringPrintf(
267 "{ this=%p, user_state=%p, state=%p, size=%d, "
268 "constant=%d, index=%d, state_offset=%d, "
269 "delta_offset=%d }",
270 this,
271 user_state_,
272 state_,
273 size_,
274 is_set_constant_,
275 index_,
276 state_offset_,
277 delta_offset_);
278 }
279
280 void EnableResidualBlockDependencies() {
281 CHECK(residual_blocks_.get() == nullptr)
282 << "Ceres bug: There is already a residual block collection "
283 << "for parameter block: " << ToString();
284 residual_blocks_.reset(new ResidualBlockSet);
285 }
286
287 void AddResidualBlock(ResidualBlock* residual_block) {
288 CHECK(residual_blocks_.get() != nullptr)
289 << "Ceres bug: The residual block collection is null for parameter "
290 << "block: " << ToString();
291 residual_blocks_->insert(residual_block);
292 }
293
294 void RemoveResidualBlock(ResidualBlock* residual_block) {
295 CHECK(residual_blocks_.get() != nullptr)
296 << "Ceres bug: The residual block collection is null for parameter "
297 << "block: " << ToString();
298 CHECK(residual_blocks_->find(residual_block) != residual_blocks_->end())
299 << "Ceres bug: Missing residual for parameter block: " << ToString();
300 residual_blocks_->erase(residual_block);
301 }
302
303 // This is only intended for iterating; perhaps this should only expose
304 // .begin() and .end().
305 ResidualBlockSet* mutable_residual_blocks() { return residual_blocks_.get(); }
306
307 double LowerBoundForParameter(int index) const {
308 if (lower_bounds_.get() == nullptr) {
309 return -std::numeric_limits<double>::max();
310 } else {
311 return lower_bounds_[index];
312 }
313 }
314
315 double UpperBoundForParameter(int index) const {
316 if (upper_bounds_.get() == nullptr) {
317 return std::numeric_limits<double>::max();
318 } else {
319 return upper_bounds_[index];
320 }
321 }
322
323 private:
324 bool UpdateLocalParameterizationJacobian() {
325 if (local_parameterization_ == nullptr) {
326 return true;
327 }
328
329 // Update the local to global Jacobian. In some cases this is
330 // wasted effort; if this is a bottleneck, we will find a solution
331 // at that time.
332
333 const int jacobian_size = Size() * LocalSize();
334 InvalidateArray(jacobian_size, local_parameterization_jacobian_.get());
335 if (!local_parameterization_->ComputeJacobian(
336 state_, local_parameterization_jacobian_.get())) {
337 LOG(WARNING) << "Local parameterization Jacobian computation failed"
338 "for x: "
339 << ConstVectorRef(state_, Size()).transpose();
340 return false;
341 }
342
343 if (!IsArrayValid(jacobian_size, local_parameterization_jacobian_.get())) {
344 LOG(WARNING) << "Local parameterization Jacobian computation returned"
345 << "an invalid matrix for x: "
346 << ConstVectorRef(state_, Size()).transpose()
347 << "\n Jacobian matrix : "
348 << ConstMatrixRef(local_parameterization_jacobian_.get(),
349 Size(),
350 LocalSize());
351 return false;
352 }
353 return true;
354 }
355
356 double* user_state_ = nullptr;
357 int size_ = -1;
358 bool is_set_constant_ = false;
359 LocalParameterization* local_parameterization_ = nullptr;
360
361 // The "state" of the parameter. These fields are only needed while the
362 // solver is running. While at first glance using mutable is a bad idea, this
363 // ends up simplifying the internals of Ceres enough to justify the potential
364 // pitfalls of using "mutable."
365 mutable const double* state_ = nullptr;
366 mutable std::unique_ptr<double[]> local_parameterization_jacobian_;
367
368 // The index of the parameter. This is used by various other parts of Ceres to
369 // permit switching from a ParameterBlock* to an index in another array.
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800370 int index_ = -1;
Austin Schuh70cc9552019-01-21 19:46:48 -0800371
372 // The offset of this parameter block inside a larger state vector.
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800373 int state_offset_ = -1;
Austin Schuh70cc9552019-01-21 19:46:48 -0800374
375 // The offset of this parameter block inside a larger delta vector.
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800376 int delta_offset_ = -1;
Austin Schuh70cc9552019-01-21 19:46:48 -0800377
378 // If non-null, contains the residual blocks this parameter block is in.
379 std::unique_ptr<ResidualBlockSet> residual_blocks_;
380
381 // Upper and lower bounds for the parameter block. SetUpperBound
382 // and SetLowerBound lazily initialize the upper_bounds_ and
383 // lower_bounds_ arrays. If they are never called, then memory for
384 // these arrays is never allocated. Thus for problems where there
385 // are no bounds, or only one sided bounds we do not pay the cost of
386 // allocating memory for the inactive bounds constraints.
387 //
388 // Upon initialization these arrays are initialized to
389 // std::numeric_limits<double>::max() and
390 // -std::numeric_limits<double>::max() respectively which correspond
391 // to the parameter block being unconstrained.
392 std::unique_ptr<double[]> upper_bounds_;
393 std::unique_ptr<double[]> lower_bounds_;
394
395 // Necessary so ProblemImpl can clean up the parameterizations.
396 friend class ProblemImpl;
397};
398
399} // namespace internal
400} // namespace ceres
401
402#endif // CERES_INTERNAL_PARAMETER_BLOCK_H_