blob: 925d1c42cd24339536bea2726a1152bf48880070 [file] [log] [blame]
Austin Schuh70cc9552019-01-21 19:46:48 -08001// Ceres Solver - A fast non-linear least squares minimizer
Austin Schuh3de38b02024-06-25 18:25:10 -07002// Copyright 2023 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"
Austin Schuh3de38b02024-06-25 18:25:10 -070043#include "ceres/internal/disable_warnings.h"
Austin Schuh70cc9552019-01-21 19:46:48 -080044#include "ceres/internal/eigen.h"
Austin Schuh3de38b02024-06-25 18:25:10 -070045#include "ceres/internal/export.h"
46#include "ceres/manifold.h"
Austin Schuh70cc9552019-01-21 19:46:48 -080047#include "ceres/stringprintf.h"
48#include "glog/logging.h"
49
Austin Schuh3de38b02024-06-25 18:25:10 -070050namespace ceres::internal {
Austin Schuh70cc9552019-01-21 19:46:48 -080051
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
Austin Schuh3de38b02024-06-25 18:25:10 -070061// parameter block may also hold a pointer to a manifold; the parameter block
62// does not take ownership of this pointer, so the user is responsible for the
63// proper disposal of the manifold.
64class CERES_NO_EXPORT ParameterBlock {
Austin Schuh70cc9552019-01-21 19:46:48 -080065 public:
Austin Schuh3de38b02024-06-25 18:25:10 -070066 using ResidualBlockSet = std::unordered_set<ResidualBlock*>;
Austin Schuh70cc9552019-01-21 19:46:48 -080067
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
Austin Schuh3de38b02024-06-25 18:25:10 -070077 ParameterBlock(double* user_state, int size, int index, Manifold* manifold)
Austin Schuh70cc9552019-01-21 19:46:48 -080078 : user_state_(user_state),
Austin Schuh70cc9552019-01-21 19:46:48 -080079 size_(size),
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080080 state_(user_state),
Austin Schuh70cc9552019-01-21 19:46:48 -080081 index_(index) {
Austin Schuh3de38b02024-06-25 18:25:10 -070082 if (manifold != nullptr) {
83 SetManifold(manifold);
Austin Schuh70cc9552019-01-21 19:46:48 -080084 }
85 }
86
87 // The size of the parameter block.
88 int Size() const { return size_; }
89
90 // Manipulate the parameter state.
91 bool SetState(const double* x) {
92 CHECK(x != nullptr) << "Tried to set the state of constant parameter "
93 << "with user location " << user_state_;
94 CHECK(!IsConstant()) << "Tried to set the state of constant parameter "
95 << "with user location " << user_state_;
96
97 state_ = x;
Austin Schuh3de38b02024-06-25 18:25:10 -070098 return UpdatePlusJacobian();
Austin Schuh70cc9552019-01-21 19:46:48 -080099 }
100
101 // Copy the current parameter state out to x. This is "GetState()" rather than
102 // simply "state()" since it is actively copying the data into the passed
103 // pointer.
104 void GetState(double* x) const {
105 if (x != state_) {
106 std::copy(state_, state_ + size_, x);
107 }
108 }
109
110 // Direct pointers to the current state.
111 const double* state() const { return state_; }
112 const double* user_state() const { return user_state_; }
113 double* mutable_user_state() { return user_state_; }
Austin Schuh3de38b02024-06-25 18:25:10 -0700114 const Manifold* manifold() const { return manifold_; }
115 Manifold* mutable_manifold() { return manifold_; }
Austin Schuh70cc9552019-01-21 19:46:48 -0800116
117 // Set this parameter block to vary or not.
118 void SetConstant() { is_set_constant_ = true; }
119 void SetVarying() { is_set_constant_ = false; }
Austin Schuh3de38b02024-06-25 18:25:10 -0700120 bool IsConstant() const { return (is_set_constant_ || TangentSize() == 0); }
Austin Schuh70cc9552019-01-21 19:46:48 -0800121
122 double UpperBound(int index) const {
123 return (upper_bounds_ ? upper_bounds_[index]
124 : std::numeric_limits<double>::max());
125 }
126
127 double LowerBound(int index) const {
128 return (lower_bounds_ ? lower_bounds_[index]
129 : -std::numeric_limits<double>::max());
130 }
131
132 bool IsUpperBounded() const { return (upper_bounds_ == nullptr); }
133 bool IsLowerBounded() const { return (lower_bounds_ == nullptr); }
134
135 // This parameter block's index in an array.
136 int index() const { return index_; }
137 void set_index(int index) { index_ = index; }
138
139 // This parameter offset inside a larger state vector.
140 int state_offset() const { return state_offset_; }
141 void set_state_offset(int state_offset) { state_offset_ = state_offset; }
142
143 // This parameter offset inside a larger delta vector.
144 int delta_offset() const { return delta_offset_; }
145 void set_delta_offset(int delta_offset) { delta_offset_ = delta_offset; }
146
Austin Schuh3de38b02024-06-25 18:25:10 -0700147 // Methods relating to the parameter block's manifold.
Austin Schuh70cc9552019-01-21 19:46:48 -0800148
Austin Schuh3de38b02024-06-25 18:25:10 -0700149 // The local to global jacobian. Returns nullptr if there is no manifold for
150 // this parameter block. The returned matrix is row-major and has Size() rows
151 // and TangentSize() columns.
152 const double* PlusJacobian() const { return plus_jacobian_.get(); }
153
154 int TangentSize() const {
155 return (manifold_ == nullptr) ? size_ : manifold_->TangentSize();
Austin Schuh70cc9552019-01-21 19:46:48 -0800156 }
157
Austin Schuh3de38b02024-06-25 18:25:10 -0700158 // Set the manifold. The parameter block does not take ownership of
159 // the manifold.
160 void SetManifold(Manifold* new_manifold) {
161 // Nothing to do if the new manifold is the same as the old
162 // manifold.
163 if (new_manifold == manifold_) {
Austin Schuh70cc9552019-01-21 19:46:48 -0800164 return;
165 }
166
Austin Schuh3de38b02024-06-25 18:25:10 -0700167 if (new_manifold == nullptr) {
168 manifold_ = nullptr;
169 plus_jacobian_ = nullptr;
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800170 return;
171 }
Austin Schuh70cc9552019-01-21 19:46:48 -0800172
Austin Schuh3de38b02024-06-25 18:25:10 -0700173 CHECK_EQ(new_manifold->AmbientSize(), size_)
174 << "The parameter block has size = " << size_
175 << " while the manifold has ambient size = "
176 << new_manifold->AmbientSize();
Austin Schuh70cc9552019-01-21 19:46:48 -0800177
Austin Schuh3de38b02024-06-25 18:25:10 -0700178 CHECK_GE(new_manifold->TangentSize(), 0)
179 << "Invalid Manifold. Manifolds must have a "
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800180 << "non-negative dimensional tangent space.";
Austin Schuh70cc9552019-01-21 19:46:48 -0800181
Austin Schuh3de38b02024-06-25 18:25:10 -0700182 manifold_ = new_manifold;
183 plus_jacobian_ = std::make_unique<double[]>(manifold_->AmbientSize() *
184 manifold_->TangentSize());
185 CHECK(UpdatePlusJacobian())
186 << "Manifold::PlusJacobian computation failed for x: "
Austin Schuh70cc9552019-01-21 19:46:48 -0800187 << ConstVectorRef(state_, Size()).transpose();
188 }
189
190 void SetUpperBound(int index, double upper_bound) {
191 CHECK_LT(index, size_);
192
193 if (upper_bound >= std::numeric_limits<double>::max() && !upper_bounds_) {
194 return;
195 }
196
197 if (!upper_bounds_) {
Austin Schuh3de38b02024-06-25 18:25:10 -0700198 upper_bounds_ = std::make_unique<double[]>(size_);
Austin Schuh70cc9552019-01-21 19:46:48 -0800199 std::fill(upper_bounds_.get(),
200 upper_bounds_.get() + size_,
201 std::numeric_limits<double>::max());
202 }
203
204 upper_bounds_[index] = upper_bound;
205 }
206
207 void SetLowerBound(int index, double lower_bound) {
208 CHECK_LT(index, size_);
209
210 if (lower_bound <= -std::numeric_limits<double>::max() && !lower_bounds_) {
211 return;
212 }
213
214 if (!lower_bounds_) {
Austin Schuh3de38b02024-06-25 18:25:10 -0700215 lower_bounds_ = std::make_unique<double[]>(size_);
Austin Schuh70cc9552019-01-21 19:46:48 -0800216 std::fill(lower_bounds_.get(),
217 lower_bounds_.get() + size_,
218 -std::numeric_limits<double>::max());
219 }
220
221 lower_bounds_[index] = lower_bound;
222 }
223
224 // Generalization of the addition operation. This is the same as
Austin Schuh3de38b02024-06-25 18:25:10 -0700225 // Manifold::Plus() followed by projection onto the
Austin Schuh70cc9552019-01-21 19:46:48 -0800226 // hyper cube implied by the bounds constraints.
227 bool Plus(const double* x, const double* delta, double* x_plus_delta) {
Austin Schuh3de38b02024-06-25 18:25:10 -0700228 if (manifold_ != nullptr) {
229 if (!manifold_->Plus(x, delta, x_plus_delta)) {
Austin Schuh70cc9552019-01-21 19:46:48 -0800230 return false;
231 }
232 } else {
233 VectorRef(x_plus_delta, size_) =
234 ConstVectorRef(x, size_) + ConstVectorRef(delta, size_);
235 }
236
237 // Project onto the box constraints.
238 if (lower_bounds_.get() != nullptr) {
239 for (int i = 0; i < size_; ++i) {
240 x_plus_delta[i] = std::max(x_plus_delta[i], lower_bounds_[i]);
241 }
242 }
243
244 if (upper_bounds_.get() != nullptr) {
245 for (int i = 0; i < size_; ++i) {
246 x_plus_delta[i] = std::min(x_plus_delta[i], upper_bounds_[i]);
247 }
248 }
249
250 return true;
251 }
252
253 std::string ToString() const {
254 return StringPrintf(
255 "{ this=%p, user_state=%p, state=%p, size=%d, "
256 "constant=%d, index=%d, state_offset=%d, "
257 "delta_offset=%d }",
258 this,
259 user_state_,
260 state_,
261 size_,
262 is_set_constant_,
263 index_,
264 state_offset_,
265 delta_offset_);
266 }
267
268 void EnableResidualBlockDependencies() {
269 CHECK(residual_blocks_.get() == nullptr)
270 << "Ceres bug: There is already a residual block collection "
271 << "for parameter block: " << ToString();
Austin Schuh3de38b02024-06-25 18:25:10 -0700272 residual_blocks_ = std::make_unique<ResidualBlockSet>();
Austin Schuh70cc9552019-01-21 19:46:48 -0800273 }
274
275 void AddResidualBlock(ResidualBlock* residual_block) {
276 CHECK(residual_blocks_.get() != nullptr)
277 << "Ceres bug: The residual block collection is null for parameter "
278 << "block: " << ToString();
279 residual_blocks_->insert(residual_block);
280 }
281
282 void RemoveResidualBlock(ResidualBlock* residual_block) {
283 CHECK(residual_blocks_.get() != nullptr)
284 << "Ceres bug: The residual block collection is null for parameter "
285 << "block: " << ToString();
286 CHECK(residual_blocks_->find(residual_block) != residual_blocks_->end())
287 << "Ceres bug: Missing residual for parameter block: " << ToString();
288 residual_blocks_->erase(residual_block);
289 }
290
291 // This is only intended for iterating; perhaps this should only expose
292 // .begin() and .end().
293 ResidualBlockSet* mutable_residual_blocks() { return residual_blocks_.get(); }
294
295 double LowerBoundForParameter(int index) const {
296 if (lower_bounds_.get() == nullptr) {
297 return -std::numeric_limits<double>::max();
298 } else {
299 return lower_bounds_[index];
300 }
301 }
302
303 double UpperBoundForParameter(int index) const {
304 if (upper_bounds_.get() == nullptr) {
305 return std::numeric_limits<double>::max();
306 } else {
307 return upper_bounds_[index];
308 }
309 }
310
311 private:
Austin Schuh3de38b02024-06-25 18:25:10 -0700312 bool UpdatePlusJacobian() {
313 if (manifold_ == nullptr) {
Austin Schuh70cc9552019-01-21 19:46:48 -0800314 return true;
315 }
316
Austin Schuh3de38b02024-06-25 18:25:10 -0700317 // Update the Plus Jacobian. In some cases this is
Austin Schuh70cc9552019-01-21 19:46:48 -0800318 // wasted effort; if this is a bottleneck, we will find a solution
319 // at that time.
Austin Schuh3de38b02024-06-25 18:25:10 -0700320 const int jacobian_size = Size() * TangentSize();
321 InvalidateArray(jacobian_size, plus_jacobian_.get());
322 if (!manifold_->PlusJacobian(state_, plus_jacobian_.get())) {
323 LOG(WARNING) << "Manifold::PlusJacobian computation failed"
Austin Schuh70cc9552019-01-21 19:46:48 -0800324 "for x: "
325 << ConstVectorRef(state_, Size()).transpose();
326 return false;
327 }
328
Austin Schuh3de38b02024-06-25 18:25:10 -0700329 if (!IsArrayValid(jacobian_size, plus_jacobian_.get())) {
330 LOG(WARNING) << "Manifold::PlusJacobian computation returned "
Austin Schuh70cc9552019-01-21 19:46:48 -0800331 << "an invalid matrix for x: "
332 << ConstVectorRef(state_, Size()).transpose()
333 << "\n Jacobian matrix : "
Austin Schuh3de38b02024-06-25 18:25:10 -0700334 << ConstMatrixRef(
335 plus_jacobian_.get(), Size(), TangentSize());
Austin Schuh70cc9552019-01-21 19:46:48 -0800336 return false;
337 }
338 return true;
339 }
340
341 double* user_state_ = nullptr;
342 int size_ = -1;
343 bool is_set_constant_ = false;
Austin Schuh3de38b02024-06-25 18:25:10 -0700344 Manifold* manifold_ = nullptr;
Austin Schuh70cc9552019-01-21 19:46:48 -0800345
346 // The "state" of the parameter. These fields are only needed while the
347 // solver is running. While at first glance using mutable is a bad idea, this
348 // ends up simplifying the internals of Ceres enough to justify the potential
349 // pitfalls of using "mutable."
350 mutable const double* state_ = nullptr;
Austin Schuh3de38b02024-06-25 18:25:10 -0700351 mutable std::unique_ptr<double[]> plus_jacobian_;
Austin Schuh70cc9552019-01-21 19:46:48 -0800352
353 // The index of the parameter. This is used by various other parts of Ceres to
354 // permit switching from a ParameterBlock* to an index in another array.
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800355 int index_ = -1;
Austin Schuh70cc9552019-01-21 19:46:48 -0800356
357 // The offset of this parameter block inside a larger state vector.
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800358 int state_offset_ = -1;
Austin Schuh70cc9552019-01-21 19:46:48 -0800359
360 // The offset of this parameter block inside a larger delta vector.
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800361 int delta_offset_ = -1;
Austin Schuh70cc9552019-01-21 19:46:48 -0800362
363 // If non-null, contains the residual blocks this parameter block is in.
364 std::unique_ptr<ResidualBlockSet> residual_blocks_;
365
366 // Upper and lower bounds for the parameter block. SetUpperBound
367 // and SetLowerBound lazily initialize the upper_bounds_ and
368 // lower_bounds_ arrays. If they are never called, then memory for
369 // these arrays is never allocated. Thus for problems where there
370 // are no bounds, or only one sided bounds we do not pay the cost of
371 // allocating memory for the inactive bounds constraints.
372 //
373 // Upon initialization these arrays are initialized to
374 // std::numeric_limits<double>::max() and
375 // -std::numeric_limits<double>::max() respectively which correspond
376 // to the parameter block being unconstrained.
377 std::unique_ptr<double[]> upper_bounds_;
378 std::unique_ptr<double[]> lower_bounds_;
379
Austin Schuh3de38b02024-06-25 18:25:10 -0700380 // Necessary so ProblemImpl can clean up the manifolds.
Austin Schuh70cc9552019-01-21 19:46:48 -0800381 friend class ProblemImpl;
382};
383
Austin Schuh3de38b02024-06-25 18:25:10 -0700384} // namespace ceres::internal
385
386#include "ceres/internal/reenable_warnings.h"
Austin Schuh70cc9552019-01-21 19:46:48 -0800387
388#endif // CERES_INTERNAL_PARAMETER_BLOCK_H_