1 // Ceres Solver - A fast non-linear least squares minimizer 2 // Copyright 2010, 2011, 2012 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 (at) google.com (Sameer Agarwal) 30 31 #include <algorithm> 32 #include <ctime> 33 #include <set> 34 #include <vector> 35 36 #include "Eigen/Dense" 37 #include "ceres/block_random_access_dense_matrix.h" 38 #include "ceres/block_random_access_matrix.h" 39 #include "ceres/block_random_access_sparse_matrix.h" 40 #include "ceres/block_sparse_matrix.h" 41 #include "ceres/block_structure.h" 42 #include "ceres/cxsparse.h" 43 #include "ceres/detect_structure.h" 44 #include "ceres/internal/eigen.h" 45 #include "ceres/internal/port.h" 46 #include "ceres/internal/scoped_ptr.h" 47 #include "ceres/lapack.h" 48 #include "ceres/linear_solver.h" 49 #include "ceres/schur_complement_solver.h" 50 #include "ceres/suitesparse.h" 51 #include "ceres/triplet_sparse_matrix.h" 52 #include "ceres/types.h" 53 #include "ceres/wall_time.h" 54 55 namespace ceres { 56 namespace internal { 57 58 LinearSolver::Summary SchurComplementSolver::SolveImpl( 59 BlockSparseMatrix* A, 60 const double* b, 61 const LinearSolver::PerSolveOptions& per_solve_options, 62 double* x) { 63 EventLogger event_logger("SchurComplementSolver::Solve"); 64 65 if (eliminator_.get() == NULL) { 66 InitStorage(A->block_structure()); 67 DetectStructure(*A->block_structure(), 68 options_.elimination_groups[0], 69 &options_.row_block_size, 70 &options_.e_block_size, 71 &options_.f_block_size); 72 eliminator_.reset(CHECK_NOTNULL(SchurEliminatorBase::Create(options_))); 73 eliminator_->Init(options_.elimination_groups[0], A->block_structure()); 74 }; 75 fill(x, x + A->num_cols(), 0.0); 76 event_logger.AddEvent("Setup"); 77 78 LinearSolver::Summary summary; 79 summary.num_iterations = 1; 80 summary.termination_type = FAILURE; 81 eliminator_->Eliminate(A, b, per_solve_options.D, lhs_.get(), rhs_.get()); 82 event_logger.AddEvent("Eliminate"); 83 84 double* reduced_solution = x + A->num_cols() - lhs_->num_cols(); 85 const bool status = SolveReducedLinearSystem(reduced_solution); 86 event_logger.AddEvent("ReducedSolve"); 87 88 if (!status) { 89 return summary; 90 } 91 92 eliminator_->BackSubstitute(A, b, per_solve_options.D, reduced_solution, x); 93 summary.termination_type = TOLERANCE; 94 95 event_logger.AddEvent("BackSubstitute"); 96 return summary; 97 } 98 99 // Initialize a BlockRandomAccessDenseMatrix to store the Schur 100 // complement. 101 void DenseSchurComplementSolver::InitStorage( 102 const CompressedRowBlockStructure* bs) { 103 const int num_eliminate_blocks = options().elimination_groups[0]; 104 const int num_col_blocks = bs->cols.size(); 105 106 vector<int> blocks(num_col_blocks - num_eliminate_blocks, 0); 107 for (int i = num_eliminate_blocks, j = 0; 108 i < num_col_blocks; 109 ++i, ++j) { 110 blocks[j] = bs->cols[i].size; 111 } 112 113 set_lhs(new BlockRandomAccessDenseMatrix(blocks)); 114 set_rhs(new double[lhs()->num_rows()]); 115 } 116 117 // Solve the system Sx = r, assuming that the matrix S is stored in a 118 // BlockRandomAccessDenseMatrix. The linear system is solved using 119 // Eigen's Cholesky factorization. 120 bool DenseSchurComplementSolver::SolveReducedLinearSystem(double* solution) { 121 const BlockRandomAccessDenseMatrix* m = 122 down_cast<const BlockRandomAccessDenseMatrix*>(lhs()); 123 const int num_rows = m->num_rows(); 124 125 // The case where there are no f blocks, and the system is block 126 // diagonal. 127 if (num_rows == 0) { 128 return true; 129 } 130 131 if (options().dense_linear_algebra_library_type == EIGEN) { 132 // TODO(sameeragarwal): Add proper error handling; this completely ignores 133 // the quality of the solution to the solve. 134 VectorRef(solution, num_rows) = 135 ConstMatrixRef(m->values(), num_rows, num_rows) 136 .selfadjointView<Eigen::Upper>() 137 .llt() 138 .solve(ConstVectorRef(rhs(), num_rows)); 139 return true; 140 } 141 142 VectorRef(solution, num_rows) = ConstVectorRef(rhs(), num_rows); 143 const int info = LAPACK::SolveInPlaceUsingCholesky(num_rows, 144 m->values(), 145 solution); 146 return (info == 0); 147 } 148 149 #if !defined(CERES_NO_SUITESPARSE) || !defined(CERES_NO_CXSPARE) 150 151 SparseSchurComplementSolver::SparseSchurComplementSolver( 152 const LinearSolver::Options& options) 153 : SchurComplementSolver(options), 154 factor_(NULL), 155 cxsparse_factor_(NULL) { 156 } 157 158 SparseSchurComplementSolver::~SparseSchurComplementSolver() { 159 #ifndef CERES_NO_SUITESPARSE 160 if (factor_ != NULL) { 161 ss_.Free(factor_); 162 factor_ = NULL; 163 } 164 #endif // CERES_NO_SUITESPARSE 165 166 #ifndef CERES_NO_CXSPARSE 167 if (cxsparse_factor_ != NULL) { 168 cxsparse_.Free(cxsparse_factor_); 169 cxsparse_factor_ = NULL; 170 } 171 #endif // CERES_NO_CXSPARSE 172 } 173 174 // Determine the non-zero blocks in the Schur Complement matrix, and 175 // initialize a BlockRandomAccessSparseMatrix object. 176 void SparseSchurComplementSolver::InitStorage( 177 const CompressedRowBlockStructure* bs) { 178 const int num_eliminate_blocks = options().elimination_groups[0]; 179 const int num_col_blocks = bs->cols.size(); 180 const int num_row_blocks = bs->rows.size(); 181 182 blocks_.resize(num_col_blocks - num_eliminate_blocks, 0); 183 for (int i = num_eliminate_blocks; i < num_col_blocks; ++i) { 184 blocks_[i - num_eliminate_blocks] = bs->cols[i].size; 185 } 186 187 set<pair<int, int> > block_pairs; 188 for (int i = 0; i < blocks_.size(); ++i) { 189 block_pairs.insert(make_pair(i, i)); 190 } 191 192 int r = 0; 193 while (r < num_row_blocks) { 194 int e_block_id = bs->rows[r].cells.front().block_id; 195 if (e_block_id >= num_eliminate_blocks) { 196 break; 197 } 198 vector<int> f_blocks; 199 200 // Add to the chunk until the first block in the row is 201 // different than the one in the first row for the chunk. 202 for (; r < num_row_blocks; ++r) { 203 const CompressedRow& row = bs->rows[r]; 204 if (row.cells.front().block_id != e_block_id) { 205 break; 206 } 207 208 // Iterate over the blocks in the row, ignoring the first 209 // block since it is the one to be eliminated. 210 for (int c = 1; c < row.cells.size(); ++c) { 211 const Cell& cell = row.cells[c]; 212 f_blocks.push_back(cell.block_id - num_eliminate_blocks); 213 } 214 } 215 216 sort(f_blocks.begin(), f_blocks.end()); 217 f_blocks.erase(unique(f_blocks.begin(), f_blocks.end()), f_blocks.end()); 218 for (int i = 0; i < f_blocks.size(); ++i) { 219 for (int j = i + 1; j < f_blocks.size(); ++j) { 220 block_pairs.insert(make_pair(f_blocks[i], f_blocks[j])); 221 } 222 } 223 } 224 225 // Remaing rows do not contribute to the chunks and directly go 226 // into the schur complement via an outer product. 227 for (; r < num_row_blocks; ++r) { 228 const CompressedRow& row = bs->rows[r]; 229 CHECK_GE(row.cells.front().block_id, num_eliminate_blocks); 230 for (int i = 0; i < row.cells.size(); ++i) { 231 int r_block1_id = row.cells[i].block_id - num_eliminate_blocks; 232 for (int j = 0; j < row.cells.size(); ++j) { 233 int r_block2_id = row.cells[j].block_id - num_eliminate_blocks; 234 if (r_block1_id <= r_block2_id) { 235 block_pairs.insert(make_pair(r_block1_id, r_block2_id)); 236 } 237 } 238 } 239 } 240 241 set_lhs(new BlockRandomAccessSparseMatrix(blocks_, block_pairs)); 242 set_rhs(new double[lhs()->num_rows()]); 243 } 244 245 bool SparseSchurComplementSolver::SolveReducedLinearSystem(double* solution) { 246 switch (options().sparse_linear_algebra_library_type) { 247 case SUITE_SPARSE: 248 return SolveReducedLinearSystemUsingSuiteSparse(solution); 249 case CX_SPARSE: 250 return SolveReducedLinearSystemUsingCXSparse(solution); 251 default: 252 LOG(FATAL) << "Unknown sparse linear algebra library : " 253 << options().sparse_linear_algebra_library_type; 254 } 255 256 LOG(FATAL) << "Unknown sparse linear algebra library : " 257 << options().sparse_linear_algebra_library_type; 258 return false; 259 } 260 261 #ifndef CERES_NO_SUITESPARSE 262 // Solve the system Sx = r, assuming that the matrix S is stored in a 263 // BlockRandomAccessSparseMatrix. The linear system is solved using 264 // CHOLMOD's sparse cholesky factorization routines. 265 bool SparseSchurComplementSolver::SolveReducedLinearSystemUsingSuiteSparse( 266 double* solution) { 267 TripletSparseMatrix* tsm = 268 const_cast<TripletSparseMatrix*>( 269 down_cast<const BlockRandomAccessSparseMatrix*>(lhs())->matrix()); 270 271 const int num_rows = tsm->num_rows(); 272 273 // The case where there are no f blocks, and the system is block 274 // diagonal. 275 if (num_rows == 0) { 276 return true; 277 } 278 279 cholmod_sparse* cholmod_lhs = NULL; 280 if (options().use_postordering) { 281 // If we are going to do a full symbolic analysis of the schur 282 // complement matrix from scratch and not rely on the 283 // pre-ordering, then the fastest path in cholmod_factorize is the 284 // one corresponding to upper triangular matrices. 285 286 // Create a upper triangular symmetric matrix. 287 cholmod_lhs = ss_.CreateSparseMatrix(tsm); 288 cholmod_lhs->stype = 1; 289 290 if (factor_ == NULL) { 291 factor_ = ss_.BlockAnalyzeCholesky(cholmod_lhs, blocks_, blocks_); 292 } 293 } else { 294 // If we are going to use the natural ordering (i.e. rely on the 295 // pre-ordering computed by solver_impl.cc), then the fastest 296 // path in cholmod_factorize is the one corresponding to lower 297 // triangular matrices. 298 299 // Create a upper triangular symmetric matrix. 300 cholmod_lhs = ss_.CreateSparseMatrixTranspose(tsm); 301 cholmod_lhs->stype = -1; 302 303 if (factor_ == NULL) { 304 factor_ = ss_.AnalyzeCholeskyWithNaturalOrdering(cholmod_lhs); 305 } 306 } 307 308 cholmod_dense* cholmod_rhs = 309 ss_.CreateDenseVector(const_cast<double*>(rhs()), num_rows, num_rows); 310 cholmod_dense* cholmod_solution = 311 ss_.SolveCholesky(cholmod_lhs, factor_, cholmod_rhs); 312 313 ss_.Free(cholmod_lhs); 314 ss_.Free(cholmod_rhs); 315 316 if (cholmod_solution == NULL) { 317 LOG(WARNING) << "CHOLMOD solve failed."; 318 return false; 319 } 320 321 VectorRef(solution, num_rows) 322 = VectorRef(static_cast<double*>(cholmod_solution->x), num_rows); 323 ss_.Free(cholmod_solution); 324 return true; 325 } 326 #else 327 bool SparseSchurComplementSolver::SolveReducedLinearSystemUsingSuiteSparse( 328 double* solution) { 329 LOG(FATAL) << "No SuiteSparse support in Ceres."; 330 return false; 331 } 332 #endif // CERES_NO_SUITESPARSE 333 334 #ifndef CERES_NO_CXSPARSE 335 // Solve the system Sx = r, assuming that the matrix S is stored in a 336 // BlockRandomAccessSparseMatrix. The linear system is solved using 337 // CXSparse's sparse cholesky factorization routines. 338 bool SparseSchurComplementSolver::SolveReducedLinearSystemUsingCXSparse( 339 double* solution) { 340 // Extract the TripletSparseMatrix that is used for actually storing S. 341 TripletSparseMatrix* tsm = 342 const_cast<TripletSparseMatrix*>( 343 down_cast<const BlockRandomAccessSparseMatrix*>(lhs())->matrix()); 344 345 const int num_rows = tsm->num_rows(); 346 347 // The case where there are no f blocks, and the system is block 348 // diagonal. 349 if (num_rows == 0) { 350 return true; 351 } 352 353 cs_di* lhs = CHECK_NOTNULL(cxsparse_.CreateSparseMatrix(tsm)); 354 VectorRef(solution, num_rows) = ConstVectorRef(rhs(), num_rows); 355 356 // Compute symbolic factorization if not available. 357 if (cxsparse_factor_ == NULL) { 358 cxsparse_factor_ = 359 CHECK_NOTNULL(cxsparse_.BlockAnalyzeCholesky(lhs, blocks_, blocks_)); 360 } 361 362 // Solve the linear system. 363 bool ok = cxsparse_.SolveCholesky(lhs, cxsparse_factor_, solution); 364 365 cxsparse_.Free(lhs); 366 return ok; 367 } 368 #else 369 bool SparseSchurComplementSolver::SolveReducedLinearSystemUsingCXSparse( 370 double* solution) { 371 LOG(FATAL) << "No CXSparse support in Ceres."; 372 return false; 373 } 374 #endif // CERES_NO_CXPARSE 375 376 #endif // !defined(CERES_NO_SUITESPARSE) || !defined(CERES_NO_CXSPARE) 377 } // namespace internal 378 } // namespace ceres 379