1// Ceres Solver - A fast non-linear least squares minimizer
2// Copyright 2010, 2011, 2012, 2013 Google Inc. All rights reserved.
3// http://code.google.com/p/ceres-solver/
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 <cstdlib>
36#include <limits>
37#include <string>
38#include "ceres/array_utils.h"
39#include "ceres/collections_port.h"
40#include "ceres/integral_types.h"
41#include "ceres/internal/eigen.h"
42#include "ceres/internal/port.h"
43#include "ceres/internal/scoped_ptr.h"
44#include "ceres/local_parameterization.h"
45#include "ceres/stringprintf.h"
46#include "glog/logging.h"
47
48namespace ceres {
49namespace internal {
50
51class ProblemImpl;
52class ResidualBlock;
53
54// The parameter block encodes the location of the user's original value, and
55// also the "current state" of the parameter. The evaluator uses whatever is in
56// the current state of the parameter when evaluating. This is inlined since the
57// methods are performance sensitive.
58//
59// The class is not thread-safe, unless only const methods are called. The
60// parameter block may also hold a pointer to a local parameterization; the
61// parameter block does not take ownership of this pointer, so the user is
62// responsible for the proper disposal of the local parameterization.
63class ParameterBlock {
64 public:
65  // TODO(keir): Decide what data structure is best here. Should this be a set?
66  // Probably not, because sets are memory inefficient. However, if it's a
67  // vector, you can get into pathological linear performance when removing a
68  // residual block from a problem where all the residual blocks depend on one
69  // parameter; for example, shared focal length in a bundle adjustment
70  // problem. It might be worth making a custom structure that is just an array
71  // when it is small, but transitions to a hash set when it has more elements.
72  //
73  // For now, use a hash set.
74  typedef HashSet<ResidualBlock*> ResidualBlockSet;
75
76  // Create a parameter block with the user state, size, and index specified.
77  // The size is the size of the parameter block and the index is the position
78  // of the parameter block inside a Program (if any).
79  ParameterBlock(double* user_state, int size, int index) {
80    Init(user_state, size, index, NULL);
81  }
82
83  ParameterBlock(double* user_state,
84                 int size,
85                 int index,
86                 LocalParameterization* local_parameterization) {
87    Init(user_state, size, index, local_parameterization);
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 != NULL)
96        << "Tried to set the state of constant parameter "
97        << "with user location " << user_state_;
98    CHECK(!is_constant_)
99        << "Tried to set the state of constant parameter "
100        << "with user location " << user_state_;
101
102    state_ = x;
103    return UpdateLocalParameterizationJacobian();
104  }
105
106  // Copy the current parameter state out to x. This is "GetState()" rather than
107  // simply "state()" since it is actively copying the data into the passed
108  // pointer.
109  void GetState(double *x) const {
110    if (x != state_) {
111      memcpy(x, state_, sizeof(*state_) * size_);
112    }
113  }
114
115  // Direct pointers to the current state.
116  const double* state() const { return state_; }
117  const double* user_state() const { return user_state_; }
118  double* mutable_user_state() { return user_state_; }
119  LocalParameterization* local_parameterization() const {
120    return local_parameterization_;
121  }
122  LocalParameterization* mutable_local_parameterization() {
123    return local_parameterization_;
124  }
125
126  // Set this parameter block to vary or not.
127  void SetConstant() { is_constant_ = true; }
128  void SetVarying() { is_constant_ = false; }
129  bool IsConstant() const { return is_constant_; }
130
131  // This parameter block's index in an array.
132  int index() const { return index_; }
133  void set_index(int index) { index_ = index; }
134
135  // This parameter offset inside a larger state vector.
136  int state_offset() const { return state_offset_; }
137  void set_state_offset(int state_offset) { state_offset_ = state_offset; }
138
139  // This parameter offset inside a larger delta vector.
140  int delta_offset() const { return delta_offset_; }
141  void set_delta_offset(int delta_offset) { delta_offset_ = delta_offset; }
142
143  // Methods relating to the parameter block's parameterization.
144
145  // The local to global jacobian. Returns NULL if there is no local
146  // parameterization for this parameter block. The returned matrix is row-major
147  // and has Size() rows and  LocalSize() columns.
148  const double* LocalParameterizationJacobian() const {
149    return local_parameterization_jacobian_.get();
150  }
151
152  int LocalSize() const {
153    return (local_parameterization_ == NULL)
154        ? size_
155        : local_parameterization_->LocalSize();
156  }
157
158  // Set the parameterization. The parameterization can be set exactly once;
159  // multiple calls to set the parameterization to different values will crash.
160  // It is an error to pass NULL for the parameterization. The parameter block
161  // does not take ownership of the parameterization.
162  void SetParameterization(LocalParameterization* new_parameterization) {
163    CHECK(new_parameterization != NULL) << "NULL parameterization invalid.";
164    CHECK(new_parameterization->GlobalSize() == size_)
165        << "Invalid parameterization for parameter block. The parameter block "
166        << "has size " << size_ << " while the parameterization has a global "
167        << "size of " << new_parameterization->GlobalSize() << ". Did you "
168        << "accidentally use the wrong parameter block or parameterization?";
169    if (new_parameterization != local_parameterization_) {
170      CHECK(local_parameterization_ == NULL)
171          << "Can't re-set the local parameterization; it leads to "
172          << "ambiguous ownership.";
173      local_parameterization_ = new_parameterization;
174      local_parameterization_jacobian_.reset(
175          new double[local_parameterization_->GlobalSize() *
176                     local_parameterization_->LocalSize()]);
177      CHECK(UpdateLocalParameterizationJacobian())
178          << "Local parameterization Jacobian computation failed for x: "
179          << ConstVectorRef(state_, Size()).transpose();
180    } else {
181      // Ignore the case that the parameterizations match.
182    }
183  }
184
185  void SetUpperBound(int index, double upper_bound) {
186    CHECK_LT(index, size_);
187
188    if (upper_bounds_.get() == NULL) {
189      upper_bounds_.reset(new double[size_]);
190      std::fill(upper_bounds_.get(),
191                upper_bounds_.get() + size_,
192                std::numeric_limits<double>::max());
193    }
194
195    upper_bounds_[index] = upper_bound;
196  };
197
198  void SetLowerBound(int index, double lower_bound) {
199    CHECK_LT(index, size_);
200
201    if (lower_bounds_.get() == NULL) {
202      lower_bounds_.reset(new double[size_]);
203      std::fill(lower_bounds_.get(),
204                lower_bounds_.get() + size_,
205                -std::numeric_limits<double>::max());
206    }
207
208    lower_bounds_[index] = lower_bound;
209  }
210
211  // Generalization of the addition operation. This is the same as
212  // LocalParameterization::Plus() followed by projection onto the
213  // hyper cube implied by the bounds constraints.
214  bool Plus(const double *x, const double* delta, double* x_plus_delta) {
215    if (local_parameterization_ != NULL) {
216      if (!local_parameterization_->Plus(x, delta, x_plus_delta)) {
217        return false;
218      }
219    } else {
220      VectorRef(x_plus_delta, size_) = ConstVectorRef(x, size_) +
221                                       ConstVectorRef(delta,  size_);
222    }
223
224    // Project onto the box constraints.
225    if (lower_bounds_.get() != NULL) {
226      for (int i = 0; i < size_; ++i) {
227        x_plus_delta[i] = std::max(x_plus_delta[i], lower_bounds_[i]);
228      }
229    }
230
231    if (upper_bounds_.get() != NULL) {
232      for (int i = 0; i < size_; ++i) {
233        x_plus_delta[i] = std::min(x_plus_delta[i], upper_bounds_[i]);
234      }
235    }
236
237    return true;
238  }
239
240  string ToString() const {
241    return StringPrintf("{ user_state=%p, state=%p, size=%d, "
242                        "constant=%d, index=%d, state_offset=%d, "
243                        "delta_offset=%d }",
244                        user_state_,
245                        state_,
246                        size_,
247                        is_constant_,
248                        index_,
249                        state_offset_,
250                        delta_offset_);
251  }
252
253  void EnableResidualBlockDependencies() {
254    CHECK(residual_blocks_.get() == NULL)
255        << "Ceres bug: There is already a residual block collection "
256        << "for parameter block: " << ToString();
257    residual_blocks_.reset(new ResidualBlockSet);
258  }
259
260  void AddResidualBlock(ResidualBlock* residual_block) {
261    CHECK(residual_blocks_.get() != NULL)
262        << "Ceres bug: The residual block collection is null for parameter "
263        << "block: " << ToString();
264    residual_blocks_->insert(residual_block);
265  }
266
267  void RemoveResidualBlock(ResidualBlock* residual_block) {
268    CHECK(residual_blocks_.get() != NULL)
269        << "Ceres bug: The residual block collection is null for parameter "
270        << "block: " << ToString();
271    CHECK(residual_blocks_->find(residual_block) != residual_blocks_->end())
272        << "Ceres bug: Missing residual for parameter block: " << ToString();
273    residual_blocks_->erase(residual_block);
274  }
275
276  // This is only intended for iterating; perhaps this should only expose
277  // .begin() and .end().
278  ResidualBlockSet* mutable_residual_blocks() {
279    return residual_blocks_.get();
280  }
281
282  double LowerBoundForParameter(int index) const {
283    if (lower_bounds_.get() == NULL) {
284      return -std::numeric_limits<double>::max();
285    } else {
286      return lower_bounds_[index];
287    }
288  }
289
290  double UpperBoundForParameter(int index) const {
291    if (upper_bounds_.get() == NULL) {
292      return std::numeric_limits<double>::max();
293    } else {
294      return upper_bounds_[index];
295    }
296  }
297
298 private:
299  void Init(double* user_state,
300            int size,
301            int index,
302            LocalParameterization* local_parameterization) {
303    user_state_ = user_state;
304    size_ = size;
305    index_ = index;
306    is_constant_ = false;
307    state_ = user_state_;
308
309    local_parameterization_ = NULL;
310    if (local_parameterization != NULL) {
311      SetParameterization(local_parameterization);
312    }
313
314    state_offset_ = -1;
315    delta_offset_ = -1;
316  }
317
318  bool UpdateLocalParameterizationJacobian() {
319    if (local_parameterization_ == NULL) {
320      return true;
321    }
322
323    // Update the local to global Jacobian. In some cases this is
324    // wasted effort; if this is a bottleneck, we will find a solution
325    // at that time.
326
327    const int jacobian_size = Size() * LocalSize();
328    InvalidateArray(jacobian_size,
329                    local_parameterization_jacobian_.get());
330    if (!local_parameterization_->ComputeJacobian(
331            state_,
332            local_parameterization_jacobian_.get())) {
333      LOG(WARNING) << "Local parameterization Jacobian computation failed"
334          "for x: " << ConstVectorRef(state_, Size()).transpose();
335      return false;
336    }
337
338    if (!IsArrayValid(jacobian_size, local_parameterization_jacobian_.get())) {
339      LOG(WARNING) << "Local parameterization Jacobian computation returned"
340                   << "an invalid matrix for x: "
341                   << ConstVectorRef(state_, Size()).transpose()
342                   << "\n Jacobian matrix : "
343                   << ConstMatrixRef(local_parameterization_jacobian_.get(),
344                                     Size(),
345                                     LocalSize());
346      return false;
347    }
348    return true;
349  }
350
351  double* user_state_;
352  int size_;
353  bool is_constant_;
354  LocalParameterization* local_parameterization_;
355
356  // The "state" of the parameter. These fields are only needed while the
357  // solver is running. While at first glance using mutable is a bad idea, this
358  // ends up simplifying the internals of Ceres enough to justify the potential
359  // pitfalls of using "mutable."
360  mutable const double* state_;
361  mutable scoped_array<double> local_parameterization_jacobian_;
362
363  // The index of the parameter. This is used by various other parts of Ceres to
364  // permit switching from a ParameterBlock* to an index in another array.
365  int32 index_;
366
367  // The offset of this parameter block inside a larger state vector.
368  int32 state_offset_;
369
370  // The offset of this parameter block inside a larger delta vector.
371  int32 delta_offset_;
372
373  // If non-null, contains the residual blocks this parameter block is in.
374  scoped_ptr<ResidualBlockSet> residual_blocks_;
375
376  // Upper and lower bounds for the parameter block.  SetUpperBound
377  // and SetLowerBound lazily initialize the upper_bounds_ and
378  // lower_bounds_ arrays. If they are never called, then memory for
379  // these arrays is never allocated. Thus for problems where there
380  // are no bounds, or only one sided bounds we do not pay the cost of
381  // allocating memory for the inactive bounds constraints.
382  //
383  // Upon initialization these arrays are initialized to
384  // std::numeric_limits<double>::max() and
385  // -std::numeric_limits<double>::max() respectively which correspond
386  // to the parameter block being unconstrained.
387  scoped_array<double> upper_bounds_;
388  scoped_array<double> lower_bounds_;
389
390  // Necessary so ProblemImpl can clean up the parameterizations.
391  friend class ProblemImpl;
392};
393
394}  // namespace internal
395}  // namespace ceres
396
397#endif  // CERES_INTERNAL_PARAMETER_BLOCK_H_
398