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: sameeragarwal@google.com (Sameer Agarwal) |
| 30 | |
| 31 | #include "ceres/covariance_impl.h" |
| 32 | |
| 33 | #include <algorithm> |
| 34 | #include <cstdlib> |
| 35 | #include <memory> |
| 36 | #include <numeric> |
| 37 | #include <sstream> |
| 38 | #include <unordered_set> |
| 39 | #include <utility> |
| 40 | #include <vector> |
| 41 | |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 42 | #include "Eigen/SVD" |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 43 | #include "Eigen/SparseCore" |
| 44 | #include "Eigen/SparseQR" |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 45 | #include "ceres/compressed_col_sparse_matrix_utils.h" |
| 46 | #include "ceres/compressed_row_sparse_matrix.h" |
| 47 | #include "ceres/covariance.h" |
| 48 | #include "ceres/crs_matrix.h" |
| 49 | #include "ceres/internal/eigen.h" |
| 50 | #include "ceres/map_util.h" |
| 51 | #include "ceres/parallel_for.h" |
| 52 | #include "ceres/parallel_utils.h" |
| 53 | #include "ceres/parameter_block.h" |
| 54 | #include "ceres/problem_impl.h" |
| 55 | #include "ceres/residual_block.h" |
| 56 | #include "ceres/suitesparse.h" |
| 57 | #include "ceres/wall_time.h" |
| 58 | #include "glog/logging.h" |
| 59 | |
| 60 | namespace ceres { |
| 61 | namespace internal { |
| 62 | |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 63 | using std::swap; |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 64 | |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 65 | using CovarianceBlocks = std::vector<std::pair<const double*, const double*>>; |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 66 | |
| 67 | CovarianceImpl::CovarianceImpl(const Covariance::Options& options) |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 68 | : options_(options), is_computed_(false), is_valid_(false) { |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 69 | #ifdef CERES_NO_THREADS |
| 70 | if (options_.num_threads > 1) { |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 71 | LOG(WARNING) << "No threading support is compiled into this binary; " |
| 72 | << "only options.num_threads = 1 is supported. Switching " |
| 73 | << "to single threaded mode."; |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 74 | options_.num_threads = 1; |
| 75 | } |
| 76 | #endif |
| 77 | |
| 78 | evaluate_options_.num_threads = options_.num_threads; |
| 79 | evaluate_options_.apply_loss_function = options_.apply_loss_function; |
| 80 | } |
| 81 | |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 82 | CovarianceImpl::~CovarianceImpl() {} |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 83 | |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 84 | template <typename T> |
| 85 | void CheckForDuplicates(std::vector<T> blocks) { |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 86 | sort(blocks.begin(), blocks.end()); |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 87 | typename std::vector<T>::iterator it = |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 88 | std::adjacent_find(blocks.begin(), blocks.end()); |
| 89 | if (it != blocks.end()) { |
| 90 | // In case there are duplicates, we search for their location. |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 91 | std::map<T, std::vector<int>> blocks_map; |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 92 | for (int i = 0; i < blocks.size(); ++i) { |
| 93 | blocks_map[blocks[i]].push_back(i); |
| 94 | } |
| 95 | |
| 96 | std::ostringstream duplicates; |
| 97 | while (it != blocks.end()) { |
| 98 | duplicates << "("; |
| 99 | for (int i = 0; i < blocks_map[*it].size() - 1; ++i) { |
| 100 | duplicates << blocks_map[*it][i] << ", "; |
| 101 | } |
| 102 | duplicates << blocks_map[*it].back() << ")"; |
| 103 | it = std::adjacent_find(it + 1, blocks.end()); |
| 104 | if (it < blocks.end()) { |
| 105 | duplicates << " and "; |
| 106 | } |
| 107 | } |
| 108 | |
| 109 | LOG(FATAL) << "Covariance::Compute called with duplicate blocks at " |
| 110 | << "indices " << duplicates.str(); |
| 111 | } |
| 112 | } |
| 113 | |
| 114 | bool CovarianceImpl::Compute(const CovarianceBlocks& covariance_blocks, |
| 115 | ProblemImpl* problem) { |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 116 | CheckForDuplicates<std::pair<const double*, const double*>>( |
| 117 | covariance_blocks); |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 118 | problem_ = problem; |
| 119 | parameter_block_to_row_index_.clear(); |
| 120 | covariance_matrix_.reset(NULL); |
| 121 | is_valid_ = (ComputeCovarianceSparsity(covariance_blocks, problem) && |
| 122 | ComputeCovarianceValues()); |
| 123 | is_computed_ = true; |
| 124 | return is_valid_; |
| 125 | } |
| 126 | |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 127 | bool CovarianceImpl::Compute(const std::vector<const double*>& parameter_blocks, |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 128 | ProblemImpl* problem) { |
| 129 | CheckForDuplicates<const double*>(parameter_blocks); |
| 130 | CovarianceBlocks covariance_blocks; |
| 131 | for (int i = 0; i < parameter_blocks.size(); ++i) { |
| 132 | for (int j = i; j < parameter_blocks.size(); ++j) { |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 133 | covariance_blocks.push_back( |
| 134 | std::make_pair(parameter_blocks[i], parameter_blocks[j])); |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 135 | } |
| 136 | } |
| 137 | |
| 138 | return Compute(covariance_blocks, problem); |
| 139 | } |
| 140 | |
| 141 | bool CovarianceImpl::GetCovarianceBlockInTangentOrAmbientSpace( |
| 142 | const double* original_parameter_block1, |
| 143 | const double* original_parameter_block2, |
| 144 | bool lift_covariance_to_ambient_space, |
| 145 | double* covariance_block) const { |
| 146 | CHECK(is_computed_) |
| 147 | << "Covariance::GetCovarianceBlock called before Covariance::Compute"; |
| 148 | CHECK(is_valid_) |
| 149 | << "Covariance::GetCovarianceBlock called when Covariance::Compute " |
| 150 | << "returned false."; |
| 151 | |
| 152 | // If either of the two parameter blocks is constant, then the |
| 153 | // covariance block is also zero. |
| 154 | if (constant_parameter_blocks_.count(original_parameter_block1) > 0 || |
| 155 | constant_parameter_blocks_.count(original_parameter_block2) > 0) { |
| 156 | const ProblemImpl::ParameterMap& parameter_map = problem_->parameter_map(); |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 157 | ParameterBlock* block1 = FindOrDie( |
| 158 | parameter_map, const_cast<double*>(original_parameter_block1)); |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 159 | |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 160 | ParameterBlock* block2 = FindOrDie( |
| 161 | parameter_map, const_cast<double*>(original_parameter_block2)); |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 162 | |
| 163 | const int block1_size = block1->Size(); |
| 164 | const int block2_size = block2->Size(); |
| 165 | const int block1_local_size = block1->LocalSize(); |
| 166 | const int block2_local_size = block2->LocalSize(); |
| 167 | if (!lift_covariance_to_ambient_space) { |
| 168 | MatrixRef(covariance_block, block1_local_size, block2_local_size) |
| 169 | .setZero(); |
| 170 | } else { |
| 171 | MatrixRef(covariance_block, block1_size, block2_size).setZero(); |
| 172 | } |
| 173 | return true; |
| 174 | } |
| 175 | |
| 176 | const double* parameter_block1 = original_parameter_block1; |
| 177 | const double* parameter_block2 = original_parameter_block2; |
| 178 | const bool transpose = parameter_block1 > parameter_block2; |
| 179 | if (transpose) { |
| 180 | swap(parameter_block1, parameter_block2); |
| 181 | } |
| 182 | |
| 183 | // Find where in the covariance matrix the block is located. |
| 184 | const int row_begin = |
| 185 | FindOrDie(parameter_block_to_row_index_, parameter_block1); |
| 186 | const int col_begin = |
| 187 | FindOrDie(parameter_block_to_row_index_, parameter_block2); |
| 188 | const int* rows = covariance_matrix_->rows(); |
| 189 | const int* cols = covariance_matrix_->cols(); |
| 190 | const int row_size = rows[row_begin + 1] - rows[row_begin]; |
| 191 | const int* cols_begin = cols + rows[row_begin]; |
| 192 | |
| 193 | // The only part that requires work is walking the compressed column |
| 194 | // vector to determine where the set of columns correspnding to the |
| 195 | // covariance block begin. |
| 196 | int offset = 0; |
| 197 | while (cols_begin[offset] != col_begin && offset < row_size) { |
| 198 | ++offset; |
| 199 | } |
| 200 | |
| 201 | if (offset == row_size) { |
| 202 | LOG(ERROR) << "Unable to find covariance block for " |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 203 | << original_parameter_block1 << " " << original_parameter_block2; |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 204 | return false; |
| 205 | } |
| 206 | |
| 207 | const ProblemImpl::ParameterMap& parameter_map = problem_->parameter_map(); |
| 208 | ParameterBlock* block1 = |
| 209 | FindOrDie(parameter_map, const_cast<double*>(parameter_block1)); |
| 210 | ParameterBlock* block2 = |
| 211 | FindOrDie(parameter_map, const_cast<double*>(parameter_block2)); |
| 212 | const LocalParameterization* local_param1 = block1->local_parameterization(); |
| 213 | const LocalParameterization* local_param2 = block2->local_parameterization(); |
| 214 | const int block1_size = block1->Size(); |
| 215 | const int block1_local_size = block1->LocalSize(); |
| 216 | const int block2_size = block2->Size(); |
| 217 | const int block2_local_size = block2->LocalSize(); |
| 218 | |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 219 | ConstMatrixRef cov( |
| 220 | covariance_matrix_->values() + rows[row_begin], block1_size, row_size); |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 221 | |
| 222 | // Fast path when there are no local parameterizations or if the |
| 223 | // user does not want it lifted to the ambient space. |
| 224 | if ((local_param1 == NULL && local_param2 == NULL) || |
| 225 | !lift_covariance_to_ambient_space) { |
| 226 | if (transpose) { |
| 227 | MatrixRef(covariance_block, block2_local_size, block1_local_size) = |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 228 | cov.block(0, offset, block1_local_size, block2_local_size) |
| 229 | .transpose(); |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 230 | } else { |
| 231 | MatrixRef(covariance_block, block1_local_size, block2_local_size) = |
| 232 | cov.block(0, offset, block1_local_size, block2_local_size); |
| 233 | } |
| 234 | return true; |
| 235 | } |
| 236 | |
| 237 | // If local parameterizations are used then the covariance that has |
| 238 | // been computed is in the tangent space and it needs to be lifted |
| 239 | // back to the ambient space. |
| 240 | // |
| 241 | // This is given by the formula |
| 242 | // |
| 243 | // C'_12 = J_1 C_12 J_2' |
| 244 | // |
| 245 | // Where C_12 is the local tangent space covariance for parameter |
| 246 | // blocks 1 and 2. J_1 and J_2 are respectively the local to global |
| 247 | // jacobians for parameter blocks 1 and 2. |
| 248 | // |
| 249 | // See Result 5.11 on page 142 of Hartley & Zisserman (2nd Edition) |
| 250 | // for a proof. |
| 251 | // |
| 252 | // TODO(sameeragarwal): Add caching of local parameterization, so |
| 253 | // that they are computed just once per parameter block. |
| 254 | Matrix block1_jacobian(block1_size, block1_local_size); |
| 255 | if (local_param1 == NULL) { |
| 256 | block1_jacobian.setIdentity(); |
| 257 | } else { |
| 258 | local_param1->ComputeJacobian(parameter_block1, block1_jacobian.data()); |
| 259 | } |
| 260 | |
| 261 | Matrix block2_jacobian(block2_size, block2_local_size); |
| 262 | // Fast path if the user is requesting a diagonal block. |
| 263 | if (parameter_block1 == parameter_block2) { |
| 264 | block2_jacobian = block1_jacobian; |
| 265 | } else { |
| 266 | if (local_param2 == NULL) { |
| 267 | block2_jacobian.setIdentity(); |
| 268 | } else { |
| 269 | local_param2->ComputeJacobian(parameter_block2, block2_jacobian.data()); |
| 270 | } |
| 271 | } |
| 272 | |
| 273 | if (transpose) { |
| 274 | MatrixRef(covariance_block, block2_size, block1_size) = |
| 275 | block2_jacobian * |
| 276 | cov.block(0, offset, block1_local_size, block2_local_size).transpose() * |
| 277 | block1_jacobian.transpose(); |
| 278 | } else { |
| 279 | MatrixRef(covariance_block, block1_size, block2_size) = |
| 280 | block1_jacobian * |
| 281 | cov.block(0, offset, block1_local_size, block2_local_size) * |
| 282 | block2_jacobian.transpose(); |
| 283 | } |
| 284 | |
| 285 | return true; |
| 286 | } |
| 287 | |
| 288 | bool CovarianceImpl::GetCovarianceMatrixInTangentOrAmbientSpace( |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 289 | const std::vector<const double*>& parameters, |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 290 | bool lift_covariance_to_ambient_space, |
| 291 | double* covariance_matrix) const { |
| 292 | CHECK(is_computed_) |
| 293 | << "Covariance::GetCovarianceMatrix called before Covariance::Compute"; |
| 294 | CHECK(is_valid_) |
| 295 | << "Covariance::GetCovarianceMatrix called when Covariance::Compute " |
| 296 | << "returned false."; |
| 297 | |
| 298 | const ProblemImpl::ParameterMap& parameter_map = problem_->parameter_map(); |
| 299 | // For OpenMP compatibility we need to define these vectors in advance |
| 300 | const int num_parameters = parameters.size(); |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 301 | std::vector<int> parameter_sizes; |
| 302 | std::vector<int> cum_parameter_size; |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 303 | parameter_sizes.reserve(num_parameters); |
| 304 | cum_parameter_size.resize(num_parameters + 1); |
| 305 | cum_parameter_size[0] = 0; |
| 306 | for (int i = 0; i < num_parameters; ++i) { |
| 307 | ParameterBlock* block = |
| 308 | FindOrDie(parameter_map, const_cast<double*>(parameters[i])); |
| 309 | if (lift_covariance_to_ambient_space) { |
| 310 | parameter_sizes.push_back(block->Size()); |
| 311 | } else { |
| 312 | parameter_sizes.push_back(block->LocalSize()); |
| 313 | } |
| 314 | } |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 315 | std::partial_sum(parameter_sizes.begin(), |
| 316 | parameter_sizes.end(), |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 317 | cum_parameter_size.begin() + 1); |
| 318 | const int max_covariance_block_size = |
| 319 | *std::max_element(parameter_sizes.begin(), parameter_sizes.end()); |
| 320 | const int covariance_size = cum_parameter_size.back(); |
| 321 | |
| 322 | // Assemble the blocks in the covariance matrix. |
| 323 | MatrixRef covariance(covariance_matrix, covariance_size, covariance_size); |
| 324 | const int num_threads = options_.num_threads; |
| 325 | std::unique_ptr<double[]> workspace( |
| 326 | new double[num_threads * max_covariance_block_size * |
| 327 | max_covariance_block_size]); |
| 328 | |
| 329 | bool success = true; |
| 330 | |
| 331 | // Technically the following code is a double nested loop where |
| 332 | // i = 1:n, j = i:n. |
| 333 | int iteration_count = (num_parameters * (num_parameters + 1)) / 2; |
| 334 | problem_->context()->EnsureMinimumThreads(num_threads); |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 335 | ParallelFor(problem_->context(), |
| 336 | 0, |
| 337 | iteration_count, |
| 338 | num_threads, |
| 339 | [&](int thread_id, int k) { |
| 340 | int i, j; |
| 341 | LinearIndexToUpperTriangularIndex(k, num_parameters, &i, &j); |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 342 | |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 343 | int covariance_row_idx = cum_parameter_size[i]; |
| 344 | int covariance_col_idx = cum_parameter_size[j]; |
| 345 | int size_i = parameter_sizes[i]; |
| 346 | int size_j = parameter_sizes[j]; |
| 347 | double* covariance_block = |
| 348 | workspace.get() + thread_id * max_covariance_block_size * |
| 349 | max_covariance_block_size; |
| 350 | if (!GetCovarianceBlockInTangentOrAmbientSpace( |
| 351 | parameters[i], |
| 352 | parameters[j], |
| 353 | lift_covariance_to_ambient_space, |
| 354 | covariance_block)) { |
| 355 | success = false; |
| 356 | } |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 357 | |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 358 | covariance.block( |
| 359 | covariance_row_idx, covariance_col_idx, size_i, size_j) = |
| 360 | MatrixRef(covariance_block, size_i, size_j); |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 361 | |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 362 | if (i != j) { |
| 363 | covariance.block( |
| 364 | covariance_col_idx, covariance_row_idx, size_j, size_i) = |
| 365 | MatrixRef(covariance_block, size_i, size_j).transpose(); |
| 366 | } |
| 367 | }); |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 368 | return success; |
| 369 | } |
| 370 | |
| 371 | // Determine the sparsity pattern of the covariance matrix based on |
| 372 | // the block pairs requested by the user. |
| 373 | bool CovarianceImpl::ComputeCovarianceSparsity( |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 374 | const CovarianceBlocks& original_covariance_blocks, ProblemImpl* problem) { |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 375 | EventLogger event_logger("CovarianceImpl::ComputeCovarianceSparsity"); |
| 376 | |
| 377 | // Determine an ordering for the parameter block, by sorting the |
| 378 | // parameter blocks by their pointers. |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 379 | std::vector<double*> all_parameter_blocks; |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 380 | problem->GetParameterBlocks(&all_parameter_blocks); |
| 381 | const ProblemImpl::ParameterMap& parameter_map = problem->parameter_map(); |
| 382 | std::unordered_set<ParameterBlock*> parameter_blocks_in_use; |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 383 | std::vector<ResidualBlock*> residual_blocks; |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 384 | problem->GetResidualBlocks(&residual_blocks); |
| 385 | |
| 386 | for (int i = 0; i < residual_blocks.size(); ++i) { |
| 387 | ResidualBlock* residual_block = residual_blocks[i]; |
| 388 | parameter_blocks_in_use.insert(residual_block->parameter_blocks(), |
| 389 | residual_block->parameter_blocks() + |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 390 | residual_block->NumParameterBlocks()); |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 391 | } |
| 392 | |
| 393 | constant_parameter_blocks_.clear(); |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 394 | std::vector<double*>& active_parameter_blocks = |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 395 | evaluate_options_.parameter_blocks; |
| 396 | active_parameter_blocks.clear(); |
| 397 | for (int i = 0; i < all_parameter_blocks.size(); ++i) { |
| 398 | double* parameter_block = all_parameter_blocks[i]; |
| 399 | ParameterBlock* block = FindOrDie(parameter_map, parameter_block); |
| 400 | if (!block->IsConstant() && (parameter_blocks_in_use.count(block) > 0)) { |
| 401 | active_parameter_blocks.push_back(parameter_block); |
| 402 | } else { |
| 403 | constant_parameter_blocks_.insert(parameter_block); |
| 404 | } |
| 405 | } |
| 406 | |
| 407 | std::sort(active_parameter_blocks.begin(), active_parameter_blocks.end()); |
| 408 | |
| 409 | // Compute the number of rows. Map each parameter block to the |
| 410 | // first row corresponding to it in the covariance matrix using the |
| 411 | // ordering of parameter blocks just constructed. |
| 412 | int num_rows = 0; |
| 413 | parameter_block_to_row_index_.clear(); |
| 414 | for (int i = 0; i < active_parameter_blocks.size(); ++i) { |
| 415 | double* parameter_block = active_parameter_blocks[i]; |
| 416 | const int parameter_block_size = |
| 417 | problem->ParameterBlockLocalSize(parameter_block); |
| 418 | parameter_block_to_row_index_[parameter_block] = num_rows; |
| 419 | num_rows += parameter_block_size; |
| 420 | } |
| 421 | |
| 422 | // Compute the number of non-zeros in the covariance matrix. Along |
| 423 | // the way flip any covariance blocks which are in the lower |
| 424 | // triangular part of the matrix. |
| 425 | int num_nonzeros = 0; |
| 426 | CovarianceBlocks covariance_blocks; |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 427 | for (int i = 0; i < original_covariance_blocks.size(); ++i) { |
| 428 | const std::pair<const double*, const double*>& block_pair = |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 429 | original_covariance_blocks[i]; |
| 430 | if (constant_parameter_blocks_.count(block_pair.first) > 0 || |
| 431 | constant_parameter_blocks_.count(block_pair.second) > 0) { |
| 432 | continue; |
| 433 | } |
| 434 | |
| 435 | int index1 = FindOrDie(parameter_block_to_row_index_, block_pair.first); |
| 436 | int index2 = FindOrDie(parameter_block_to_row_index_, block_pair.second); |
| 437 | const int size1 = problem->ParameterBlockLocalSize(block_pair.first); |
| 438 | const int size2 = problem->ParameterBlockLocalSize(block_pair.second); |
| 439 | num_nonzeros += size1 * size2; |
| 440 | |
| 441 | // Make sure we are constructing a block upper triangular matrix. |
| 442 | if (index1 > index2) { |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 443 | covariance_blocks.push_back( |
| 444 | std::make_pair(block_pair.second, block_pair.first)); |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 445 | } else { |
| 446 | covariance_blocks.push_back(block_pair); |
| 447 | } |
| 448 | } |
| 449 | |
| 450 | if (covariance_blocks.size() == 0) { |
| 451 | VLOG(2) << "No non-zero covariance blocks found"; |
| 452 | covariance_matrix_.reset(NULL); |
| 453 | return true; |
| 454 | } |
| 455 | |
| 456 | // Sort the block pairs. As a consequence we get the covariance |
| 457 | // blocks as they will occur in the CompressedRowSparseMatrix that |
| 458 | // will store the covariance. |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 459 | std::sort(covariance_blocks.begin(), covariance_blocks.end()); |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 460 | |
| 461 | // Fill the sparsity pattern of the covariance matrix. |
| 462 | covariance_matrix_.reset( |
| 463 | new CompressedRowSparseMatrix(num_rows, num_rows, num_nonzeros)); |
| 464 | |
| 465 | int* rows = covariance_matrix_->mutable_rows(); |
| 466 | int* cols = covariance_matrix_->mutable_cols(); |
| 467 | |
| 468 | // Iterate over parameter blocks and in turn over the rows of the |
| 469 | // covariance matrix. For each parameter block, look in the upper |
| 470 | // triangular part of the covariance matrix to see if there are any |
| 471 | // blocks requested by the user. If this is the case then fill out a |
| 472 | // set of compressed rows corresponding to this parameter block. |
| 473 | // |
| 474 | // The key thing that makes this loop work is the fact that the |
| 475 | // row/columns of the covariance matrix are ordered by the pointer |
| 476 | // values of the parameter blocks. Thus iterating over the keys of |
| 477 | // parameter_block_to_row_index_ corresponds to iterating over the |
| 478 | // rows of the covariance matrix in order. |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 479 | int i = 0; // index into covariance_blocks. |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 480 | int cursor = 0; // index into the covariance matrix. |
| 481 | for (const auto& entry : parameter_block_to_row_index_) { |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 482 | const double* row_block = entry.first; |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 483 | const int row_block_size = problem->ParameterBlockLocalSize(row_block); |
| 484 | int row_begin = entry.second; |
| 485 | |
| 486 | // Iterate over the covariance blocks contained in this row block |
| 487 | // and count the number of columns in this row block. |
| 488 | int num_col_blocks = 0; |
| 489 | int num_columns = 0; |
| 490 | for (int j = i; j < covariance_blocks.size(); ++j, ++num_col_blocks) { |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 491 | const std::pair<const double*, const double*>& block_pair = |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 492 | covariance_blocks[j]; |
| 493 | if (block_pair.first != row_block) { |
| 494 | break; |
| 495 | } |
| 496 | num_columns += problem->ParameterBlockLocalSize(block_pair.second); |
| 497 | } |
| 498 | |
| 499 | // Fill out all the compressed rows for this parameter block. |
| 500 | for (int r = 0; r < row_block_size; ++r) { |
| 501 | rows[row_begin + r] = cursor; |
| 502 | for (int c = 0; c < num_col_blocks; ++c) { |
| 503 | const double* col_block = covariance_blocks[i + c].second; |
| 504 | const int col_block_size = problem->ParameterBlockLocalSize(col_block); |
| 505 | int col_begin = FindOrDie(parameter_block_to_row_index_, col_block); |
| 506 | for (int k = 0; k < col_block_size; ++k) { |
| 507 | cols[cursor++] = col_begin++; |
| 508 | } |
| 509 | } |
| 510 | } |
| 511 | |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 512 | i += num_col_blocks; |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 513 | } |
| 514 | |
| 515 | rows[num_rows] = cursor; |
| 516 | return true; |
| 517 | } |
| 518 | |
| 519 | bool CovarianceImpl::ComputeCovarianceValues() { |
| 520 | if (options_.algorithm_type == DENSE_SVD) { |
| 521 | return ComputeCovarianceValuesUsingDenseSVD(); |
| 522 | } |
| 523 | |
| 524 | if (options_.algorithm_type == SPARSE_QR) { |
| 525 | if (options_.sparse_linear_algebra_library_type == EIGEN_SPARSE) { |
| 526 | return ComputeCovarianceValuesUsingEigenSparseQR(); |
| 527 | } |
| 528 | |
| 529 | if (options_.sparse_linear_algebra_library_type == SUITE_SPARSE) { |
| 530 | #if !defined(CERES_NO_SUITESPARSE) |
| 531 | return ComputeCovarianceValuesUsingSuiteSparseQR(); |
| 532 | #else |
| 533 | LOG(ERROR) << "SuiteSparse is required to use the SPARSE_QR algorithm " |
| 534 | << "with " |
| 535 | << "Covariance::Options::sparse_linear_algebra_library_type " |
| 536 | << "= SUITE_SPARSE."; |
| 537 | return false; |
| 538 | #endif |
| 539 | } |
| 540 | |
| 541 | LOG(ERROR) << "Unsupported " |
| 542 | << "Covariance::Options::sparse_linear_algebra_library_type " |
| 543 | << "= " |
| 544 | << SparseLinearAlgebraLibraryTypeToString( |
| 545 | options_.sparse_linear_algebra_library_type); |
| 546 | return false; |
| 547 | } |
| 548 | |
| 549 | LOG(ERROR) << "Unsupported Covariance::Options::algorithm_type = " |
| 550 | << CovarianceAlgorithmTypeToString(options_.algorithm_type); |
| 551 | return false; |
| 552 | } |
| 553 | |
| 554 | bool CovarianceImpl::ComputeCovarianceValuesUsingSuiteSparseQR() { |
| 555 | EventLogger event_logger( |
| 556 | "CovarianceImpl::ComputeCovarianceValuesUsingSparseQR"); |
| 557 | |
| 558 | #ifndef CERES_NO_SUITESPARSE |
| 559 | if (covariance_matrix_.get() == NULL) { |
| 560 | // Nothing to do, all zeros covariance matrix. |
| 561 | return true; |
| 562 | } |
| 563 | |
| 564 | CRSMatrix jacobian; |
| 565 | problem_->Evaluate(evaluate_options_, NULL, NULL, NULL, &jacobian); |
| 566 | event_logger.AddEvent("Evaluate"); |
| 567 | |
| 568 | // Construct a compressed column form of the Jacobian. |
| 569 | const int num_rows = jacobian.num_rows; |
| 570 | const int num_cols = jacobian.num_cols; |
| 571 | const int num_nonzeros = jacobian.values.size(); |
| 572 | |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 573 | std::vector<SuiteSparse_long> transpose_rows(num_cols + 1, 0); |
| 574 | std::vector<SuiteSparse_long> transpose_cols(num_nonzeros, 0); |
| 575 | std::vector<double> transpose_values(num_nonzeros, 0); |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 576 | |
| 577 | for (int idx = 0; idx < num_nonzeros; ++idx) { |
| 578 | transpose_rows[jacobian.cols[idx] + 1] += 1; |
| 579 | } |
| 580 | |
| 581 | for (int i = 1; i < transpose_rows.size(); ++i) { |
| 582 | transpose_rows[i] += transpose_rows[i - 1]; |
| 583 | } |
| 584 | |
| 585 | for (int r = 0; r < num_rows; ++r) { |
| 586 | for (int idx = jacobian.rows[r]; idx < jacobian.rows[r + 1]; ++idx) { |
| 587 | const int c = jacobian.cols[idx]; |
| 588 | const int transpose_idx = transpose_rows[c]; |
| 589 | transpose_cols[transpose_idx] = r; |
| 590 | transpose_values[transpose_idx] = jacobian.values[idx]; |
| 591 | ++transpose_rows[c]; |
| 592 | } |
| 593 | } |
| 594 | |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 595 | for (int i = transpose_rows.size() - 1; i > 0; --i) { |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 596 | transpose_rows[i] = transpose_rows[i - 1]; |
| 597 | } |
| 598 | transpose_rows[0] = 0; |
| 599 | |
| 600 | cholmod_sparse cholmod_jacobian; |
| 601 | cholmod_jacobian.nrow = num_rows; |
| 602 | cholmod_jacobian.ncol = num_cols; |
| 603 | cholmod_jacobian.nzmax = num_nonzeros; |
| 604 | cholmod_jacobian.nz = NULL; |
| 605 | cholmod_jacobian.p = reinterpret_cast<void*>(&transpose_rows[0]); |
| 606 | cholmod_jacobian.i = reinterpret_cast<void*>(&transpose_cols[0]); |
| 607 | cholmod_jacobian.x = reinterpret_cast<void*>(&transpose_values[0]); |
| 608 | cholmod_jacobian.z = NULL; |
| 609 | cholmod_jacobian.stype = 0; // Matrix is not symmetric. |
| 610 | cholmod_jacobian.itype = CHOLMOD_LONG; |
| 611 | cholmod_jacobian.xtype = CHOLMOD_REAL; |
| 612 | cholmod_jacobian.dtype = CHOLMOD_DOUBLE; |
| 613 | cholmod_jacobian.sorted = 1; |
| 614 | cholmod_jacobian.packed = 1; |
| 615 | |
| 616 | cholmod_common cc; |
| 617 | cholmod_l_start(&cc); |
| 618 | |
| 619 | cholmod_sparse* R = NULL; |
| 620 | SuiteSparse_long* permutation = NULL; |
| 621 | |
| 622 | // Compute a Q-less QR factorization of the Jacobian. Since we are |
| 623 | // only interested in inverting J'J = R'R, we do not need Q. This |
| 624 | // saves memory and gives us R as a permuted compressed column |
| 625 | // sparse matrix. |
| 626 | // |
| 627 | // TODO(sameeragarwal): Currently the symbolic factorization and the |
| 628 | // numeric factorization is done at the same time, and this does not |
| 629 | // explicitly account for the block column and row structure in the |
| 630 | // matrix. When using AMD, we have observed in the past that |
| 631 | // computing the ordering with the block matrix is significantly |
| 632 | // more efficient, both in runtime as well as the quality of |
| 633 | // ordering computed. So, it maybe worth doing that analysis |
| 634 | // separately. |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 635 | const SuiteSparse_long rank = SuiteSparseQR<double>(SPQR_ORDERING_BESTAMD, |
| 636 | SPQR_DEFAULT_TOL, |
| 637 | cholmod_jacobian.ncol, |
| 638 | &cholmod_jacobian, |
| 639 | &R, |
| 640 | &permutation, |
| 641 | &cc); |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 642 | event_logger.AddEvent("Numeric Factorization"); |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 643 | if (R == nullptr) { |
| 644 | LOG(ERROR) << "Something is wrong. SuiteSparseQR returned R = nullptr."; |
| 645 | free(permutation); |
| 646 | cholmod_l_finish(&cc); |
| 647 | return false; |
| 648 | } |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 649 | |
| 650 | if (rank < cholmod_jacobian.ncol) { |
| 651 | LOG(ERROR) << "Jacobian matrix is rank deficient. " |
| 652 | << "Number of columns: " << cholmod_jacobian.ncol |
| 653 | << " rank: " << rank; |
| 654 | free(permutation); |
| 655 | cholmod_l_free_sparse(&R, &cc); |
| 656 | cholmod_l_finish(&cc); |
| 657 | return false; |
| 658 | } |
| 659 | |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 660 | std::vector<int> inverse_permutation(num_cols); |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 661 | if (permutation) { |
| 662 | for (SuiteSparse_long i = 0; i < num_cols; ++i) { |
| 663 | inverse_permutation[permutation[i]] = i; |
| 664 | } |
| 665 | } else { |
| 666 | for (SuiteSparse_long i = 0; i < num_cols; ++i) { |
| 667 | inverse_permutation[i] = i; |
| 668 | } |
| 669 | } |
| 670 | |
| 671 | const int* rows = covariance_matrix_->rows(); |
| 672 | const int* cols = covariance_matrix_->cols(); |
| 673 | double* values = covariance_matrix_->mutable_values(); |
| 674 | |
| 675 | // The following loop exploits the fact that the i^th column of A^{-1} |
| 676 | // is given by the solution to the linear system |
| 677 | // |
| 678 | // A x = e_i |
| 679 | // |
| 680 | // where e_i is a vector with e(i) = 1 and all other entries zero. |
| 681 | // |
| 682 | // Since the covariance matrix is symmetric, the i^th row and column |
| 683 | // are equal. |
| 684 | const int num_threads = options_.num_threads; |
| 685 | std::unique_ptr<double[]> workspace(new double[num_threads * num_cols]); |
| 686 | |
| 687 | problem_->context()->EnsureMinimumThreads(num_threads); |
| 688 | ParallelFor( |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 689 | problem_->context(), 0, num_cols, num_threads, [&](int thread_id, int r) { |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 690 | const int row_begin = rows[r]; |
| 691 | const int row_end = rows[r + 1]; |
| 692 | if (row_end != row_begin) { |
| 693 | double* solution = workspace.get() + thread_id * num_cols; |
| 694 | SolveRTRWithSparseRHS<SuiteSparse_long>( |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 695 | num_cols, |
| 696 | static_cast<SuiteSparse_long*>(R->i), |
| 697 | static_cast<SuiteSparse_long*>(R->p), |
| 698 | static_cast<double*>(R->x), |
| 699 | inverse_permutation[r], |
| 700 | solution); |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 701 | for (int idx = row_begin; idx < row_end; ++idx) { |
| 702 | const int c = cols[idx]; |
| 703 | values[idx] = solution[inverse_permutation[c]]; |
| 704 | } |
| 705 | } |
| 706 | }); |
| 707 | |
| 708 | free(permutation); |
| 709 | cholmod_l_free_sparse(&R, &cc); |
| 710 | cholmod_l_finish(&cc); |
| 711 | event_logger.AddEvent("Inversion"); |
| 712 | return true; |
| 713 | |
| 714 | #else // CERES_NO_SUITESPARSE |
| 715 | |
| 716 | return false; |
| 717 | |
| 718 | #endif // CERES_NO_SUITESPARSE |
| 719 | } |
| 720 | |
| 721 | bool CovarianceImpl::ComputeCovarianceValuesUsingDenseSVD() { |
| 722 | EventLogger event_logger( |
| 723 | "CovarianceImpl::ComputeCovarianceValuesUsingDenseSVD"); |
| 724 | if (covariance_matrix_.get() == NULL) { |
| 725 | // Nothing to do, all zeros covariance matrix. |
| 726 | return true; |
| 727 | } |
| 728 | |
| 729 | CRSMatrix jacobian; |
| 730 | problem_->Evaluate(evaluate_options_, NULL, NULL, NULL, &jacobian); |
| 731 | event_logger.AddEvent("Evaluate"); |
| 732 | |
| 733 | Matrix dense_jacobian(jacobian.num_rows, jacobian.num_cols); |
| 734 | dense_jacobian.setZero(); |
| 735 | for (int r = 0; r < jacobian.num_rows; ++r) { |
| 736 | for (int idx = jacobian.rows[r]; idx < jacobian.rows[r + 1]; ++idx) { |
| 737 | const int c = jacobian.cols[idx]; |
| 738 | dense_jacobian(r, c) = jacobian.values[idx]; |
| 739 | } |
| 740 | } |
| 741 | event_logger.AddEvent("ConvertToDenseMatrix"); |
| 742 | |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 743 | Eigen::BDCSVD<Matrix> svd(dense_jacobian, |
| 744 | Eigen::ComputeThinU | Eigen::ComputeThinV); |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 745 | |
| 746 | event_logger.AddEvent("SingularValueDecomposition"); |
| 747 | |
| 748 | const Vector singular_values = svd.singularValues(); |
| 749 | const int num_singular_values = singular_values.rows(); |
| 750 | Vector inverse_squared_singular_values(num_singular_values); |
| 751 | inverse_squared_singular_values.setZero(); |
| 752 | |
| 753 | const double max_singular_value = singular_values[0]; |
| 754 | const double min_singular_value_ratio = |
| 755 | sqrt(options_.min_reciprocal_condition_number); |
| 756 | |
| 757 | const bool automatic_truncation = (options_.null_space_rank < 0); |
| 758 | const int max_rank = std::min(num_singular_values, |
| 759 | num_singular_values - options_.null_space_rank); |
| 760 | |
| 761 | // Compute the squared inverse of the singular values. Truncate the |
| 762 | // computation based on min_singular_value_ratio and |
| 763 | // null_space_rank. When either of these two quantities are active, |
| 764 | // the resulting covariance matrix is a Moore-Penrose inverse |
| 765 | // instead of a regular inverse. |
| 766 | for (int i = 0; i < max_rank; ++i) { |
| 767 | const double singular_value_ratio = singular_values[i] / max_singular_value; |
| 768 | if (singular_value_ratio < min_singular_value_ratio) { |
| 769 | // Since the singular values are in decreasing order, if |
| 770 | // automatic truncation is enabled, then from this point on |
| 771 | // all values will fail the ratio test and there is nothing to |
| 772 | // do in this loop. |
| 773 | if (automatic_truncation) { |
| 774 | break; |
| 775 | } else { |
| 776 | LOG(ERROR) << "Error: Covariance matrix is near rank deficient " |
| 777 | << "and the user did not specify a non-zero" |
| 778 | << "Covariance::Options::null_space_rank " |
| 779 | << "to enable the computation of a Pseudo-Inverse. " |
| 780 | << "Reciprocal condition number: " |
| 781 | << singular_value_ratio * singular_value_ratio << " " |
| 782 | << "min_reciprocal_condition_number: " |
| 783 | << options_.min_reciprocal_condition_number; |
| 784 | return false; |
| 785 | } |
| 786 | } |
| 787 | |
| 788 | inverse_squared_singular_values[i] = |
| 789 | 1.0 / (singular_values[i] * singular_values[i]); |
| 790 | } |
| 791 | |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 792 | Matrix dense_covariance = svd.matrixV() * |
| 793 | inverse_squared_singular_values.asDiagonal() * |
| 794 | svd.matrixV().transpose(); |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 795 | event_logger.AddEvent("PseudoInverse"); |
| 796 | |
| 797 | const int num_rows = covariance_matrix_->num_rows(); |
| 798 | const int* rows = covariance_matrix_->rows(); |
| 799 | const int* cols = covariance_matrix_->cols(); |
| 800 | double* values = covariance_matrix_->mutable_values(); |
| 801 | |
| 802 | for (int r = 0; r < num_rows; ++r) { |
| 803 | for (int idx = rows[r]; idx < rows[r + 1]; ++idx) { |
| 804 | const int c = cols[idx]; |
| 805 | values[idx] = dense_covariance(r, c); |
| 806 | } |
| 807 | } |
| 808 | event_logger.AddEvent("CopyToCovarianceMatrix"); |
| 809 | return true; |
| 810 | } |
| 811 | |
| 812 | bool CovarianceImpl::ComputeCovarianceValuesUsingEigenSparseQR() { |
| 813 | EventLogger event_logger( |
| 814 | "CovarianceImpl::ComputeCovarianceValuesUsingEigenSparseQR"); |
| 815 | if (covariance_matrix_.get() == NULL) { |
| 816 | // Nothing to do, all zeros covariance matrix. |
| 817 | return true; |
| 818 | } |
| 819 | |
| 820 | CRSMatrix jacobian; |
| 821 | problem_->Evaluate(evaluate_options_, NULL, NULL, NULL, &jacobian); |
| 822 | event_logger.AddEvent("Evaluate"); |
| 823 | |
| 824 | typedef Eigen::SparseMatrix<double, Eigen::ColMajor> EigenSparseMatrix; |
| 825 | |
| 826 | // Convert the matrix to column major order as required by SparseQR. |
| 827 | EigenSparseMatrix sparse_jacobian = |
| 828 | Eigen::MappedSparseMatrix<double, Eigen::RowMajor>( |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 829 | jacobian.num_rows, |
| 830 | jacobian.num_cols, |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 831 | static_cast<int>(jacobian.values.size()), |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 832 | jacobian.rows.data(), |
| 833 | jacobian.cols.data(), |
| 834 | jacobian.values.data()); |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 835 | event_logger.AddEvent("ConvertToSparseMatrix"); |
| 836 | |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 837 | Eigen::SparseQR<EigenSparseMatrix, Eigen::COLAMDOrdering<int>> qr_solver( |
| 838 | sparse_jacobian); |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 839 | event_logger.AddEvent("QRDecomposition"); |
| 840 | |
| 841 | if (qr_solver.info() != Eigen::Success) { |
| 842 | LOG(ERROR) << "Eigen::SparseQR decomposition failed."; |
| 843 | return false; |
| 844 | } |
| 845 | |
| 846 | if (qr_solver.rank() < jacobian.num_cols) { |
| 847 | LOG(ERROR) << "Jacobian matrix is rank deficient. " |
| 848 | << "Number of columns: " << jacobian.num_cols |
| 849 | << " rank: " << qr_solver.rank(); |
| 850 | return false; |
| 851 | } |
| 852 | |
| 853 | const int* rows = covariance_matrix_->rows(); |
| 854 | const int* cols = covariance_matrix_->cols(); |
| 855 | double* values = covariance_matrix_->mutable_values(); |
| 856 | |
| 857 | // Compute the inverse column permutation used by QR factorization. |
| 858 | Eigen::PermutationMatrix<Eigen::Dynamic, Eigen::Dynamic> inverse_permutation = |
| 859 | qr_solver.colsPermutation().inverse(); |
| 860 | |
| 861 | // The following loop exploits the fact that the i^th column of A^{-1} |
| 862 | // is given by the solution to the linear system |
| 863 | // |
| 864 | // A x = e_i |
| 865 | // |
| 866 | // where e_i is a vector with e(i) = 1 and all other entries zero. |
| 867 | // |
| 868 | // Since the covariance matrix is symmetric, the i^th row and column |
| 869 | // are equal. |
| 870 | const int num_cols = jacobian.num_cols; |
| 871 | const int num_threads = options_.num_threads; |
| 872 | std::unique_ptr<double[]> workspace(new double[num_threads * num_cols]); |
| 873 | |
| 874 | problem_->context()->EnsureMinimumThreads(num_threads); |
| 875 | ParallelFor( |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 876 | problem_->context(), 0, num_cols, num_threads, [&](int thread_id, int r) { |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 877 | const int row_begin = rows[r]; |
| 878 | const int row_end = rows[r + 1]; |
| 879 | if (row_end != row_begin) { |
| 880 | double* solution = workspace.get() + thread_id * num_cols; |
Austin Schuh | 1d1e6ea | 2020-12-23 21:56:30 -0800 | [diff] [blame^] | 881 | SolveRTRWithSparseRHS<int>(num_cols, |
| 882 | qr_solver.matrixR().innerIndexPtr(), |
| 883 | qr_solver.matrixR().outerIndexPtr(), |
| 884 | &qr_solver.matrixR().data().value(0), |
| 885 | inverse_permutation.indices().coeff(r), |
| 886 | solution); |
Austin Schuh | 70cc955 | 2019-01-21 19:46:48 -0800 | [diff] [blame] | 887 | |
| 888 | // Assign the values of the computed covariance using the |
| 889 | // inverse permutation used in the QR factorization. |
| 890 | for (int idx = row_begin; idx < row_end; ++idx) { |
| 891 | const int c = cols[idx]; |
| 892 | values[idx] = solution[inverse_permutation.indices().coeff(c)]; |
| 893 | } |
| 894 | } |
| 895 | }); |
| 896 | |
| 897 | event_logger.AddEvent("Inverse"); |
| 898 | |
| 899 | return true; |
| 900 | } |
| 901 | |
| 902 | } // namespace internal |
| 903 | } // namespace ceres |