Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame^] | 1 | // Ceres Solver - A fast non-linear least squares minimizer |
| 2 | // Copyright 2016 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: wjr@google.com (William Rucklidge) |
| 30 | // |
| 31 | // This file contains tests for the GradientChecker class. |
| 32 | |
| 33 | #include "ceres/gradient_checker.h" |
| 34 | |
| 35 | #include <cmath> |
| 36 | #include <cstdlib> |
| 37 | #include <vector> |
| 38 | |
| 39 | #include "ceres/cost_function.h" |
| 40 | #include "ceres/problem.h" |
| 41 | #include "ceres/random.h" |
| 42 | #include "ceres/solver.h" |
| 43 | #include "ceres/test_util.h" |
| 44 | #include "glog/logging.h" |
| 45 | #include "gtest/gtest.h" |
| 46 | |
| 47 | namespace ceres { |
| 48 | namespace internal { |
| 49 | |
| 50 | using std::vector; |
| 51 | |
| 52 | // We pick a (non-quadratic) function whose derivative are easy: |
| 53 | // |
| 54 | // f = exp(- a' x). |
| 55 | // df = - f a. |
| 56 | // |
| 57 | // where 'a' is a vector of the same size as 'x'. In the block |
| 58 | // version, they are both block vectors, of course. |
| 59 | class GoodTestTerm : public CostFunction { |
| 60 | public: |
| 61 | GoodTestTerm(int arity, int const* dim) : arity_(arity), return_value_(true) { |
| 62 | // Make 'arity' random vectors. |
| 63 | a_.resize(arity_); |
| 64 | for (int j = 0; j < arity_; ++j) { |
| 65 | a_[j].resize(dim[j]); |
| 66 | for (int u = 0; u < dim[j]; ++u) { |
| 67 | a_[j][u] = 2.0 * RandDouble() - 1.0; |
| 68 | } |
| 69 | } |
| 70 | |
| 71 | for (int i = 0; i < arity_; i++) { |
| 72 | mutable_parameter_block_sizes()->push_back(dim[i]); |
| 73 | } |
| 74 | set_num_residuals(1); |
| 75 | } |
| 76 | |
| 77 | bool Evaluate(double const* const* parameters, |
| 78 | double* residuals, |
| 79 | double** jacobians) const { |
| 80 | if (!return_value_) { |
| 81 | return false; |
| 82 | } |
| 83 | // Compute a . x. |
| 84 | double ax = 0; |
| 85 | for (int j = 0; j < arity_; ++j) { |
| 86 | for (int u = 0; u < parameter_block_sizes()[j]; ++u) { |
| 87 | ax += a_[j][u] * parameters[j][u]; |
| 88 | } |
| 89 | } |
| 90 | |
| 91 | // This is the cost, but also appears as a factor |
| 92 | // in the derivatives. |
| 93 | double f = *residuals = exp(-ax); |
| 94 | |
| 95 | // Accumulate 1st order derivatives. |
| 96 | if (jacobians) { |
| 97 | for (int j = 0; j < arity_; ++j) { |
| 98 | if (jacobians[j]) { |
| 99 | for (int u = 0; u < parameter_block_sizes()[j]; ++u) { |
| 100 | // See comments before class. |
| 101 | jacobians[j][u] = -f * a_[j][u]; |
| 102 | } |
| 103 | } |
| 104 | } |
| 105 | } |
| 106 | |
| 107 | return true; |
| 108 | } |
| 109 | |
| 110 | void SetReturnValue(bool return_value) { return_value_ = return_value; } |
| 111 | |
| 112 | private: |
| 113 | int arity_; |
| 114 | bool return_value_; |
| 115 | vector<vector<double>> a_; // our vectors. |
| 116 | }; |
| 117 | |
| 118 | class BadTestTerm : public CostFunction { |
| 119 | public: |
| 120 | BadTestTerm(int arity, int const* dim) : arity_(arity) { |
| 121 | // Make 'arity' random vectors. |
| 122 | a_.resize(arity_); |
| 123 | for (int j = 0; j < arity_; ++j) { |
| 124 | a_[j].resize(dim[j]); |
| 125 | for (int u = 0; u < dim[j]; ++u) { |
| 126 | a_[j][u] = 2.0 * RandDouble() - 1.0; |
| 127 | } |
| 128 | } |
| 129 | |
| 130 | for (int i = 0; i < arity_; i++) { |
| 131 | mutable_parameter_block_sizes()->push_back(dim[i]); |
| 132 | } |
| 133 | set_num_residuals(1); |
| 134 | } |
| 135 | |
| 136 | bool Evaluate(double const* const* parameters, |
| 137 | double* residuals, |
| 138 | double** jacobians) const { |
| 139 | // Compute a . x. |
| 140 | double ax = 0; |
| 141 | for (int j = 0; j < arity_; ++j) { |
| 142 | for (int u = 0; u < parameter_block_sizes()[j]; ++u) { |
| 143 | ax += a_[j][u] * parameters[j][u]; |
| 144 | } |
| 145 | } |
| 146 | |
| 147 | // This is the cost, but also appears as a factor |
| 148 | // in the derivatives. |
| 149 | double f = *residuals = exp(-ax); |
| 150 | |
| 151 | // Accumulate 1st order derivatives. |
| 152 | if (jacobians) { |
| 153 | for (int j = 0; j < arity_; ++j) { |
| 154 | if (jacobians[j]) { |
| 155 | for (int u = 0; u < parameter_block_sizes()[j]; ++u) { |
| 156 | // See comments before class. |
| 157 | jacobians[j][u] = -f * a_[j][u] + 0.001; |
| 158 | } |
| 159 | } |
| 160 | } |
| 161 | } |
| 162 | |
| 163 | return true; |
| 164 | } |
| 165 | |
| 166 | private: |
| 167 | int arity_; |
| 168 | vector<vector<double>> a_; // our vectors. |
| 169 | }; |
| 170 | |
| 171 | const double kTolerance = 1e-6; |
| 172 | |
| 173 | void CheckDimensions(const GradientChecker::ProbeResults& results, |
| 174 | const std::vector<int>& parameter_sizes, |
| 175 | const std::vector<int>& local_parameter_sizes, |
| 176 | int residual_size) { |
| 177 | CHECK_EQ(parameter_sizes.size(), local_parameter_sizes.size()); |
| 178 | int num_parameters = parameter_sizes.size(); |
| 179 | ASSERT_EQ(residual_size, results.residuals.size()); |
| 180 | ASSERT_EQ(num_parameters, results.local_jacobians.size()); |
| 181 | ASSERT_EQ(num_parameters, results.local_numeric_jacobians.size()); |
| 182 | ASSERT_EQ(num_parameters, results.jacobians.size()); |
| 183 | ASSERT_EQ(num_parameters, results.numeric_jacobians.size()); |
| 184 | for (int i = 0; i < num_parameters; ++i) { |
| 185 | EXPECT_EQ(residual_size, results.local_jacobians.at(i).rows()); |
| 186 | EXPECT_EQ(local_parameter_sizes[i], results.local_jacobians.at(i).cols()); |
| 187 | EXPECT_EQ(residual_size, results.local_numeric_jacobians.at(i).rows()); |
| 188 | EXPECT_EQ(local_parameter_sizes[i], |
| 189 | results.local_numeric_jacobians.at(i).cols()); |
| 190 | EXPECT_EQ(residual_size, results.jacobians.at(i).rows()); |
| 191 | EXPECT_EQ(parameter_sizes[i], results.jacobians.at(i).cols()); |
| 192 | EXPECT_EQ(residual_size, results.numeric_jacobians.at(i).rows()); |
| 193 | EXPECT_EQ(parameter_sizes[i], results.numeric_jacobians.at(i).cols()); |
| 194 | } |
| 195 | } |
| 196 | |
| 197 | TEST(GradientChecker, SmokeTest) { |
| 198 | srand(5); |
| 199 | |
| 200 | // Test with 3 blocks of size 2, 3 and 4. |
| 201 | int const num_parameters = 3; |
| 202 | std::vector<int> parameter_sizes(3); |
| 203 | parameter_sizes[0] = 2; |
| 204 | parameter_sizes[1] = 3; |
| 205 | parameter_sizes[2] = 4; |
| 206 | |
| 207 | // Make a random set of blocks. |
| 208 | FixedArray<double*> parameters(num_parameters); |
| 209 | for (int j = 0; j < num_parameters; ++j) { |
| 210 | parameters[j] = new double[parameter_sizes[j]]; |
| 211 | for (int u = 0; u < parameter_sizes[j]; ++u) { |
| 212 | parameters[j][u] = 2.0 * RandDouble() - 1.0; |
| 213 | } |
| 214 | } |
| 215 | |
| 216 | NumericDiffOptions numeric_diff_options; |
| 217 | GradientChecker::ProbeResults results; |
| 218 | |
| 219 | // Test that Probe returns true for correct Jacobians. |
| 220 | GoodTestTerm good_term(num_parameters, parameter_sizes.data()); |
| 221 | GradientChecker good_gradient_checker(&good_term, NULL, numeric_diff_options); |
| 222 | EXPECT_TRUE(good_gradient_checker.Probe(parameters.get(), kTolerance, NULL)); |
| 223 | EXPECT_TRUE( |
| 224 | good_gradient_checker.Probe(parameters.get(), kTolerance, &results)) |
| 225 | << results.error_log; |
| 226 | |
| 227 | // Check that results contain sensible data. |
| 228 | ASSERT_EQ(results.return_value, true); |
| 229 | ASSERT_EQ(results.residuals.size(), 1); |
| 230 | CheckDimensions(results, parameter_sizes, parameter_sizes, 1); |
| 231 | EXPECT_GE(results.maximum_relative_error, 0.0); |
| 232 | EXPECT_TRUE(results.error_log.empty()); |
| 233 | |
| 234 | // Test that if the cost function return false, Probe should return false. |
| 235 | good_term.SetReturnValue(false); |
| 236 | EXPECT_FALSE(good_gradient_checker.Probe(parameters.get(), kTolerance, NULL)); |
| 237 | EXPECT_FALSE( |
| 238 | good_gradient_checker.Probe(parameters.get(), kTolerance, &results)) |
| 239 | << results.error_log; |
| 240 | |
| 241 | // Check that results contain sensible data. |
| 242 | ASSERT_EQ(results.return_value, false); |
| 243 | ASSERT_EQ(results.residuals.size(), 1); |
| 244 | CheckDimensions(results, parameter_sizes, parameter_sizes, 1); |
| 245 | for (int i = 0; i < num_parameters; ++i) { |
| 246 | EXPECT_EQ(results.local_jacobians.at(i).norm(), 0); |
| 247 | EXPECT_EQ(results.local_numeric_jacobians.at(i).norm(), 0); |
| 248 | } |
| 249 | EXPECT_EQ(results.maximum_relative_error, 0.0); |
| 250 | EXPECT_FALSE(results.error_log.empty()); |
| 251 | |
| 252 | // Test that Probe returns false for incorrect Jacobians. |
| 253 | BadTestTerm bad_term(num_parameters, parameter_sizes.data()); |
| 254 | GradientChecker bad_gradient_checker(&bad_term, NULL, numeric_diff_options); |
| 255 | EXPECT_FALSE(bad_gradient_checker.Probe(parameters.get(), kTolerance, NULL)); |
| 256 | EXPECT_FALSE( |
| 257 | bad_gradient_checker.Probe(parameters.get(), kTolerance, &results)); |
| 258 | |
| 259 | // Check that results contain sensible data. |
| 260 | ASSERT_EQ(results.return_value, true); |
| 261 | ASSERT_EQ(results.residuals.size(), 1); |
| 262 | CheckDimensions(results, parameter_sizes, parameter_sizes, 1); |
| 263 | EXPECT_GT(results.maximum_relative_error, kTolerance); |
| 264 | EXPECT_FALSE(results.error_log.empty()); |
| 265 | |
| 266 | // Setting a high threshold should make the test pass. |
| 267 | EXPECT_TRUE(bad_gradient_checker.Probe(parameters.get(), 1.0, &results)); |
| 268 | |
| 269 | // Check that results contain sensible data. |
| 270 | ASSERT_EQ(results.return_value, true); |
| 271 | ASSERT_EQ(results.residuals.size(), 1); |
| 272 | CheckDimensions(results, parameter_sizes, parameter_sizes, 1); |
| 273 | EXPECT_GT(results.maximum_relative_error, 0.0); |
| 274 | EXPECT_TRUE(results.error_log.empty()); |
| 275 | |
| 276 | for (int j = 0; j < num_parameters; j++) { |
| 277 | delete[] parameters[j]; |
| 278 | } |
| 279 | } |
| 280 | |
| 281 | /** |
| 282 | * Helper cost function that multiplies the parameters by the given jacobians |
| 283 | * and adds a constant offset. |
| 284 | */ |
| 285 | class LinearCostFunction : public CostFunction { |
| 286 | public: |
| 287 | explicit LinearCostFunction(const Vector& residuals_offset) |
| 288 | : residuals_offset_(residuals_offset) { |
| 289 | set_num_residuals(residuals_offset_.size()); |
| 290 | } |
| 291 | |
| 292 | virtual bool Evaluate(double const* const* parameter_ptrs, |
| 293 | double* residuals_ptr, |
| 294 | double** residual_J_params) const { |
| 295 | CHECK_GE(residual_J_params_.size(), 0.0); |
| 296 | VectorRef residuals(residuals_ptr, residual_J_params_[0].rows()); |
| 297 | residuals = residuals_offset_; |
| 298 | |
| 299 | for (size_t i = 0; i < residual_J_params_.size(); ++i) { |
| 300 | const Matrix& residual_J_param = residual_J_params_[i]; |
| 301 | int parameter_size = residual_J_param.cols(); |
| 302 | ConstVectorRef param(parameter_ptrs[i], parameter_size); |
| 303 | |
| 304 | // Compute residual. |
| 305 | residuals += residual_J_param * param; |
| 306 | |
| 307 | // Return Jacobian. |
| 308 | if (residual_J_params != NULL && residual_J_params[i] != NULL) { |
| 309 | Eigen::Map<Matrix> residual_J_param_out(residual_J_params[i], |
| 310 | residual_J_param.rows(), |
| 311 | residual_J_param.cols()); |
| 312 | if (jacobian_offsets_.count(i) != 0) { |
| 313 | residual_J_param_out = residual_J_param + jacobian_offsets_.at(i); |
| 314 | } else { |
| 315 | residual_J_param_out = residual_J_param; |
| 316 | } |
| 317 | } |
| 318 | } |
| 319 | return true; |
| 320 | } |
| 321 | |
| 322 | void AddParameter(const Matrix& residual_J_param) { |
| 323 | CHECK_EQ(num_residuals(), residual_J_param.rows()); |
| 324 | residual_J_params_.push_back(residual_J_param); |
| 325 | mutable_parameter_block_sizes()->push_back(residual_J_param.cols()); |
| 326 | } |
| 327 | |
| 328 | /// Add offset to the given Jacobian before returning it from Evaluate(), |
| 329 | /// thus introducing an error in the comutation. |
| 330 | void SetJacobianOffset(size_t index, Matrix offset) { |
| 331 | CHECK_LT(index, residual_J_params_.size()); |
| 332 | CHECK_EQ(residual_J_params_[index].rows(), offset.rows()); |
| 333 | CHECK_EQ(residual_J_params_[index].cols(), offset.cols()); |
| 334 | jacobian_offsets_[index] = offset; |
| 335 | } |
| 336 | |
| 337 | private: |
| 338 | std::vector<Matrix> residual_J_params_; |
| 339 | std::map<int, Matrix> jacobian_offsets_; |
| 340 | Vector residuals_offset_; |
| 341 | }; |
| 342 | |
| 343 | /** |
| 344 | * Helper local parameterization that multiplies the delta vector by the given |
| 345 | * jacobian and adds it to the parameter. |
| 346 | */ |
| 347 | class MatrixParameterization : public LocalParameterization { |
| 348 | public: |
| 349 | virtual bool Plus(const double* x, |
| 350 | const double* delta, |
| 351 | double* x_plus_delta) const { |
| 352 | VectorRef(x_plus_delta, GlobalSize()) = |
| 353 | ConstVectorRef(x, GlobalSize()) + |
| 354 | (global_J_local * ConstVectorRef(delta, LocalSize())); |
| 355 | return true; |
| 356 | } |
| 357 | |
| 358 | virtual bool ComputeJacobian(const double* /*x*/, double* jacobian) const { |
| 359 | MatrixRef(jacobian, GlobalSize(), LocalSize()) = global_J_local; |
| 360 | return true; |
| 361 | } |
| 362 | |
| 363 | virtual int GlobalSize() const { return global_J_local.rows(); } |
| 364 | virtual int LocalSize() const { return global_J_local.cols(); } |
| 365 | |
| 366 | Matrix global_J_local; |
| 367 | }; |
| 368 | |
| 369 | // Helper function to compare two Eigen matrices (used in the test below). |
| 370 | void ExpectMatricesClose(Matrix p, Matrix q, double tolerance) { |
| 371 | ASSERT_EQ(p.rows(), q.rows()); |
| 372 | ASSERT_EQ(p.cols(), q.cols()); |
| 373 | ExpectArraysClose(p.size(), p.data(), q.data(), tolerance); |
| 374 | } |
| 375 | |
| 376 | TEST(GradientChecker, TestCorrectnessWithLocalParameterizations) { |
| 377 | // Create cost function. |
| 378 | Eigen::Vector3d residual_offset(100.0, 200.0, 300.0); |
| 379 | LinearCostFunction cost_function(residual_offset); |
| 380 | Eigen::Matrix<double, 3, 3, Eigen::RowMajor> j0; |
| 381 | j0.row(0) << 1.0, 2.0, 3.0; |
| 382 | j0.row(1) << 4.0, 5.0, 6.0; |
| 383 | j0.row(2) << 7.0, 8.0, 9.0; |
| 384 | Eigen::Matrix<double, 3, 2, Eigen::RowMajor> j1; |
| 385 | j1.row(0) << 10.0, 11.0; |
| 386 | j1.row(1) << 12.0, 13.0; |
| 387 | j1.row(2) << 14.0, 15.0; |
| 388 | |
| 389 | Eigen::Vector3d param0(1.0, 2.0, 3.0); |
| 390 | Eigen::Vector2d param1(4.0, 5.0); |
| 391 | |
| 392 | cost_function.AddParameter(j0); |
| 393 | cost_function.AddParameter(j1); |
| 394 | |
| 395 | std::vector<int> parameter_sizes(2); |
| 396 | parameter_sizes[0] = 3; |
| 397 | parameter_sizes[1] = 2; |
| 398 | std::vector<int> local_parameter_sizes(2); |
| 399 | local_parameter_sizes[0] = 2; |
| 400 | local_parameter_sizes[1] = 2; |
| 401 | |
| 402 | // Test cost function for correctness. |
| 403 | Eigen::Matrix<double, 3, 3, Eigen::RowMajor> j1_out; |
| 404 | Eigen::Matrix<double, 3, 2, Eigen::RowMajor> j2_out; |
| 405 | Eigen::Vector3d residual; |
| 406 | std::vector<const double*> parameters(2); |
| 407 | parameters[0] = param0.data(); |
| 408 | parameters[1] = param1.data(); |
| 409 | std::vector<double*> jacobians(2); |
| 410 | jacobians[0] = j1_out.data(); |
| 411 | jacobians[1] = j2_out.data(); |
| 412 | cost_function.Evaluate(parameters.data(), residual.data(), jacobians.data()); |
| 413 | |
| 414 | Matrix residual_expected = residual_offset + j0 * param0 + j1 * param1; |
| 415 | |
| 416 | ExpectMatricesClose(j1_out, j0, std::numeric_limits<double>::epsilon()); |
| 417 | ExpectMatricesClose(j2_out, j1, std::numeric_limits<double>::epsilon()); |
| 418 | ExpectMatricesClose(residual, residual_expected, kTolerance); |
| 419 | |
| 420 | // Create local parameterization. |
| 421 | Eigen::Matrix<double, 3, 2, Eigen::RowMajor> global_J_local; |
| 422 | global_J_local.row(0) << 1.5, 2.5; |
| 423 | global_J_local.row(1) << 3.5, 4.5; |
| 424 | global_J_local.row(2) << 5.5, 6.5; |
| 425 | |
| 426 | MatrixParameterization parameterization; |
| 427 | parameterization.global_J_local = global_J_local; |
| 428 | |
| 429 | // Test local parameterization for correctness. |
| 430 | Eigen::Vector3d x(7.0, 8.0, 9.0); |
| 431 | Eigen::Vector2d delta(10.0, 11.0); |
| 432 | |
| 433 | Eigen::Matrix<double, 3, 2, Eigen::RowMajor> global_J_local_out; |
| 434 | parameterization.ComputeJacobian(x.data(), global_J_local_out.data()); |
| 435 | ExpectMatricesClose(global_J_local_out, |
| 436 | global_J_local, |
| 437 | std::numeric_limits<double>::epsilon()); |
| 438 | |
| 439 | Eigen::Vector3d x_plus_delta; |
| 440 | parameterization.Plus(x.data(), delta.data(), x_plus_delta.data()); |
| 441 | Eigen::Vector3d x_plus_delta_expected = x + (global_J_local * delta); |
| 442 | ExpectMatricesClose(x_plus_delta, x_plus_delta_expected, kTolerance); |
| 443 | |
| 444 | // Now test GradientChecker. |
| 445 | std::vector<const LocalParameterization*> parameterizations(2); |
| 446 | parameterizations[0] = ¶meterization; |
| 447 | parameterizations[1] = NULL; |
| 448 | NumericDiffOptions numeric_diff_options; |
| 449 | GradientChecker::ProbeResults results; |
| 450 | GradientChecker gradient_checker( |
| 451 | &cost_function, ¶meterizations, numeric_diff_options); |
| 452 | |
| 453 | Problem::Options problem_options; |
| 454 | problem_options.cost_function_ownership = DO_NOT_TAKE_OWNERSHIP; |
| 455 | problem_options.local_parameterization_ownership = DO_NOT_TAKE_OWNERSHIP; |
| 456 | Problem problem(problem_options); |
| 457 | Eigen::Vector3d param0_solver; |
| 458 | Eigen::Vector2d param1_solver; |
| 459 | problem.AddParameterBlock(param0_solver.data(), 3, ¶meterization); |
| 460 | problem.AddParameterBlock(param1_solver.data(), 2); |
| 461 | problem.AddResidualBlock( |
| 462 | &cost_function, NULL, param0_solver.data(), param1_solver.data()); |
| 463 | Solver::Options solver_options; |
| 464 | solver_options.check_gradients = true; |
| 465 | solver_options.initial_trust_region_radius = 1e10; |
| 466 | Solver solver; |
| 467 | Solver::Summary summary; |
| 468 | |
| 469 | // First test case: everything is correct. |
| 470 | EXPECT_TRUE(gradient_checker.Probe(parameters.data(), kTolerance, NULL)); |
| 471 | EXPECT_TRUE(gradient_checker.Probe(parameters.data(), kTolerance, &results)) |
| 472 | << results.error_log; |
| 473 | |
| 474 | // Check that results contain correct data. |
| 475 | ASSERT_EQ(results.return_value, true); |
| 476 | ExpectMatricesClose( |
| 477 | results.residuals, residual, std::numeric_limits<double>::epsilon()); |
| 478 | CheckDimensions(results, parameter_sizes, local_parameter_sizes, 3); |
| 479 | ExpectMatricesClose( |
| 480 | results.local_jacobians.at(0), j0 * global_J_local, kTolerance); |
| 481 | ExpectMatricesClose(results.local_jacobians.at(1), |
| 482 | j1, |
| 483 | std::numeric_limits<double>::epsilon()); |
| 484 | ExpectMatricesClose( |
| 485 | results.local_numeric_jacobians.at(0), j0 * global_J_local, kTolerance); |
| 486 | ExpectMatricesClose(results.local_numeric_jacobians.at(1), j1, kTolerance); |
| 487 | ExpectMatricesClose( |
| 488 | results.jacobians.at(0), j0, std::numeric_limits<double>::epsilon()); |
| 489 | ExpectMatricesClose( |
| 490 | results.jacobians.at(1), j1, std::numeric_limits<double>::epsilon()); |
| 491 | ExpectMatricesClose(results.numeric_jacobians.at(0), j0, kTolerance); |
| 492 | ExpectMatricesClose(results.numeric_jacobians.at(1), j1, kTolerance); |
| 493 | EXPECT_GE(results.maximum_relative_error, 0.0); |
| 494 | EXPECT_TRUE(results.error_log.empty()); |
| 495 | |
| 496 | // Test interaction with the 'check_gradients' option in Solver. |
| 497 | param0_solver = param0; |
| 498 | param1_solver = param1; |
| 499 | solver.Solve(solver_options, &problem, &summary); |
| 500 | EXPECT_EQ(CONVERGENCE, summary.termination_type); |
| 501 | EXPECT_LE(summary.final_cost, 1e-12); |
| 502 | |
| 503 | // Second test case: Mess up reported derivatives with respect to 3rd |
| 504 | // component of 1st parameter. Check should fail. |
| 505 | Eigen::Matrix<double, 3, 3, Eigen::RowMajor> j0_offset; |
| 506 | j0_offset.setZero(); |
| 507 | j0_offset.col(2).setConstant(0.001); |
| 508 | cost_function.SetJacobianOffset(0, j0_offset); |
| 509 | EXPECT_FALSE(gradient_checker.Probe(parameters.data(), kTolerance, NULL)); |
| 510 | EXPECT_FALSE(gradient_checker.Probe(parameters.data(), kTolerance, &results)) |
| 511 | << results.error_log; |
| 512 | |
| 513 | // Check that results contain correct data. |
| 514 | ASSERT_EQ(results.return_value, true); |
| 515 | ExpectMatricesClose( |
| 516 | results.residuals, residual, std::numeric_limits<double>::epsilon()); |
| 517 | CheckDimensions(results, parameter_sizes, local_parameter_sizes, 3); |
| 518 | ASSERT_EQ(results.local_jacobians.size(), 2); |
| 519 | ASSERT_EQ(results.local_numeric_jacobians.size(), 2); |
| 520 | ExpectMatricesClose(results.local_jacobians.at(0), |
| 521 | (j0 + j0_offset) * global_J_local, |
| 522 | kTolerance); |
| 523 | ExpectMatricesClose(results.local_jacobians.at(1), |
| 524 | j1, |
| 525 | std::numeric_limits<double>::epsilon()); |
| 526 | ExpectMatricesClose( |
| 527 | results.local_numeric_jacobians.at(0), j0 * global_J_local, kTolerance); |
| 528 | ExpectMatricesClose(results.local_numeric_jacobians.at(1), j1, kTolerance); |
| 529 | ExpectMatricesClose(results.jacobians.at(0), j0 + j0_offset, kTolerance); |
| 530 | ExpectMatricesClose( |
| 531 | results.jacobians.at(1), j1, std::numeric_limits<double>::epsilon()); |
| 532 | ExpectMatricesClose(results.numeric_jacobians.at(0), j0, kTolerance); |
| 533 | ExpectMatricesClose(results.numeric_jacobians.at(1), j1, kTolerance); |
| 534 | EXPECT_GT(results.maximum_relative_error, 0.0); |
| 535 | EXPECT_FALSE(results.error_log.empty()); |
| 536 | |
| 537 | // Test interaction with the 'check_gradients' option in Solver. |
| 538 | param0_solver = param0; |
| 539 | param1_solver = param1; |
| 540 | solver.Solve(solver_options, &problem, &summary); |
| 541 | EXPECT_EQ(FAILURE, summary.termination_type); |
| 542 | |
| 543 | // Now, zero out the local parameterization Jacobian of the 1st parameter |
| 544 | // with respect to the 3rd component. This makes the combination of |
| 545 | // cost function and local parameterization return correct values again. |
| 546 | parameterization.global_J_local.row(2).setZero(); |
| 547 | |
| 548 | // Verify that the gradient checker does not treat this as an error. |
| 549 | EXPECT_TRUE(gradient_checker.Probe(parameters.data(), kTolerance, &results)) |
| 550 | << results.error_log; |
| 551 | |
| 552 | // Check that results contain correct data. |
| 553 | ASSERT_EQ(results.return_value, true); |
| 554 | ExpectMatricesClose( |
| 555 | results.residuals, residual, std::numeric_limits<double>::epsilon()); |
| 556 | CheckDimensions(results, parameter_sizes, local_parameter_sizes, 3); |
| 557 | ASSERT_EQ(results.local_jacobians.size(), 2); |
| 558 | ASSERT_EQ(results.local_numeric_jacobians.size(), 2); |
| 559 | ExpectMatricesClose(results.local_jacobians.at(0), |
| 560 | (j0 + j0_offset) * parameterization.global_J_local, |
| 561 | kTolerance); |
| 562 | ExpectMatricesClose(results.local_jacobians.at(1), |
| 563 | j1, |
| 564 | std::numeric_limits<double>::epsilon()); |
| 565 | ExpectMatricesClose(results.local_numeric_jacobians.at(0), |
| 566 | j0 * parameterization.global_J_local, |
| 567 | kTolerance); |
| 568 | ExpectMatricesClose(results.local_numeric_jacobians.at(1), j1, kTolerance); |
| 569 | ExpectMatricesClose(results.jacobians.at(0), j0 + j0_offset, kTolerance); |
| 570 | ExpectMatricesClose( |
| 571 | results.jacobians.at(1), j1, std::numeric_limits<double>::epsilon()); |
| 572 | ExpectMatricesClose(results.numeric_jacobians.at(0), j0, kTolerance); |
| 573 | ExpectMatricesClose(results.numeric_jacobians.at(1), j1, kTolerance); |
| 574 | EXPECT_GE(results.maximum_relative_error, 0.0); |
| 575 | EXPECT_TRUE(results.error_log.empty()); |
| 576 | |
| 577 | // Test interaction with the 'check_gradients' option in Solver. |
| 578 | param0_solver = param0; |
| 579 | param1_solver = param1; |
| 580 | solver.Solve(solver_options, &problem, &summary); |
| 581 | EXPECT_EQ(CONVERGENCE, summary.termination_type); |
| 582 | EXPECT_LE(summary.final_cost, 1e-12); |
| 583 | } |
| 584 | |
| 585 | } // namespace internal |
| 586 | } // namespace ceres |