Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame^] | 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 | // sameeragarwal@google.com (Sameer Agarwal) |
| 31 | |
| 32 | #ifndef CERES_PUBLIC_LOCAL_PARAMETERIZATION_H_ |
| 33 | #define CERES_PUBLIC_LOCAL_PARAMETERIZATION_H_ |
| 34 | |
| 35 | #include <vector> |
| 36 | #include "ceres/internal/port.h" |
| 37 | #include "ceres/internal/disable_warnings.h" |
| 38 | |
| 39 | namespace ceres { |
| 40 | |
| 41 | // Purpose: Sometimes parameter blocks x can overparameterize a problem |
| 42 | // |
| 43 | // min f(x) |
| 44 | // x |
| 45 | // |
| 46 | // In that case it is desirable to choose a parameterization for the |
| 47 | // block itself to remove the null directions of the cost. More |
| 48 | // generally, if x lies on a manifold of a smaller dimension than the |
| 49 | // ambient space that it is embedded in, then it is numerically and |
| 50 | // computationally more effective to optimize it using a |
| 51 | // parameterization that lives in the tangent space of that manifold |
| 52 | // at each point. |
| 53 | // |
| 54 | // For example, a sphere in three dimensions is a 2 dimensional |
| 55 | // manifold, embedded in a three dimensional space. At each point on |
| 56 | // the sphere, the plane tangent to it defines a two dimensional |
| 57 | // tangent space. For a cost function defined on this sphere, given a |
| 58 | // point x, moving in the direction normal to the sphere at that point |
| 59 | // is not useful. Thus a better way to do a local optimization is to |
| 60 | // optimize over two dimensional vector delta in the tangent space at |
| 61 | // that point and then "move" to the point x + delta, where the move |
| 62 | // operation involves projecting back onto the sphere. Doing so |
| 63 | // removes a redundant dimension from the optimization, making it |
| 64 | // numerically more robust and efficient. |
| 65 | // |
| 66 | // More generally we can define a function |
| 67 | // |
| 68 | // x_plus_delta = Plus(x, delta), |
| 69 | // |
| 70 | // where x_plus_delta has the same size as x, and delta is of size |
| 71 | // less than or equal to x. The function Plus, generalizes the |
| 72 | // definition of vector addition. Thus it satisfies the identify |
| 73 | // |
| 74 | // Plus(x, 0) = x, for all x. |
| 75 | // |
| 76 | // A trivial version of Plus is when delta is of the same size as x |
| 77 | // and |
| 78 | // |
| 79 | // Plus(x, delta) = x + delta |
| 80 | // |
| 81 | // A more interesting case if x is two dimensional vector, and the |
| 82 | // user wishes to hold the first coordinate constant. Then, delta is a |
| 83 | // scalar and Plus is defined as |
| 84 | // |
| 85 | // Plus(x, delta) = x + [0] * delta |
| 86 | // [1] |
| 87 | // |
| 88 | // An example that occurs commonly in Structure from Motion problems |
| 89 | // is when camera rotations are parameterized using Quaternion. There, |
| 90 | // it is useful only make updates orthogonal to that 4-vector defining |
| 91 | // the quaternion. One way to do this is to let delta be a 3 |
| 92 | // dimensional vector and define Plus to be |
| 93 | // |
| 94 | // Plus(x, delta) = [cos(|delta|), sin(|delta|) delta / |delta|] * x |
| 95 | // |
| 96 | // The multiplication between the two 4-vectors on the RHS is the |
| 97 | // standard quaternion product. |
| 98 | // |
| 99 | // Given g and a point x, optimizing f can now be restated as |
| 100 | // |
| 101 | // min f(Plus(x, delta)) |
| 102 | // delta |
| 103 | // |
| 104 | // Given a solution delta to this problem, the optimal value is then |
| 105 | // given by |
| 106 | // |
| 107 | // x* = Plus(x, delta) |
| 108 | // |
| 109 | // The class LocalParameterization defines the function Plus and its |
| 110 | // Jacobian which is needed to compute the Jacobian of f w.r.t delta. |
| 111 | class CERES_EXPORT LocalParameterization { |
| 112 | public: |
| 113 | virtual ~LocalParameterization(); |
| 114 | |
| 115 | // Generalization of the addition operation, |
| 116 | // |
| 117 | // x_plus_delta = Plus(x, delta) |
| 118 | // |
| 119 | // with the condition that Plus(x, 0) = x. |
| 120 | virtual bool Plus(const double* x, |
| 121 | const double* delta, |
| 122 | double* x_plus_delta) const = 0; |
| 123 | |
| 124 | // The jacobian of Plus(x, delta) w.r.t delta at delta = 0. |
| 125 | // |
| 126 | // jacobian is a row-major GlobalSize() x LocalSize() matrix. |
| 127 | virtual bool ComputeJacobian(const double* x, double* jacobian) const = 0; |
| 128 | |
| 129 | // local_matrix = global_matrix * jacobian |
| 130 | // |
| 131 | // global_matrix is a num_rows x GlobalSize row major matrix. |
| 132 | // local_matrix is a num_rows x LocalSize row major matrix. |
| 133 | // jacobian(x) is the matrix returned by ComputeJacobian at x. |
| 134 | // |
| 135 | // This is only used by GradientProblem. For most normal uses, it is |
| 136 | // okay to use the default implementation. |
| 137 | virtual bool MultiplyByJacobian(const double* x, |
| 138 | const int num_rows, |
| 139 | const double* global_matrix, |
| 140 | double* local_matrix) const; |
| 141 | |
| 142 | // Size of x. |
| 143 | virtual int GlobalSize() const = 0; |
| 144 | |
| 145 | // Size of delta. |
| 146 | virtual int LocalSize() const = 0; |
| 147 | }; |
| 148 | |
| 149 | // Some basic parameterizations |
| 150 | |
| 151 | // Identity Parameterization: Plus(x, delta) = x + delta |
| 152 | class CERES_EXPORT IdentityParameterization : public LocalParameterization { |
| 153 | public: |
| 154 | explicit IdentityParameterization(int size); |
| 155 | virtual ~IdentityParameterization() {} |
| 156 | virtual bool Plus(const double* x, |
| 157 | const double* delta, |
| 158 | double* x_plus_delta) const; |
| 159 | virtual bool ComputeJacobian(const double* x, |
| 160 | double* jacobian) const; |
| 161 | virtual bool MultiplyByJacobian(const double* x, |
| 162 | const int num_cols, |
| 163 | const double* global_matrix, |
| 164 | double* local_matrix) const; |
| 165 | virtual int GlobalSize() const { return size_; } |
| 166 | virtual int LocalSize() const { return size_; } |
| 167 | |
| 168 | private: |
| 169 | const int size_; |
| 170 | }; |
| 171 | |
| 172 | // Hold a subset of the parameters inside a parameter block constant. |
| 173 | class CERES_EXPORT SubsetParameterization : public LocalParameterization { |
| 174 | public: |
| 175 | explicit SubsetParameterization(int size, |
| 176 | const std::vector<int>& constant_parameters); |
| 177 | virtual ~SubsetParameterization() {} |
| 178 | virtual bool Plus(const double* x, |
| 179 | const double* delta, |
| 180 | double* x_plus_delta) const; |
| 181 | virtual bool ComputeJacobian(const double* x, |
| 182 | double* jacobian) const; |
| 183 | virtual bool MultiplyByJacobian(const double* x, |
| 184 | const int num_cols, |
| 185 | const double* global_matrix, |
| 186 | double* local_matrix) const; |
| 187 | virtual int GlobalSize() const { |
| 188 | return static_cast<int>(constancy_mask_.size()); |
| 189 | } |
| 190 | virtual int LocalSize() const { return local_size_; } |
| 191 | |
| 192 | private: |
| 193 | const int local_size_; |
| 194 | std::vector<char> constancy_mask_; |
| 195 | }; |
| 196 | |
| 197 | // Plus(x, delta) = [cos(|delta|), sin(|delta|) delta / |delta|] * x |
| 198 | // with * being the quaternion multiplication operator. Here we assume |
| 199 | // that the first element of the quaternion vector is the real (cos |
| 200 | // theta) part. |
| 201 | class CERES_EXPORT QuaternionParameterization : public LocalParameterization { |
| 202 | public: |
| 203 | virtual ~QuaternionParameterization() {} |
| 204 | virtual bool Plus(const double* x, |
| 205 | const double* delta, |
| 206 | double* x_plus_delta) const; |
| 207 | virtual bool ComputeJacobian(const double* x, |
| 208 | double* jacobian) const; |
| 209 | virtual int GlobalSize() const { return 4; } |
| 210 | virtual int LocalSize() const { return 3; } |
| 211 | }; |
| 212 | |
| 213 | // Implements the quaternion local parameterization for Eigen's representation |
| 214 | // of the quaternion. Eigen uses a different internal memory layout for the |
| 215 | // elements of the quaternion than what is commonly used. Specifically, Eigen |
| 216 | // stores the elements in memory as [x, y, z, w] where the real part is last |
| 217 | // whereas it is typically stored first. Note, when creating an Eigen quaternion |
| 218 | // through the constructor the elements are accepted in w, x, y, z order. Since |
| 219 | // Ceres operates on parameter blocks which are raw double pointers this |
| 220 | // difference is important and requires a different parameterization. |
| 221 | // |
| 222 | // Plus(x, delta) = [sin(|delta|) delta / |delta|, cos(|delta|)] * x |
| 223 | // with * being the quaternion multiplication operator. |
| 224 | class CERES_EXPORT EigenQuaternionParameterization |
| 225 | : public ceres::LocalParameterization { |
| 226 | public: |
| 227 | virtual ~EigenQuaternionParameterization() {} |
| 228 | virtual bool Plus(const double* x, |
| 229 | const double* delta, |
| 230 | double* x_plus_delta) const; |
| 231 | virtual bool ComputeJacobian(const double* x, double* jacobian) const; |
| 232 | virtual int GlobalSize() const { return 4; } |
| 233 | virtual int LocalSize() const { return 3; } |
| 234 | }; |
| 235 | |
| 236 | // This provides a parameterization for homogeneous vectors which are commonly |
| 237 | // used in Structure for Motion problems. One example where they are used is |
| 238 | // in representing points whose triangulation is ill-conditioned. Here |
| 239 | // it is advantageous to use an over-parameterization since homogeneous vectors |
| 240 | // can represent points at infinity. |
| 241 | // |
| 242 | // The plus operator is defined as |
| 243 | // Plus(x, delta) = |
| 244 | // [sin(0.5 * |delta|) * delta / |delta|, cos(0.5 * |delta|)] * x |
| 245 | // with * defined as an operator which applies the update orthogonal to x to |
| 246 | // remain on the sphere. We assume that the last element of x is the scalar |
| 247 | // component. The size of the homogeneous vector is required to be greater than |
| 248 | // 1. |
| 249 | class CERES_EXPORT HomogeneousVectorParameterization : |
| 250 | public LocalParameterization { |
| 251 | public: |
| 252 | explicit HomogeneousVectorParameterization(int size); |
| 253 | virtual ~HomogeneousVectorParameterization() {} |
| 254 | virtual bool Plus(const double* x, |
| 255 | const double* delta, |
| 256 | double* x_plus_delta) const; |
| 257 | virtual bool ComputeJacobian(const double* x, |
| 258 | double* jacobian) const; |
| 259 | virtual int GlobalSize() const { return size_; } |
| 260 | virtual int LocalSize() const { return size_ - 1; } |
| 261 | |
| 262 | private: |
| 263 | const int size_; |
| 264 | }; |
| 265 | |
| 266 | // Construct a local parameterization by taking the Cartesian product |
| 267 | // of a number of other local parameterizations. This is useful, when |
| 268 | // a parameter block is the cartesian product of two or more |
| 269 | // manifolds. For example the parameters of a camera consist of a |
| 270 | // rotation and a translation, i.e., SO(3) x R^3. |
| 271 | // |
| 272 | // Currently this class supports taking the cartesian product of up to |
| 273 | // four local parameterizations. |
| 274 | // |
| 275 | // Example usage: |
| 276 | // |
| 277 | // ProductParameterization product_param(new QuaterionionParameterization(), |
| 278 | // new IdentityParameterization(3)); |
| 279 | // |
| 280 | // is the local parameterization for a rigid transformation, where the |
| 281 | // rotation is represented using a quaternion. |
| 282 | class CERES_EXPORT ProductParameterization : public LocalParameterization { |
| 283 | public: |
| 284 | // |
| 285 | // NOTE: All the constructors take ownership of the input local |
| 286 | // parameterizations. |
| 287 | // |
| 288 | ProductParameterization(LocalParameterization* local_param1, |
| 289 | LocalParameterization* local_param2); |
| 290 | |
| 291 | ProductParameterization(LocalParameterization* local_param1, |
| 292 | LocalParameterization* local_param2, |
| 293 | LocalParameterization* local_param3); |
| 294 | |
| 295 | ProductParameterization(LocalParameterization* local_param1, |
| 296 | LocalParameterization* local_param2, |
| 297 | LocalParameterization* local_param3, |
| 298 | LocalParameterization* local_param4); |
| 299 | |
| 300 | virtual ~ProductParameterization(); |
| 301 | virtual bool Plus(const double* x, |
| 302 | const double* delta, |
| 303 | double* x_plus_delta) const; |
| 304 | virtual bool ComputeJacobian(const double* x, |
| 305 | double* jacobian) const; |
| 306 | virtual int GlobalSize() const { return global_size_; } |
| 307 | virtual int LocalSize() const { return local_size_; } |
| 308 | |
| 309 | private: |
| 310 | void Init(); |
| 311 | |
| 312 | std::vector<LocalParameterization*> local_params_; |
| 313 | int local_size_; |
| 314 | int global_size_; |
| 315 | int buffer_size_; |
| 316 | }; |
| 317 | |
| 318 | } // namespace ceres |
| 319 | |
| 320 | #include "ceres/internal/reenable_warnings.h" |
| 321 | |
| 322 | #endif // CERES_PUBLIC_LOCAL_PARAMETERIZATION_H_ |