Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 1 | // Ceres Solver - A fast non-linear least squares minimizer |
| 2 | // Copyright 2013 Google Inc. All rights reserved. |
| 3 | // http://code.google.com/p/ceres-solver/ |
| 4 | // |
| 5 | // Redistribution and use in source and binary forms, with or without |
| 6 | // modification, are permitted provided that the following conditions are met: |
| 7 | // |
| 8 | // * Redistributions of source code must retain the above copyright notice, |
| 9 | // this list of conditions and the following disclaimer. |
| 10 | // * Redistributions in binary form must reproduce the above copyright notice, |
| 11 | // this list of conditions and the following disclaimer in the documentation |
| 12 | // and/or other materials provided with the distribution. |
| 13 | // * Neither the name of Google Inc. nor the names of its contributors may be |
| 14 | // used to endorse or promote products derived from this software without |
| 15 | // specific prior written permission. |
| 16 | // |
| 17 | // THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" |
| 18 | // AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE |
| 19 | // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE |
| 20 | // ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE |
| 21 | // LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR |
| 22 | // CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF |
| 23 | // SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS |
| 24 | // INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN |
| 25 | // CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) |
| 26 | // ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE |
| 27 | // POSSIBILITY OF SUCH DAMAGE. |
| 28 | // |
| 29 | // Author: sameeragarwal@google.com (Sameer Agarwal) |
| 30 | |
| 31 | #include "ceres/covariance_impl.h" |
| 32 | |
Sameer Agarwal | 8f7e896 | 2013-06-03 13:07:39 -0700 | [diff] [blame] | 33 | #ifdef CERES_USE_OPENMP |
| 34 | #include <omp.h> |
| 35 | #endif |
| 36 | |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 37 | #include <algorithm> |
Sameer Agarwal | 29310f0 | 2013-09-23 10:59:45 -0700 | [diff] [blame] | 38 | #include <cstdlib> |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 39 | #include <utility> |
| 40 | #include <vector> |
| 41 | #include "Eigen/SVD" |
Sameer Agarwal | b22d063 | 2013-08-15 22:55:23 -0700 | [diff] [blame] | 42 | #include "ceres/compressed_col_sparse_matrix_utils.h" |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 43 | #include "ceres/compressed_row_sparse_matrix.h" |
| 44 | #include "ceres/covariance.h" |
| 45 | #include "ceres/crs_matrix.h" |
| 46 | #include "ceres/internal/eigen.h" |
| 47 | #include "ceres/map_util.h" |
| 48 | #include "ceres/parameter_block.h" |
| 49 | #include "ceres/problem_impl.h" |
| 50 | #include "ceres/suitesparse.h" |
Sameer Agarwal | df01256 | 2013-05-21 10:59:54 -0700 | [diff] [blame] | 51 | #include "ceres/wall_time.h" |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 52 | #include "glog/logging.h" |
Sergey Sharybin | fb465a0 | 2013-08-05 22:35:14 -0700 | [diff] [blame] | 53 | |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 54 | namespace ceres { |
| 55 | namespace internal { |
Sameer Agarwal | 8f7e896 | 2013-06-03 13:07:39 -0700 | [diff] [blame] | 56 | namespace { |
| 57 | |
| 58 | // Per thread storage for SuiteSparse. |
Sergey Sharybin | 8b26cc7 | 2013-06-19 21:04:02 +0600 | [diff] [blame] | 59 | #ifndef CERES_NO_SUITESPARSE |
Sameer Agarwal | b22d063 | 2013-08-15 22:55:23 -0700 | [diff] [blame] | 60 | |
Sameer Agarwal | 8f7e896 | 2013-06-03 13:07:39 -0700 | [diff] [blame] | 61 | struct PerThreadContext { |
Sameer Agarwal | d48feb8 | 2013-06-04 16:47:50 -0700 | [diff] [blame] | 62 | explicit PerThreadContext(int num_rows) |
Sameer Agarwal | 8f7e896 | 2013-06-03 13:07:39 -0700 | [diff] [blame] | 63 | : solution(NULL), |
| 64 | solution_set(NULL), |
| 65 | y_workspace(NULL), |
| 66 | e_workspace(NULL), |
| 67 | rhs(NULL) { |
| 68 | rhs = ss.CreateDenseVector(NULL, num_rows, num_rows); |
| 69 | } |
| 70 | |
| 71 | ~PerThreadContext() { |
| 72 | ss.Free(solution); |
| 73 | ss.Free(solution_set); |
| 74 | ss.Free(y_workspace); |
| 75 | ss.Free(e_workspace); |
| 76 | ss.Free(rhs); |
| 77 | } |
| 78 | |
| 79 | cholmod_dense* solution; |
| 80 | cholmod_sparse* solution_set; |
| 81 | cholmod_dense* y_workspace; |
| 82 | cholmod_dense* e_workspace; |
| 83 | cholmod_dense* rhs; |
| 84 | SuiteSparse ss; |
| 85 | }; |
Sameer Agarwal | b22d063 | 2013-08-15 22:55:23 -0700 | [diff] [blame] | 86 | |
Sergey Sharybin | 8b26cc7 | 2013-06-19 21:04:02 +0600 | [diff] [blame] | 87 | #endif |
Sameer Agarwal | 8f7e896 | 2013-06-03 13:07:39 -0700 | [diff] [blame] | 88 | |
Sameer Agarwal | d48feb8 | 2013-06-04 16:47:50 -0700 | [diff] [blame] | 89 | } // namespace |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 90 | |
| 91 | typedef vector<pair<const double*, const double*> > CovarianceBlocks; |
| 92 | |
| 93 | CovarianceImpl::CovarianceImpl(const Covariance::Options& options) |
Sameer Agarwal | df01256 | 2013-05-21 10:59:54 -0700 | [diff] [blame] | 94 | : options_(options), |
| 95 | is_computed_(false), |
| 96 | is_valid_(false) { |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 97 | evaluate_options_.num_threads = options.num_threads; |
| 98 | evaluate_options_.apply_loss_function = options.apply_loss_function; |
| 99 | } |
| 100 | |
| 101 | CovarianceImpl::~CovarianceImpl() { |
| 102 | } |
| 103 | |
| 104 | bool CovarianceImpl::Compute(const CovarianceBlocks& covariance_blocks, |
| 105 | ProblemImpl* problem) { |
| 106 | problem_ = problem; |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 107 | parameter_block_to_row_index_.clear(); |
| 108 | covariance_matrix_.reset(NULL); |
Sameer Agarwal | df01256 | 2013-05-21 10:59:54 -0700 | [diff] [blame] | 109 | is_valid_ = (ComputeCovarianceSparsity(covariance_blocks, problem) && |
| 110 | ComputeCovarianceValues()); |
| 111 | is_computed_ = true; |
| 112 | return is_valid_; |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 113 | } |
| 114 | |
| 115 | bool CovarianceImpl::GetCovarianceBlock(const double* original_parameter_block1, |
| 116 | const double* original_parameter_block2, |
| 117 | double* covariance_block) const { |
Sameer Agarwal | df01256 | 2013-05-21 10:59:54 -0700 | [diff] [blame] | 118 | CHECK(is_computed_) |
| 119 | << "Covariance::GetCovarianceBlock called before Covariance::Compute"; |
| 120 | CHECK(is_valid_) |
| 121 | << "Covariance::GetCovarianceBlock called when Covariance::Compute " |
| 122 | << "returned false."; |
| 123 | |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 124 | // If either of the two parameter blocks is constant, then the |
| 125 | // covariance block is also zero. |
| 126 | if (constant_parameter_blocks_.count(original_parameter_block1) > 0 || |
| 127 | constant_parameter_blocks_.count(original_parameter_block2) > 0) { |
| 128 | const ProblemImpl::ParameterMap& parameter_map = problem_->parameter_map(); |
| 129 | ParameterBlock* block1 = |
| 130 | FindOrDie(parameter_map, |
| 131 | const_cast<double*>(original_parameter_block1)); |
| 132 | |
| 133 | ParameterBlock* block2 = |
| 134 | FindOrDie(parameter_map, |
| 135 | const_cast<double*>(original_parameter_block2)); |
| 136 | const int block1_size = block1->Size(); |
| 137 | const int block2_size = block2->Size(); |
| 138 | MatrixRef(covariance_block, block1_size, block2_size).setZero(); |
| 139 | return true; |
| 140 | } |
| 141 | |
| 142 | const double* parameter_block1 = original_parameter_block1; |
| 143 | const double* parameter_block2 = original_parameter_block2; |
| 144 | const bool transpose = parameter_block1 > parameter_block2; |
| 145 | if (transpose) { |
| 146 | std::swap(parameter_block1, parameter_block2); |
| 147 | } |
| 148 | |
| 149 | // Find where in the covariance matrix the block is located. |
| 150 | const int row_begin = |
| 151 | FindOrDie(parameter_block_to_row_index_, parameter_block1); |
| 152 | const int col_begin = |
| 153 | FindOrDie(parameter_block_to_row_index_, parameter_block2); |
| 154 | const int* rows = covariance_matrix_->rows(); |
| 155 | const int* cols = covariance_matrix_->cols(); |
| 156 | const int row_size = rows[row_begin + 1] - rows[row_begin]; |
| 157 | const int* cols_begin = cols + rows[row_begin]; |
| 158 | |
| 159 | // The only part that requires work is walking the compressed column |
| 160 | // vector to determine where the set of columns correspnding to the |
| 161 | // covariance block begin. |
| 162 | int offset = 0; |
| 163 | while (cols_begin[offset] != col_begin && offset < row_size) { |
| 164 | ++offset; |
| 165 | } |
| 166 | |
| 167 | if (offset == row_size) { |
Sameer Agarwal | b16e118 | 2013-11-25 05:47:43 -0800 | [diff] [blame] | 168 | LOG(ERROR) << "Unable to find covariance block for " |
| 169 | << original_parameter_block1 << " " |
| 170 | << original_parameter_block2; |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 171 | return false; |
| 172 | } |
| 173 | |
| 174 | const ProblemImpl::ParameterMap& parameter_map = problem_->parameter_map(); |
| 175 | ParameterBlock* block1 = |
| 176 | FindOrDie(parameter_map, const_cast<double*>(parameter_block1)); |
| 177 | ParameterBlock* block2 = |
| 178 | FindOrDie(parameter_map, const_cast<double*>(parameter_block2)); |
| 179 | const LocalParameterization* local_param1 = block1->local_parameterization(); |
| 180 | const LocalParameterization* local_param2 = block2->local_parameterization(); |
| 181 | const int block1_size = block1->Size(); |
| 182 | const int block1_local_size = block1->LocalSize(); |
| 183 | const int block2_size = block2->Size(); |
| 184 | const int block2_local_size = block2->LocalSize(); |
| 185 | |
| 186 | ConstMatrixRef cov(covariance_matrix_->values() + rows[row_begin], |
| 187 | block1_size, |
| 188 | row_size); |
| 189 | |
| 190 | // Fast path when there are no local parameterizations. |
| 191 | if (local_param1 == NULL && local_param2 == NULL) { |
| 192 | if (transpose) { |
| 193 | MatrixRef(covariance_block, block2_size, block1_size) = |
| 194 | cov.block(0, offset, block1_size, block2_size).transpose(); |
| 195 | } else { |
| 196 | MatrixRef(covariance_block, block1_size, block2_size) = |
| 197 | cov.block(0, offset, block1_size, block2_size); |
| 198 | } |
| 199 | return true; |
| 200 | } |
| 201 | |
| 202 | // If local parameterizations are used then the covariance that has |
| 203 | // been computed is in the tangent space and it needs to be lifted |
| 204 | // back to the ambient space. |
| 205 | // |
| 206 | // This is given by the formula |
| 207 | // |
| 208 | // C'_12 = J_1 C_12 J_2' |
| 209 | // |
| 210 | // Where C_12 is the local tangent space covariance for parameter |
| 211 | // blocks 1 and 2. J_1 and J_2 are respectively the local to global |
| 212 | // jacobians for parameter blocks 1 and 2. |
| 213 | // |
| 214 | // See Result 5.11 on page 142 of Hartley & Zisserman (2nd Edition) |
| 215 | // for a proof. |
| 216 | // |
| 217 | // TODO(sameeragarwal): Add caching of local parameterization, so |
| 218 | // that they are computed just once per parameter block. |
| 219 | Matrix block1_jacobian(block1_size, block1_local_size); |
| 220 | if (local_param1 == NULL) { |
| 221 | block1_jacobian.setIdentity(); |
| 222 | } else { |
| 223 | local_param1->ComputeJacobian(parameter_block1, block1_jacobian.data()); |
| 224 | } |
| 225 | |
| 226 | Matrix block2_jacobian(block2_size, block2_local_size); |
| 227 | // Fast path if the user is requesting a diagonal block. |
| 228 | if (parameter_block1 == parameter_block2) { |
| 229 | block2_jacobian = block1_jacobian; |
| 230 | } else { |
| 231 | if (local_param2 == NULL) { |
| 232 | block2_jacobian.setIdentity(); |
| 233 | } else { |
| 234 | local_param2->ComputeJacobian(parameter_block2, block2_jacobian.data()); |
| 235 | } |
| 236 | } |
| 237 | |
| 238 | if (transpose) { |
| 239 | MatrixRef(covariance_block, block2_size, block1_size) = |
| 240 | block2_jacobian * |
| 241 | cov.block(0, offset, block1_local_size, block2_local_size).transpose() * |
| 242 | block1_jacobian.transpose(); |
| 243 | } else { |
| 244 | MatrixRef(covariance_block, block1_size, block2_size) = |
| 245 | block1_jacobian * |
| 246 | cov.block(0, offset, block1_local_size, block2_local_size) * |
| 247 | block2_jacobian.transpose(); |
| 248 | } |
| 249 | |
| 250 | return true; |
| 251 | } |
| 252 | |
| 253 | // Determine the sparsity pattern of the covariance matrix based on |
| 254 | // the block pairs requested by the user. |
| 255 | bool CovarianceImpl::ComputeCovarianceSparsity( |
| 256 | const CovarianceBlocks& original_covariance_blocks, |
| 257 | ProblemImpl* problem) { |
Sameer Agarwal | df01256 | 2013-05-21 10:59:54 -0700 | [diff] [blame] | 258 | EventLogger event_logger("CovarianceImpl::ComputeCovarianceSparsity"); |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 259 | |
| 260 | // Determine an ordering for the parameter block, by sorting the |
| 261 | // parameter blocks by their pointers. |
| 262 | vector<double*> all_parameter_blocks; |
| 263 | problem->GetParameterBlocks(&all_parameter_blocks); |
| 264 | const ProblemImpl::ParameterMap& parameter_map = problem->parameter_map(); |
| 265 | constant_parameter_blocks_.clear(); |
| 266 | vector<double*>& active_parameter_blocks = evaluate_options_.parameter_blocks; |
| 267 | active_parameter_blocks.clear(); |
| 268 | for (int i = 0; i < all_parameter_blocks.size(); ++i) { |
| 269 | double* parameter_block = all_parameter_blocks[i]; |
| 270 | |
| 271 | ParameterBlock* block = FindOrDie(parameter_map, parameter_block); |
| 272 | if (block->IsConstant()) { |
| 273 | constant_parameter_blocks_.insert(parameter_block); |
| 274 | } else { |
| 275 | active_parameter_blocks.push_back(parameter_block); |
| 276 | } |
| 277 | } |
| 278 | |
| 279 | sort(active_parameter_blocks.begin(), active_parameter_blocks.end()); |
| 280 | |
| 281 | // Compute the number of rows. Map each parameter block to the |
| 282 | // first row corresponding to it in the covariance matrix using the |
| 283 | // ordering of parameter blocks just constructed. |
| 284 | int num_rows = 0; |
| 285 | parameter_block_to_row_index_.clear(); |
| 286 | for (int i = 0; i < active_parameter_blocks.size(); ++i) { |
| 287 | double* parameter_block = active_parameter_blocks[i]; |
| 288 | const int parameter_block_size = |
| 289 | problem->ParameterBlockLocalSize(parameter_block); |
| 290 | parameter_block_to_row_index_[parameter_block] = num_rows; |
| 291 | num_rows += parameter_block_size; |
| 292 | } |
| 293 | |
| 294 | // Compute the number of non-zeros in the covariance matrix. Along |
| 295 | // the way flip any covariance blocks which are in the lower |
| 296 | // triangular part of the matrix. |
| 297 | int num_nonzeros = 0; |
| 298 | CovarianceBlocks covariance_blocks; |
| 299 | for (int i = 0; i < original_covariance_blocks.size(); ++i) { |
| 300 | const pair<const double*, const double*>& block_pair = |
| 301 | original_covariance_blocks[i]; |
| 302 | if (constant_parameter_blocks_.count(block_pair.first) > 0 || |
| 303 | constant_parameter_blocks_.count(block_pair.second) > 0) { |
| 304 | continue; |
| 305 | } |
| 306 | |
| 307 | int index1 = FindOrDie(parameter_block_to_row_index_, block_pair.first); |
| 308 | int index2 = FindOrDie(parameter_block_to_row_index_, block_pair.second); |
| 309 | const int size1 = problem->ParameterBlockLocalSize(block_pair.first); |
| 310 | const int size2 = problem->ParameterBlockLocalSize(block_pair.second); |
| 311 | num_nonzeros += size1 * size2; |
| 312 | |
| 313 | // Make sure we are constructing a block upper triangular matrix. |
| 314 | if (index1 > index2) { |
| 315 | covariance_blocks.push_back(make_pair(block_pair.second, |
| 316 | block_pair.first)); |
| 317 | } else { |
| 318 | covariance_blocks.push_back(block_pair); |
| 319 | } |
| 320 | } |
| 321 | |
| 322 | if (covariance_blocks.size() == 0) { |
| 323 | VLOG(2) << "No non-zero covariance blocks found"; |
| 324 | covariance_matrix_.reset(NULL); |
| 325 | return true; |
| 326 | } |
| 327 | |
| 328 | // Sort the block pairs. As a consequence we get the covariance |
| 329 | // blocks as they will occur in the CompressedRowSparseMatrix that |
| 330 | // will store the covariance. |
| 331 | sort(covariance_blocks.begin(), covariance_blocks.end()); |
| 332 | |
| 333 | // Fill the sparsity pattern of the covariance matrix. |
| 334 | covariance_matrix_.reset( |
| 335 | new CompressedRowSparseMatrix(num_rows, num_rows, num_nonzeros)); |
| 336 | |
| 337 | int* rows = covariance_matrix_->mutable_rows(); |
| 338 | int* cols = covariance_matrix_->mutable_cols(); |
| 339 | |
| 340 | // Iterate over parameter blocks and in turn over the rows of the |
| 341 | // covariance matrix. For each parameter block, look in the upper |
| 342 | // triangular part of the covariance matrix to see if there are any |
| 343 | // blocks requested by the user. If this is the case then fill out a |
| 344 | // set of compressed rows corresponding to this parameter block. |
| 345 | // |
| 346 | // The key thing that makes this loop work is the fact that the |
| 347 | // row/columns of the covariance matrix are ordered by the pointer |
| 348 | // values of the parameter blocks. Thus iterating over the keys of |
| 349 | // parameter_block_to_row_index_ corresponds to iterating over the |
| 350 | // rows of the covariance matrix in order. |
Sameer Agarwal | 89a592f | 2013-11-26 11:35:49 -0800 | [diff] [blame] | 351 | int i = 0; // index into covariance_blocks. |
| 352 | int cursor = 0; // index into the covariance matrix. |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 353 | for (map<const double*, int>::const_iterator it = |
| 354 | parameter_block_to_row_index_.begin(); |
| 355 | it != parameter_block_to_row_index_.end(); |
| 356 | ++it) { |
| 357 | const double* row_block = it->first; |
| 358 | const int row_block_size = problem->ParameterBlockLocalSize(row_block); |
| 359 | int row_begin = it->second; |
| 360 | |
| 361 | // Iterate over the covariance blocks contained in this row block |
| 362 | // and count the number of columns in this row block. |
| 363 | int num_col_blocks = 0; |
| 364 | int num_columns = 0; |
| 365 | for (int j = i; j < covariance_blocks.size(); ++j, ++num_col_blocks) { |
| 366 | const pair<const double*, const double*>& block_pair = |
| 367 | covariance_blocks[j]; |
| 368 | if (block_pair.first != row_block) { |
| 369 | break; |
| 370 | } |
| 371 | num_columns += problem->ParameterBlockLocalSize(block_pair.second); |
| 372 | } |
| 373 | |
| 374 | // Fill out all the compressed rows for this parameter block. |
| 375 | for (int r = 0; r < row_block_size; ++r) { |
| 376 | rows[row_begin + r] = cursor; |
| 377 | for (int c = 0; c < num_col_blocks; ++c) { |
| 378 | const double* col_block = covariance_blocks[i + c].second; |
| 379 | const int col_block_size = problem->ParameterBlockLocalSize(col_block); |
| 380 | int col_begin = FindOrDie(parameter_block_to_row_index_, col_block); |
| 381 | for (int k = 0; k < col_block_size; ++k) { |
| 382 | cols[cursor++] = col_begin++; |
| 383 | } |
| 384 | } |
| 385 | } |
| 386 | |
| 387 | i+= num_col_blocks; |
| 388 | } |
| 389 | |
| 390 | rows[num_rows] = cursor; |
| 391 | return true; |
| 392 | } |
| 393 | |
| 394 | bool CovarianceImpl::ComputeCovarianceValues() { |
Sameer Agarwal | 5a97471 | 2013-06-07 22:38:30 -0700 | [diff] [blame] | 395 | switch (options_.algorithm_type) { |
Sameer Agarwal | 89a592f | 2013-11-26 11:35:49 -0800 | [diff] [blame] | 396 | case DENSE_SVD: |
Sameer Agarwal | 5a97471 | 2013-06-07 22:38:30 -0700 | [diff] [blame] | 397 | return ComputeCovarianceValuesUsingDenseSVD(); |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 398 | #ifndef CERES_NO_SUITESPARSE |
Sameer Agarwal | 89a592f | 2013-11-26 11:35:49 -0800 | [diff] [blame] | 399 | case SPARSE_CHOLESKY: |
Sameer Agarwal | 5a97471 | 2013-06-07 22:38:30 -0700 | [diff] [blame] | 400 | return ComputeCovarianceValuesUsingSparseCholesky(); |
Sameer Agarwal | 89a592f | 2013-11-26 11:35:49 -0800 | [diff] [blame] | 401 | case SPARSE_QR: |
Sameer Agarwal | 5a97471 | 2013-06-07 22:38:30 -0700 | [diff] [blame] | 402 | return ComputeCovarianceValuesUsingSparseQR(); |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 403 | #endif |
Sameer Agarwal | 5a97471 | 2013-06-07 22:38:30 -0700 | [diff] [blame] | 404 | default: |
| 405 | LOG(ERROR) << "Unsupported covariance estimation algorithm type: " |
| 406 | << CovarianceAlgorithmTypeToString(options_.algorithm_type); |
| 407 | return false; |
| 408 | } |
| 409 | return false; |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 410 | } |
| 411 | |
Sameer Agarwal | 5a97471 | 2013-06-07 22:38:30 -0700 | [diff] [blame] | 412 | bool CovarianceImpl::ComputeCovarianceValuesUsingSparseCholesky() { |
Sameer Agarwal | df01256 | 2013-05-21 10:59:54 -0700 | [diff] [blame] | 413 | EventLogger event_logger( |
Sameer Agarwal | 5a97471 | 2013-06-07 22:38:30 -0700 | [diff] [blame] | 414 | "CovarianceImpl::ComputeCovarianceValuesUsingSparseCholesky"); |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 415 | #ifndef CERES_NO_SUITESPARSE |
| 416 | if (covariance_matrix_.get() == NULL) { |
| 417 | // Nothing to do, all zeros covariance matrix. |
| 418 | return true; |
| 419 | } |
| 420 | |
Sameer Agarwal | 5a97471 | 2013-06-07 22:38:30 -0700 | [diff] [blame] | 421 | SuiteSparse ss; |
| 422 | |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 423 | CRSMatrix jacobian; |
| 424 | problem_->Evaluate(evaluate_options_, NULL, NULL, NULL, &jacobian); |
| 425 | |
Sameer Agarwal | df01256 | 2013-05-21 10:59:54 -0700 | [diff] [blame] | 426 | event_logger.AddEvent("Evaluate"); |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 427 | // m is a transposed view of the Jacobian. |
| 428 | cholmod_sparse cholmod_jacobian_view; |
| 429 | cholmod_jacobian_view.nrow = jacobian.num_cols; |
| 430 | cholmod_jacobian_view.ncol = jacobian.num_rows; |
| 431 | cholmod_jacobian_view.nzmax = jacobian.values.size(); |
| 432 | cholmod_jacobian_view.nz = NULL; |
| 433 | cholmod_jacobian_view.p = reinterpret_cast<void*>(&jacobian.rows[0]); |
| 434 | cholmod_jacobian_view.i = reinterpret_cast<void*>(&jacobian.cols[0]); |
| 435 | cholmod_jacobian_view.x = reinterpret_cast<void*>(&jacobian.values[0]); |
| 436 | cholmod_jacobian_view.z = NULL; |
| 437 | cholmod_jacobian_view.stype = 0; // Matrix is not symmetric. |
| 438 | cholmod_jacobian_view.itype = CHOLMOD_INT; |
| 439 | cholmod_jacobian_view.xtype = CHOLMOD_REAL; |
| 440 | cholmod_jacobian_view.dtype = CHOLMOD_DOUBLE; |
| 441 | cholmod_jacobian_view.sorted = 1; |
| 442 | cholmod_jacobian_view.packed = 1; |
| 443 | |
Sameer Agarwal | 3b1ad31 | 2013-12-02 15:43:20 -0800 | [diff] [blame] | 444 | string message; |
| 445 | cholmod_factor* factor = ss.AnalyzeCholesky(&cholmod_jacobian_view, &message); |
Sameer Agarwal | df01256 | 2013-05-21 10:59:54 -0700 | [diff] [blame] | 446 | event_logger.AddEvent("Symbolic Factorization"); |
Sameer Agarwal | 79bde35 | 2013-11-21 21:33:51 -0800 | [diff] [blame] | 447 | if (factor == NULL) { |
Sameer Agarwal | b16e118 | 2013-11-25 05:47:43 -0800 | [diff] [blame] | 448 | LOG(ERROR) << "Covariance estimation failed. " |
| 449 | << "CHOLMOD symbolic cholesky factorization returned with: " |
Sameer Agarwal | 3b1ad31 | 2013-12-02 15:43:20 -0800 | [diff] [blame] | 450 | << message; |
Sameer Agarwal | 79bde35 | 2013-11-21 21:33:51 -0800 | [diff] [blame] | 451 | return false; |
Sameer Agarwal | 7129cd3 | 2013-06-02 19:35:29 -0700 | [diff] [blame] | 452 | } |
| 453 | |
Sameer Agarwal | 66e15b4 | 2013-11-22 07:59:23 -0800 | [diff] [blame] | 454 | LinearSolverTerminationType termination_type = |
Sameer Agarwal | 3b1ad31 | 2013-12-02 15:43:20 -0800 | [diff] [blame] | 455 | ss.Cholesky(&cholmod_jacobian_view, factor, &message); |
Sameer Agarwal | df01256 | 2013-05-21 10:59:54 -0700 | [diff] [blame] | 456 | event_logger.AddEvent("Numeric Factorization"); |
Sameer Agarwal | 33e01b9 | 2013-11-27 10:24:03 -0800 | [diff] [blame] | 457 | if (termination_type != LINEAR_SOLVER_SUCCESS) { |
Sameer Agarwal | b16e118 | 2013-11-25 05:47:43 -0800 | [diff] [blame] | 458 | LOG(ERROR) << "Covariance estimation failed. " |
| 459 | << "CHOLMOD numeric cholesky factorization returned with: " |
Sameer Agarwal | 3b1ad31 | 2013-12-02 15:43:20 -0800 | [diff] [blame] | 460 | << message; |
Sameer Agarwal | b16e118 | 2013-11-25 05:47:43 -0800 | [diff] [blame] | 461 | ss.Free(factor); |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 462 | return false; |
| 463 | } |
| 464 | |
Sameer Agarwal | 79bde35 | 2013-11-21 21:33:51 -0800 | [diff] [blame] | 465 | const double reciprocal_condition_number = |
| 466 | cholmod_rcond(factor, ss.mutable_cc()); |
| 467 | |
| 468 | if (reciprocal_condition_number < |
| 469 | options_.min_reciprocal_condition_number) { |
Sameer Agarwal | b16e118 | 2013-11-25 05:47:43 -0800 | [diff] [blame] | 470 | LOG(ERROR) << "Cholesky factorization of J'J is not reliable. " |
| 471 | << "Reciprocal condition number: " |
| 472 | << reciprocal_condition_number << " " |
Sameer Agarwal | ed92366 | 2013-11-28 06:50:43 -0800 | [diff] [blame] | 473 | << "min_reciprocal_condition_number: " |
Sameer Agarwal | b16e118 | 2013-11-25 05:47:43 -0800 | [diff] [blame] | 474 | << options_.min_reciprocal_condition_number; |
Sameer Agarwal | 79bde35 | 2013-11-21 21:33:51 -0800 | [diff] [blame] | 475 | ss.Free(factor); |
| 476 | return false; |
| 477 | } |
| 478 | |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 479 | const int num_rows = covariance_matrix_->num_rows(); |
| 480 | const int* rows = covariance_matrix_->rows(); |
| 481 | const int* cols = covariance_matrix_->cols(); |
| 482 | double* values = covariance_matrix_->mutable_values(); |
| 483 | |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 484 | // The following loop exploits the fact that the i^th column of A^{-1} |
| 485 | // is given by the solution to the linear system |
| 486 | // |
| 487 | // A x = e_i |
| 488 | // |
| 489 | // where e_i is a vector with e(i) = 1 and all other entries zero. |
| 490 | // |
| 491 | // Since the covariance matrix is symmetric, the i^th row and column |
| 492 | // are equal. |
| 493 | // |
| 494 | // The ifdef separates two different version of SuiteSparse. Newer |
| 495 | // versions of SuiteSparse have the cholmod_solve2 function which |
| 496 | // re-uses memory across calls. |
| 497 | #if (SUITESPARSE_VERSION < 4002) |
Sameer Agarwal | 5a97471 | 2013-06-07 22:38:30 -0700 | [diff] [blame] | 498 | cholmod_dense* rhs = ss.CreateDenseVector(NULL, num_rows, num_rows); |
Sameer Agarwal | 8f7e896 | 2013-06-03 13:07:39 -0700 | [diff] [blame] | 499 | double* rhs_x = reinterpret_cast<double*>(rhs->x); |
| 500 | |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 501 | for (int r = 0; r < num_rows; ++r) { |
| 502 | int row_begin = rows[r]; |
| 503 | int row_end = rows[r + 1]; |
| 504 | if (row_end == row_begin) { |
| 505 | continue; |
| 506 | } |
| 507 | |
| 508 | rhs_x[r] = 1.0; |
Sameer Agarwal | 3b1ad31 | 2013-12-02 15:43:20 -0800 | [diff] [blame] | 509 | cholmod_dense* solution = ss.Solve(factor, rhs, &message); |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 510 | double* solution_x = reinterpret_cast<double*>(solution->x); |
| 511 | for (int idx = row_begin; idx < row_end; ++idx) { |
| 512 | const int c = cols[idx]; |
| 513 | values[idx] = solution_x[c]; |
| 514 | } |
Sameer Agarwal | 5a97471 | 2013-06-07 22:38:30 -0700 | [diff] [blame] | 515 | ss.Free(solution); |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 516 | rhs_x[r] = 0.0; |
| 517 | } |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 518 | |
Sameer Agarwal | 5a97471 | 2013-06-07 22:38:30 -0700 | [diff] [blame] | 519 | ss.Free(rhs); |
Sameer Agarwal | 8f7e896 | 2013-06-03 13:07:39 -0700 | [diff] [blame] | 520 | #else // SUITESPARSE_VERSION < 4002 |
| 521 | |
| 522 | const int num_threads = options_.num_threads; |
| 523 | vector<PerThreadContext*> contexts(num_threads); |
| 524 | for (int i = 0; i < num_threads; ++i) { |
| 525 | contexts[i] = new PerThreadContext(num_rows); |
| 526 | } |
| 527 | |
| 528 | // The first call to cholmod_solve2 is not thread safe, since it |
| 529 | // changes the factorization from supernodal to simplicial etc. |
| 530 | { |
| 531 | PerThreadContext* context = contexts[0]; |
| 532 | double* context_rhs_x = reinterpret_cast<double*>(context->rhs->x); |
| 533 | context_rhs_x[0] = 1.0; |
| 534 | cholmod_solve2(CHOLMOD_A, |
| 535 | factor, |
| 536 | context->rhs, |
| 537 | NULL, |
| 538 | &context->solution, |
| 539 | &context->solution_set, |
| 540 | &context->y_workspace, |
| 541 | &context->e_workspace, |
| 542 | context->ss.mutable_cc()); |
| 543 | context_rhs_x[0] = 0.0; |
| 544 | } |
| 545 | |
| 546 | #pragma omp parallel for num_threads(num_threads) schedule(dynamic) |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 547 | for (int r = 0; r < num_rows; ++r) { |
| 548 | int row_begin = rows[r]; |
| 549 | int row_end = rows[r + 1]; |
| 550 | if (row_end == row_begin) { |
| 551 | continue; |
| 552 | } |
| 553 | |
Sameer Agarwal | 8f7e896 | 2013-06-03 13:07:39 -0700 | [diff] [blame] | 554 | # ifdef CERES_USE_OPENMP |
| 555 | int thread_id = omp_get_thread_num(); |
| 556 | # else |
| 557 | int thread_id = 0; |
| 558 | # endif |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 559 | |
Sameer Agarwal | 8f7e896 | 2013-06-03 13:07:39 -0700 | [diff] [blame] | 560 | PerThreadContext* context = contexts[thread_id]; |
| 561 | double* context_rhs_x = reinterpret_cast<double*>(context->rhs->x); |
| 562 | context_rhs_x[r] = 1.0; |
| 563 | |
| 564 | // TODO(sameeragarwal) There should be a more efficient way |
| 565 | // involving the use of Bset but I am unable to make it work right |
| 566 | // now. |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 567 | cholmod_solve2(CHOLMOD_A, |
| 568 | factor, |
Sameer Agarwal | 8f7e896 | 2013-06-03 13:07:39 -0700 | [diff] [blame] | 569 | context->rhs, |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 570 | NULL, |
Sameer Agarwal | 8f7e896 | 2013-06-03 13:07:39 -0700 | [diff] [blame] | 571 | &context->solution, |
| 572 | &context->solution_set, |
| 573 | &context->y_workspace, |
| 574 | &context->e_workspace, |
| 575 | context->ss.mutable_cc()); |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 576 | |
Sameer Agarwal | 8f7e896 | 2013-06-03 13:07:39 -0700 | [diff] [blame] | 577 | double* solution_x = reinterpret_cast<double*>(context->solution->x); |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 578 | for (int idx = row_begin; idx < row_end; ++idx) { |
| 579 | const int c = cols[idx]; |
| 580 | values[idx] = solution_x[c]; |
| 581 | } |
Sameer Agarwal | 8f7e896 | 2013-06-03 13:07:39 -0700 | [diff] [blame] | 582 | context_rhs_x[r] = 0.0; |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 583 | } |
| 584 | |
Sameer Agarwal | 8f7e896 | 2013-06-03 13:07:39 -0700 | [diff] [blame] | 585 | for (int i = 0; i < num_threads; ++i) { |
| 586 | delete contexts[i]; |
| 587 | } |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 588 | |
Sameer Agarwal | 8f7e896 | 2013-06-03 13:07:39 -0700 | [diff] [blame] | 589 | #endif // SUITESPARSE_VERSION < 4002 |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 590 | |
Sameer Agarwal | 5a97471 | 2013-06-07 22:38:30 -0700 | [diff] [blame] | 591 | ss.Free(factor); |
Sameer Agarwal | df01256 | 2013-05-21 10:59:54 -0700 | [diff] [blame] | 592 | event_logger.AddEvent("Inversion"); |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 593 | return true; |
Sameer Agarwal | 8f7e896 | 2013-06-03 13:07:39 -0700 | [diff] [blame] | 594 | |
| 595 | #else // CERES_NO_SUITESPARSE |
| 596 | |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 597 | return false; |
Sameer Agarwal | 8f7e896 | 2013-06-03 13:07:39 -0700 | [diff] [blame] | 598 | |
| 599 | #endif // CERES_NO_SUITESPARSE |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 600 | }; |
| 601 | |
Sameer Agarwal | 5a97471 | 2013-06-07 22:38:30 -0700 | [diff] [blame] | 602 | bool CovarianceImpl::ComputeCovarianceValuesUsingSparseQR() { |
Sameer Agarwal | df01256 | 2013-05-21 10:59:54 -0700 | [diff] [blame] | 603 | EventLogger event_logger( |
Sameer Agarwal | 5a97471 | 2013-06-07 22:38:30 -0700 | [diff] [blame] | 604 | "CovarianceImpl::ComputeCovarianceValuesUsingSparseQR"); |
| 605 | |
| 606 | #ifndef CERES_NO_SUITESPARSE |
| 607 | if (covariance_matrix_.get() == NULL) { |
| 608 | // Nothing to do, all zeros covariance matrix. |
| 609 | return true; |
| 610 | } |
| 611 | |
| 612 | CRSMatrix jacobian; |
| 613 | problem_->Evaluate(evaluate_options_, NULL, NULL, NULL, &jacobian); |
| 614 | event_logger.AddEvent("Evaluate"); |
| 615 | |
| 616 | // Construct a compressed column form of the Jacobian. |
| 617 | const int num_rows = jacobian.num_rows; |
| 618 | const int num_cols = jacobian.num_cols; |
| 619 | const int num_nonzeros = jacobian.values.size(); |
| 620 | |
| 621 | vector<SuiteSparse_long> transpose_rows(num_cols + 1, 0); |
| 622 | vector<SuiteSparse_long> transpose_cols(num_nonzeros, 0); |
| 623 | vector<double> transpose_values(num_nonzeros, 0); |
| 624 | |
| 625 | for (int idx = 0; idx < num_nonzeros; ++idx) { |
| 626 | transpose_rows[jacobian.cols[idx] + 1] += 1; |
| 627 | } |
| 628 | |
| 629 | for (int i = 1; i < transpose_rows.size(); ++i) { |
| 630 | transpose_rows[i] += transpose_rows[i - 1]; |
| 631 | } |
| 632 | |
| 633 | for (int r = 0; r < num_rows; ++r) { |
| 634 | for (int idx = jacobian.rows[r]; idx < jacobian.rows[r + 1]; ++idx) { |
| 635 | const int c = jacobian.cols[idx]; |
| 636 | const int transpose_idx = transpose_rows[c]; |
| 637 | transpose_cols[transpose_idx] = r; |
| 638 | transpose_values[transpose_idx] = jacobian.values[idx]; |
| 639 | ++transpose_rows[c]; |
| 640 | } |
| 641 | } |
| 642 | |
| 643 | for (int i = transpose_rows.size() - 1; i > 0 ; --i) { |
| 644 | transpose_rows[i] = transpose_rows[i - 1]; |
| 645 | } |
| 646 | transpose_rows[0] = 0; |
| 647 | |
| 648 | cholmod_sparse cholmod_jacobian; |
| 649 | cholmod_jacobian.nrow = num_rows; |
| 650 | cholmod_jacobian.ncol = num_cols; |
| 651 | cholmod_jacobian.nzmax = num_nonzeros; |
| 652 | cholmod_jacobian.nz = NULL; |
| 653 | cholmod_jacobian.p = reinterpret_cast<void*>(&transpose_rows[0]); |
| 654 | cholmod_jacobian.i = reinterpret_cast<void*>(&transpose_cols[0]); |
| 655 | cholmod_jacobian.x = reinterpret_cast<void*>(&transpose_values[0]); |
| 656 | cholmod_jacobian.z = NULL; |
| 657 | cholmod_jacobian.stype = 0; // Matrix is not symmetric. |
| 658 | cholmod_jacobian.itype = CHOLMOD_LONG; |
| 659 | cholmod_jacobian.xtype = CHOLMOD_REAL; |
| 660 | cholmod_jacobian.dtype = CHOLMOD_DOUBLE; |
| 661 | cholmod_jacobian.sorted = 1; |
| 662 | cholmod_jacobian.packed = 1; |
| 663 | |
| 664 | cholmod_common cc; |
| 665 | cholmod_l_start(&cc); |
| 666 | |
Sameer Agarwal | b22d063 | 2013-08-15 22:55:23 -0700 | [diff] [blame] | 667 | cholmod_sparse* R = NULL; |
| 668 | SuiteSparse_long* permutation = NULL; |
Sameer Agarwal | 5a97471 | 2013-06-07 22:38:30 -0700 | [diff] [blame] | 669 | |
Sameer Agarwal | b22d063 | 2013-08-15 22:55:23 -0700 | [diff] [blame] | 670 | // Compute a Q-less QR factorization of the Jacobian. Since we are |
| 671 | // only interested in inverting J'J = R'R, we do not need Q. This |
| 672 | // saves memory and gives us R as a permuted compressed column |
| 673 | // sparse matrix. |
| 674 | // |
| 675 | // TODO(sameeragarwal): Currently the symbolic factorization and the |
| 676 | // numeric factorization is done at the same time, and this does not |
| 677 | // explicitly account for the block column and row structure in the |
| 678 | // matrix. When using AMD, we have observed in the past that |
| 679 | // computing the ordering with the block matrix is significantly |
| 680 | // more efficient, both in runtime as well as the quality of |
| 681 | // ordering computed. So, it maybe worth doing that analysis |
| 682 | // separately. |
| 683 | const SuiteSparse_long rank = |
| 684 | SuiteSparseQR<double>(SPQR_ORDERING_BESTAMD, |
| 685 | SPQR_DEFAULT_TOL, |
| 686 | cholmod_jacobian.ncol, |
| 687 | &cholmod_jacobian, |
| 688 | &R, |
| 689 | &permutation, |
| 690 | &cc); |
| 691 | event_logger.AddEvent("Numeric Factorization"); |
| 692 | CHECK_NOTNULL(permutation); |
| 693 | CHECK_NOTNULL(R); |
| 694 | |
Sameer Agarwal | 5a97471 | 2013-06-07 22:38:30 -0700 | [diff] [blame] | 695 | if (rank < cholmod_jacobian.ncol) { |
Sameer Agarwal | b16e118 | 2013-11-25 05:47:43 -0800 | [diff] [blame] | 696 | LOG(ERROR) << "Jacobian matrix is rank deficient. " |
| 697 | << "Number of columns: " << cholmod_jacobian.ncol |
| 698 | << " rank: " << rank; |
Sameer Agarwal | 29310f0 | 2013-09-23 10:59:45 -0700 | [diff] [blame] | 699 | free(permutation); |
Sameer Agarwal | b22d063 | 2013-08-15 22:55:23 -0700 | [diff] [blame] | 700 | cholmod_l_free_sparse(&R, &cc); |
Sameer Agarwal | 5a97471 | 2013-06-07 22:38:30 -0700 | [diff] [blame] | 701 | cholmod_l_finish(&cc); |
| 702 | return false; |
| 703 | } |
| 704 | |
Sameer Agarwal | b22d063 | 2013-08-15 22:55:23 -0700 | [diff] [blame] | 705 | vector<int> inverse_permutation(num_cols); |
| 706 | for (SuiteSparse_long i = 0; i < num_cols; ++i) { |
| 707 | inverse_permutation[permutation[i]] = i; |
| 708 | } |
| 709 | |
Sameer Agarwal | 5a97471 | 2013-06-07 22:38:30 -0700 | [diff] [blame] | 710 | const int* rows = covariance_matrix_->rows(); |
| 711 | const int* cols = covariance_matrix_->cols(); |
| 712 | double* values = covariance_matrix_->mutable_values(); |
| 713 | |
| 714 | // The following loop exploits the fact that the i^th column of A^{-1} |
| 715 | // is given by the solution to the linear system |
| 716 | // |
| 717 | // A x = e_i |
| 718 | // |
| 719 | // where e_i is a vector with e(i) = 1 and all other entries zero. |
| 720 | // |
| 721 | // Since the covariance matrix is symmetric, the i^th row and column |
| 722 | // are equal. |
Sameer Agarwal | b22d063 | 2013-08-15 22:55:23 -0700 | [diff] [blame] | 723 | const int num_threads = options_.num_threads; |
| 724 | scoped_array<double> workspace(new double[num_threads * num_cols]); |
Sameer Agarwal | 5a97471 | 2013-06-07 22:38:30 -0700 | [diff] [blame] | 725 | |
Sameer Agarwal | b22d063 | 2013-08-15 22:55:23 -0700 | [diff] [blame] | 726 | #pragma omp parallel for num_threads(num_threads) schedule(dynamic) |
Sameer Agarwal | 5a97471 | 2013-06-07 22:38:30 -0700 | [diff] [blame] | 727 | for (int r = 0; r < num_cols; ++r) { |
Sameer Agarwal | b22d063 | 2013-08-15 22:55:23 -0700 | [diff] [blame] | 728 | const int row_begin = rows[r]; |
| 729 | const int row_end = rows[r + 1]; |
Sameer Agarwal | 5a97471 | 2013-06-07 22:38:30 -0700 | [diff] [blame] | 730 | if (row_end == row_begin) { |
| 731 | continue; |
| 732 | } |
| 733 | |
Sameer Agarwal | b22d063 | 2013-08-15 22:55:23 -0700 | [diff] [blame] | 734 | # ifdef CERES_USE_OPENMP |
| 735 | int thread_id = omp_get_thread_num(); |
| 736 | # else |
| 737 | int thread_id = 0; |
| 738 | # endif |
Sameer Agarwal | 5a97471 | 2013-06-07 22:38:30 -0700 | [diff] [blame] | 739 | |
Sameer Agarwal | b22d063 | 2013-08-15 22:55:23 -0700 | [diff] [blame] | 740 | double* solution = workspace.get() + thread_id * num_cols; |
| 741 | SolveRTRWithSparseRHS<SuiteSparse_long>( |
| 742 | num_cols, |
| 743 | static_cast<SuiteSparse_long*>(R->i), |
| 744 | static_cast<SuiteSparse_long*>(R->p), |
| 745 | static_cast<double*>(R->x), |
| 746 | inverse_permutation[r], |
| 747 | solution); |
Sameer Agarwal | 5a97471 | 2013-06-07 22:38:30 -0700 | [diff] [blame] | 748 | for (int idx = row_begin; idx < row_end; ++idx) { |
Sameer Agarwal | b22d063 | 2013-08-15 22:55:23 -0700 | [diff] [blame] | 749 | const int c = cols[idx]; |
| 750 | values[idx] = solution[inverse_permutation[c]]; |
Sameer Agarwal | 5a97471 | 2013-06-07 22:38:30 -0700 | [diff] [blame] | 751 | } |
Sameer Agarwal | 5a97471 | 2013-06-07 22:38:30 -0700 | [diff] [blame] | 752 | } |
| 753 | |
Sameer Agarwal | 29310f0 | 2013-09-23 10:59:45 -0700 | [diff] [blame] | 754 | free(permutation); |
Sameer Agarwal | b22d063 | 2013-08-15 22:55:23 -0700 | [diff] [blame] | 755 | cholmod_l_free_sparse(&R, &cc); |
Sameer Agarwal | 5a97471 | 2013-06-07 22:38:30 -0700 | [diff] [blame] | 756 | cholmod_l_finish(&cc); |
| 757 | event_logger.AddEvent("Inversion"); |
| 758 | return true; |
| 759 | |
| 760 | #else // CERES_NO_SUITESPARSE |
| 761 | |
| 762 | return false; |
| 763 | |
| 764 | #endif // CERES_NO_SUITESPARSE |
Sameer Agarwal | c5bcfc0 | 2013-07-19 15:50:27 -0700 | [diff] [blame] | 765 | } |
Sameer Agarwal | 5a97471 | 2013-06-07 22:38:30 -0700 | [diff] [blame] | 766 | |
| 767 | bool CovarianceImpl::ComputeCovarianceValuesUsingDenseSVD() { |
| 768 | EventLogger event_logger( |
| 769 | "CovarianceImpl::ComputeCovarianceValuesUsingDenseSVD"); |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 770 | if (covariance_matrix_.get() == NULL) { |
| 771 | // Nothing to do, all zeros covariance matrix. |
| 772 | return true; |
| 773 | } |
| 774 | |
| 775 | CRSMatrix jacobian; |
| 776 | problem_->Evaluate(evaluate_options_, NULL, NULL, NULL, &jacobian); |
Sameer Agarwal | df01256 | 2013-05-21 10:59:54 -0700 | [diff] [blame] | 777 | event_logger.AddEvent("Evaluate"); |
| 778 | |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 779 | Matrix dense_jacobian(jacobian.num_rows, jacobian.num_cols); |
| 780 | dense_jacobian.setZero(); |
| 781 | for (int r = 0; r < jacobian.num_rows; ++r) { |
| 782 | for (int idx = jacobian.rows[r]; idx < jacobian.rows[r + 1]; ++idx) { |
| 783 | const int c = jacobian.cols[idx]; |
| 784 | dense_jacobian(r, c) = jacobian.values[idx]; |
| 785 | } |
| 786 | } |
Sameer Agarwal | df01256 | 2013-05-21 10:59:54 -0700 | [diff] [blame] | 787 | event_logger.AddEvent("ConvertToDenseMatrix"); |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 788 | |
| 789 | Eigen::JacobiSVD<Matrix> svd(dense_jacobian, |
| 790 | Eigen::ComputeThinU | Eigen::ComputeThinV); |
Sameer Agarwal | 5a97471 | 2013-06-07 22:38:30 -0700 | [diff] [blame] | 791 | |
Sameer Agarwal | 7129cd3 | 2013-06-02 19:35:29 -0700 | [diff] [blame] | 792 | event_logger.AddEvent("SingularValueDecomposition"); |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 793 | |
Sameer Agarwal | 7129cd3 | 2013-06-02 19:35:29 -0700 | [diff] [blame] | 794 | const Vector singular_values = svd.singularValues(); |
| 795 | const int num_singular_values = singular_values.rows(); |
| 796 | Vector inverse_squared_singular_values(num_singular_values); |
| 797 | inverse_squared_singular_values.setZero(); |
| 798 | |
| 799 | const double max_singular_value = singular_values[0]; |
| 800 | const double min_singular_value_ratio = |
| 801 | sqrt(options_.min_reciprocal_condition_number); |
| 802 | |
| 803 | const bool automatic_truncation = (options_.null_space_rank < 0); |
| 804 | const int max_rank = min(num_singular_values, |
| 805 | num_singular_values - options_.null_space_rank); |
| 806 | |
Sameer Agarwal | 8f7e896 | 2013-06-03 13:07:39 -0700 | [diff] [blame] | 807 | // Compute the squared inverse of the singular values. Truncate the |
| 808 | // computation based on min_singular_value_ratio and |
| 809 | // null_space_rank. When either of these two quantities are active, |
| 810 | // the resulting covariance matrix is a Moore-Penrose inverse |
| 811 | // instead of a regular inverse. |
Sameer Agarwal | 7129cd3 | 2013-06-02 19:35:29 -0700 | [diff] [blame] | 812 | for (int i = 0; i < max_rank; ++i) { |
| 813 | const double singular_value_ratio = singular_values[i] / max_singular_value; |
Sameer Agarwal | 8f7e896 | 2013-06-03 13:07:39 -0700 | [diff] [blame] | 814 | if (singular_value_ratio < min_singular_value_ratio) { |
| 815 | // Since the singular values are in decreasing order, if |
| 816 | // automatic truncation is enabled, then from this point on |
| 817 | // all values will fail the ratio test and there is nothing to |
| 818 | // do in this loop. |
| 819 | if (automatic_truncation) { |
| 820 | break; |
| 821 | } else { |
Sameer Agarwal | b16e118 | 2013-11-25 05:47:43 -0800 | [diff] [blame] | 822 | LOG(ERROR) << "Cholesky factorization of J'J is not reliable. " |
| 823 | << "Reciprocal condition number: " |
| 824 | << singular_value_ratio * singular_value_ratio << " " |
Sameer Agarwal | ed92366 | 2013-11-28 06:50:43 -0800 | [diff] [blame] | 825 | << "min_reciprocal_condition_number: " |
Sameer Agarwal | b16e118 | 2013-11-25 05:47:43 -0800 | [diff] [blame] | 826 | << options_.min_reciprocal_condition_number; |
Sameer Agarwal | 8f7e896 | 2013-06-03 13:07:39 -0700 | [diff] [blame] | 827 | return false; |
| 828 | } |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 829 | } |
Sameer Agarwal | 8f7e896 | 2013-06-03 13:07:39 -0700 | [diff] [blame] | 830 | |
| 831 | inverse_squared_singular_values[i] = |
| 832 | 1.0 / (singular_values[i] * singular_values[i]); |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 833 | } |
| 834 | |
| 835 | Matrix dense_covariance = |
| 836 | svd.matrixV() * |
Sameer Agarwal | 7129cd3 | 2013-06-02 19:35:29 -0700 | [diff] [blame] | 837 | inverse_squared_singular_values.asDiagonal() * |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 838 | svd.matrixV().transpose(); |
Sameer Agarwal | df01256 | 2013-05-21 10:59:54 -0700 | [diff] [blame] | 839 | event_logger.AddEvent("PseudoInverse"); |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 840 | |
| 841 | const int num_rows = covariance_matrix_->num_rows(); |
| 842 | const int* rows = covariance_matrix_->rows(); |
| 843 | const int* cols = covariance_matrix_->cols(); |
| 844 | double* values = covariance_matrix_->mutable_values(); |
| 845 | |
| 846 | for (int r = 0; r < num_rows; ++r) { |
| 847 | for (int idx = rows[r]; idx < rows[r + 1]; ++idx) { |
| 848 | const int c = cols[idx]; |
| 849 | values[idx] = dense_covariance(r, c); |
| 850 | } |
| 851 | } |
Sameer Agarwal | df01256 | 2013-05-21 10:59:54 -0700 | [diff] [blame] | 852 | event_logger.AddEvent("CopyToCovarianceMatrix"); |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 853 | return true; |
| 854 | }; |
| 855 | |
Sameer Agarwal | 02706c1 | 2013-05-12 22:07:55 -0700 | [diff] [blame] | 856 | } // namespace internal |
| 857 | } // namespace ceres |