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 | // Authors: wjr@google.com (William Rucklidge), |
| 30 | // keir@google.com (Keir Mierle), |
| 31 | // dgossow@google.com (David Gossow) |
| 32 | |
| 33 | #include "ceres/gradient_checker.h" |
| 34 | |
| 35 | #include <algorithm> |
| 36 | #include <cmath> |
| 37 | #include <cstdint> |
| 38 | #include <numeric> |
| 39 | #include <string> |
| 40 | #include <vector> |
| 41 | |
| 42 | #include "ceres/is_close.h" |
| 43 | #include "ceres/stringprintf.h" |
| 44 | #include "ceres/types.h" |
| 45 | |
| 46 | namespace ceres { |
| 47 | |
| 48 | using internal::IsClose; |
| 49 | using internal::StringAppendF; |
| 50 | using internal::StringPrintf; |
| 51 | using std::string; |
| 52 | using std::vector; |
| 53 | |
| 54 | namespace { |
| 55 | // Evaluate the cost function and transform the returned Jacobians to |
| 56 | // the local space of the respective local parameterizations. |
| 57 | bool EvaluateCostFunction( |
| 58 | const ceres::CostFunction* function, |
| 59 | double const* const * parameters, |
| 60 | const std::vector<const ceres::LocalParameterization*>& |
| 61 | local_parameterizations, |
| 62 | Vector* residuals, |
| 63 | std::vector<Matrix>* jacobians, |
| 64 | std::vector<Matrix>* local_jacobians) { |
| 65 | CHECK(residuals != nullptr); |
| 66 | CHECK(jacobians != nullptr); |
| 67 | CHECK(local_jacobians != nullptr); |
| 68 | |
| 69 | const vector<int32_t>& block_sizes = function->parameter_block_sizes(); |
| 70 | const int num_parameter_blocks = block_sizes.size(); |
| 71 | |
| 72 | // Allocate Jacobian matrices in local space. |
| 73 | local_jacobians->resize(num_parameter_blocks); |
| 74 | vector<double*> local_jacobian_data(num_parameter_blocks); |
| 75 | for (int i = 0; i < num_parameter_blocks; ++i) { |
| 76 | int block_size = block_sizes.at(i); |
| 77 | if (local_parameterizations.at(i) != NULL) { |
| 78 | block_size = local_parameterizations.at(i)->LocalSize(); |
| 79 | } |
| 80 | local_jacobians->at(i).resize(function->num_residuals(), block_size); |
| 81 | local_jacobians->at(i).setZero(); |
| 82 | local_jacobian_data.at(i) = local_jacobians->at(i).data(); |
| 83 | } |
| 84 | |
| 85 | // Allocate Jacobian matrices in global space. |
| 86 | jacobians->resize(num_parameter_blocks); |
| 87 | vector<double*> jacobian_data(num_parameter_blocks); |
| 88 | for (int i = 0; i < num_parameter_blocks; ++i) { |
| 89 | jacobians->at(i).resize(function->num_residuals(), block_sizes.at(i)); |
| 90 | jacobians->at(i).setZero(); |
| 91 | jacobian_data.at(i) = jacobians->at(i).data(); |
| 92 | } |
| 93 | |
| 94 | // Compute residuals & jacobians. |
| 95 | CHECK_NE(0, function->num_residuals()); |
| 96 | residuals->resize(function->num_residuals()); |
| 97 | residuals->setZero(); |
| 98 | if (!function->Evaluate(parameters, residuals->data(), |
| 99 | jacobian_data.data())) { |
| 100 | return false; |
| 101 | } |
| 102 | |
| 103 | // Convert Jacobians from global to local space. |
| 104 | for (size_t i = 0; i < local_jacobians->size(); ++i) { |
| 105 | if (local_parameterizations.at(i) == NULL) { |
| 106 | local_jacobians->at(i) = jacobians->at(i); |
| 107 | } else { |
| 108 | int global_size = local_parameterizations.at(i)->GlobalSize(); |
| 109 | int local_size = local_parameterizations.at(i)->LocalSize(); |
| 110 | CHECK_EQ(jacobians->at(i).cols(), global_size); |
| 111 | Matrix global_J_local(global_size, local_size); |
| 112 | local_parameterizations.at(i)->ComputeJacobian( |
| 113 | parameters[i], global_J_local.data()); |
| 114 | local_jacobians->at(i).noalias() = jacobians->at(i) * global_J_local; |
| 115 | } |
| 116 | } |
| 117 | return true; |
| 118 | } |
| 119 | } // namespace |
| 120 | |
| 121 | GradientChecker::GradientChecker( |
| 122 | const CostFunction* function, |
| 123 | const vector<const LocalParameterization*>* local_parameterizations, |
| 124 | const NumericDiffOptions& options) : |
| 125 | function_(function) { |
| 126 | CHECK(function != nullptr); |
| 127 | if (local_parameterizations != NULL) { |
| 128 | local_parameterizations_ = *local_parameterizations; |
| 129 | } else { |
| 130 | local_parameterizations_.resize(function->parameter_block_sizes().size(), |
| 131 | NULL); |
| 132 | } |
| 133 | DynamicNumericDiffCostFunction<CostFunction, CENTRAL>* |
| 134 | finite_diff_cost_function = |
| 135 | new DynamicNumericDiffCostFunction<CostFunction, CENTRAL>( |
| 136 | function, DO_NOT_TAKE_OWNERSHIP, options); |
| 137 | finite_diff_cost_function_.reset(finite_diff_cost_function); |
| 138 | |
| 139 | const vector<int32_t>& parameter_block_sizes = |
| 140 | function->parameter_block_sizes(); |
| 141 | const int num_parameter_blocks = parameter_block_sizes.size(); |
| 142 | for (int i = 0; i < num_parameter_blocks; ++i) { |
| 143 | finite_diff_cost_function->AddParameterBlock(parameter_block_sizes[i]); |
| 144 | } |
| 145 | finite_diff_cost_function->SetNumResiduals(function->num_residuals()); |
| 146 | } |
| 147 | |
| 148 | bool GradientChecker::Probe(double const* const * parameters, |
| 149 | double relative_precision, |
| 150 | ProbeResults* results_param) const { |
| 151 | int num_residuals = function_->num_residuals(); |
| 152 | |
| 153 | // Make sure that we have a place to store results, no matter if the user has |
| 154 | // provided an output argument. |
| 155 | ProbeResults* results; |
| 156 | ProbeResults results_local; |
| 157 | if (results_param != NULL) { |
| 158 | results = results_param; |
| 159 | results->residuals.resize(0); |
| 160 | results->jacobians.clear(); |
| 161 | results->numeric_jacobians.clear(); |
| 162 | results->local_jacobians.clear(); |
| 163 | results->local_numeric_jacobians.clear(); |
| 164 | results->error_log.clear(); |
| 165 | } else { |
| 166 | results = &results_local; |
| 167 | } |
| 168 | results->maximum_relative_error = 0.0; |
| 169 | results->return_value = true; |
| 170 | |
| 171 | // Evaluate the derivative using the user supplied code. |
| 172 | vector<Matrix>& jacobians = results->jacobians; |
| 173 | vector<Matrix>& local_jacobians = results->local_jacobians; |
| 174 | if (!EvaluateCostFunction(function_, parameters, local_parameterizations_, |
| 175 | &results->residuals, &jacobians, &local_jacobians)) { |
| 176 | results->error_log = "Function evaluation with Jacobians failed."; |
| 177 | results->return_value = false; |
| 178 | } |
| 179 | |
| 180 | // Evaluate the derivative using numeric derivatives. |
| 181 | vector<Matrix>& numeric_jacobians = results->numeric_jacobians; |
| 182 | vector<Matrix>& local_numeric_jacobians = results->local_numeric_jacobians; |
| 183 | Vector finite_diff_residuals; |
| 184 | if (!EvaluateCostFunction(finite_diff_cost_function_.get(), parameters, |
| 185 | local_parameterizations_, &finite_diff_residuals, |
| 186 | &numeric_jacobians, &local_numeric_jacobians)) { |
| 187 | results->error_log += "\nFunction evaluation with numerical " |
| 188 | "differentiation failed."; |
| 189 | results->return_value = false; |
| 190 | } |
| 191 | |
| 192 | if (!results->return_value) { |
| 193 | return false; |
| 194 | } |
| 195 | |
| 196 | for (int i = 0; i < num_residuals; ++i) { |
| 197 | if (!IsClose( |
| 198 | results->residuals[i], |
| 199 | finite_diff_residuals[i], |
| 200 | relative_precision, |
| 201 | NULL, |
| 202 | NULL)) { |
| 203 | results->error_log = "Function evaluation with and without Jacobians " |
| 204 | "resulted in different residuals."; |
| 205 | LOG(INFO) << results->residuals.transpose(); |
| 206 | LOG(INFO) << finite_diff_residuals.transpose(); |
| 207 | return false; |
| 208 | } |
| 209 | } |
| 210 | |
| 211 | // See if any elements have relative error larger than the threshold. |
| 212 | int num_bad_jacobian_components = 0; |
| 213 | double& worst_relative_error = results->maximum_relative_error; |
| 214 | worst_relative_error = 0; |
| 215 | |
| 216 | // Accumulate the error message for all the jacobians, since it won't get |
| 217 | // output if there are no bad jacobian components. |
| 218 | string error_log; |
| 219 | for (int k = 0; k < function_->parameter_block_sizes().size(); k++) { |
| 220 | StringAppendF(&error_log, |
| 221 | "========== " |
| 222 | "Jacobian for " "block %d: (%ld by %ld)) " |
| 223 | "==========\n", |
| 224 | k, |
| 225 | static_cast<long>(local_jacobians[k].rows()), |
| 226 | static_cast<long>(local_jacobians[k].cols())); |
| 227 | // The funny spacing creates appropriately aligned column headers. |
| 228 | error_log += |
| 229 | " block row col user dx/dy num diff dx/dy " |
| 230 | "abs error relative error parameter residual\n"; |
| 231 | |
| 232 | for (int i = 0; i < local_jacobians[k].rows(); i++) { |
| 233 | for (int j = 0; j < local_jacobians[k].cols(); j++) { |
| 234 | double term_jacobian = local_jacobians[k](i, j); |
| 235 | double finite_jacobian = local_numeric_jacobians[k](i, j); |
| 236 | double relative_error, absolute_error; |
| 237 | bool bad_jacobian_entry = |
| 238 | !IsClose(term_jacobian, |
| 239 | finite_jacobian, |
| 240 | relative_precision, |
| 241 | &relative_error, |
| 242 | &absolute_error); |
| 243 | worst_relative_error = std::max(worst_relative_error, relative_error); |
| 244 | |
| 245 | StringAppendF(&error_log, |
| 246 | "%6d %4d %4d %17g %17g %17g %17g %17g %17g", |
| 247 | k, i, j, |
| 248 | term_jacobian, finite_jacobian, |
| 249 | absolute_error, relative_error, |
| 250 | parameters[k][j], |
| 251 | results->residuals[i]); |
| 252 | |
| 253 | if (bad_jacobian_entry) { |
| 254 | num_bad_jacobian_components++; |
| 255 | StringAppendF( |
| 256 | &error_log, |
| 257 | " ------ (%d,%d,%d) Relative error worse than %g", |
| 258 | k, i, j, relative_precision); |
| 259 | } |
| 260 | error_log += "\n"; |
| 261 | } |
| 262 | } |
| 263 | } |
| 264 | |
| 265 | // Since there were some bad errors, dump comprehensive debug info. |
| 266 | if (num_bad_jacobian_components) { |
| 267 | string header = StringPrintf("\nDetected %d bad Jacobian component(s). " |
| 268 | "Worst relative error was %g.\n", |
| 269 | num_bad_jacobian_components, |
| 270 | worst_relative_error); |
| 271 | results->error_log = header + "\n" + error_log; |
| 272 | return false; |
| 273 | } |
| 274 | return true; |
| 275 | } |
| 276 | |
| 277 | } // namespace ceres |