Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 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@google.com (Sameer Agarwal) |
| 30 | |
| 31 | #include <algorithm> |
| 32 | #include <ctime> |
| 33 | #include <set> |
| 34 | #include <vector> |
Sameer Agarwal | b051873 | 2012-05-29 00:27:57 -0700 | [diff] [blame] | 35 | |
| 36 | #ifndef CERES_NO_CXSPARSE |
| 37 | #include "cs.h" |
| 38 | #endif // CERES_NO_CXSPARSE |
| 39 | |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 40 | #include "Eigen/Dense" |
| 41 | #include "ceres/block_random_access_dense_matrix.h" |
| 42 | #include "ceres/block_random_access_matrix.h" |
| 43 | #include "ceres/block_random_access_sparse_matrix.h" |
| 44 | #include "ceres/block_sparse_matrix.h" |
| 45 | #include "ceres/block_structure.h" |
| 46 | #include "ceres/detect_structure.h" |
| 47 | #include "ceres/linear_solver.h" |
| 48 | #include "ceres/schur_complement_solver.h" |
| 49 | #include "ceres/suitesparse.h" |
| 50 | #include "ceres/triplet_sparse_matrix.h" |
| 51 | #include "ceres/internal/eigen.h" |
| 52 | #include "ceres/internal/port.h" |
| 53 | #include "ceres/internal/scoped_ptr.h" |
| 54 | #include "ceres/types.h" |
| 55 | |
Sameer Agarwal | b051873 | 2012-05-29 00:27:57 -0700 | [diff] [blame] | 56 | |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 57 | namespace ceres { |
| 58 | namespace internal { |
| 59 | |
| 60 | LinearSolver::Summary SchurComplementSolver::SolveImpl( |
| 61 | BlockSparseMatrixBase* A, |
| 62 | const double* b, |
| 63 | const LinearSolver::PerSolveOptions& per_solve_options, |
| 64 | double* x) { |
| 65 | const time_t start_time = time(NULL); |
Sameer Agarwal | a9d8ef8 | 2012-05-14 02:28:05 -0700 | [diff] [blame] | 66 | if (eliminator_.get() == NULL) { |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 67 | InitStorage(A->block_structure()); |
| 68 | DetectStructure(*A->block_structure(), |
| 69 | options_.num_eliminate_blocks, |
| 70 | &options_.row_block_size, |
| 71 | &options_.e_block_size, |
| 72 | &options_.f_block_size); |
| 73 | eliminator_.reset(CHECK_NOTNULL(SchurEliminatorBase::Create(options_))); |
| 74 | eliminator_->Init(options_.num_eliminate_blocks, A->block_structure()); |
| 75 | }; |
| 76 | const time_t init_time = time(NULL); |
| 77 | fill(x, x + A->num_cols(), 0.0); |
| 78 | |
| 79 | LinearSolver::Summary summary; |
| 80 | summary.num_iterations = 1; |
| 81 | summary.termination_type = FAILURE; |
| 82 | eliminator_->Eliminate(A, b, per_solve_options.D, lhs_.get(), rhs_.get()); |
| 83 | const time_t eliminate_time = time(NULL); |
| 84 | |
| 85 | double* reduced_solution = x + A->num_cols() - lhs_->num_cols(); |
| 86 | const bool status = SolveReducedLinearSystem(reduced_solution); |
| 87 | const time_t solve_time = time(NULL); |
| 88 | |
| 89 | if (!status) { |
| 90 | return summary; |
| 91 | } |
| 92 | |
| 93 | eliminator_->BackSubstitute(A, b, per_solve_options.D, reduced_solution, x); |
| 94 | const time_t backsubstitute_time = time(NULL); |
| 95 | summary.termination_type = TOLERANCE; |
| 96 | |
Sameer Agarwal | 7a3c43b | 2012-06-05 23:10:59 -0700 | [diff] [blame] | 97 | VLOG(2) << "time (sec) total: " << (backsubstitute_time - start_time) |
| 98 | << " init: " << (init_time - start_time) |
| 99 | << " eliminate: " << (eliminate_time - init_time) |
| 100 | << " solve: " << (solve_time - eliminate_time) |
| 101 | << " backsubstitute: " << (backsubstitute_time - solve_time); |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 102 | return summary; |
| 103 | } |
| 104 | |
| 105 | // Initialize a BlockRandomAccessDenseMatrix to store the Schur |
| 106 | // complement. |
| 107 | void DenseSchurComplementSolver::InitStorage( |
| 108 | const CompressedRowBlockStructure* bs) { |
| 109 | const int num_eliminate_blocks = options().num_eliminate_blocks; |
| 110 | const int num_col_blocks = bs->cols.size(); |
| 111 | |
| 112 | vector<int> blocks(num_col_blocks - num_eliminate_blocks, 0); |
| 113 | for (int i = num_eliminate_blocks, j = 0; |
| 114 | i < num_col_blocks; |
| 115 | ++i, ++j) { |
| 116 | blocks[j] = bs->cols[i].size; |
| 117 | } |
| 118 | |
| 119 | set_lhs(new BlockRandomAccessDenseMatrix(blocks)); |
| 120 | set_rhs(new double[lhs()->num_rows()]); |
| 121 | } |
| 122 | |
| 123 | // Solve the system Sx = r, assuming that the matrix S is stored in a |
| 124 | // BlockRandomAccessDenseMatrix. The linear system is solved using |
| 125 | // Eigen's Cholesky factorization. |
| 126 | bool DenseSchurComplementSolver::SolveReducedLinearSystem(double* solution) { |
| 127 | const BlockRandomAccessDenseMatrix* m = |
| 128 | down_cast<const BlockRandomAccessDenseMatrix*>(lhs()); |
| 129 | const int num_rows = m->num_rows(); |
| 130 | |
| 131 | // The case where there are no f blocks, and the system is block |
| 132 | // diagonal. |
| 133 | if (num_rows == 0) { |
| 134 | return true; |
| 135 | } |
| 136 | |
| 137 | // TODO(sameeragarwal): Add proper error handling; this completely ignores |
| 138 | // the quality of the solution to the solve. |
| 139 | VectorRef(solution, num_rows) = |
| 140 | ConstMatrixRef(m->values(), num_rows, num_rows) |
| 141 | .selfadjointView<Eigen::Upper>() |
| 142 | .ldlt() |
| 143 | .solve(ConstVectorRef(rhs(), num_rows)); |
| 144 | |
| 145 | return true; |
| 146 | } |
| 147 | |
Sameer Agarwal | b051873 | 2012-05-29 00:27:57 -0700 | [diff] [blame] | 148 | |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 149 | SparseSchurComplementSolver::SparseSchurComplementSolver( |
| 150 | const LinearSolver::Options& options) |
Sameer Agarwal | b051873 | 2012-05-29 00:27:57 -0700 | [diff] [blame] | 151 | : SchurComplementSolver(options) { |
| 152 | #ifndef CERES_NO_SUITESPARSE |
Sameer Agarwal | 7a3c43b | 2012-06-05 23:10:59 -0700 | [diff] [blame] | 153 | factor_ = NULL; |
Sameer Agarwal | b051873 | 2012-05-29 00:27:57 -0700 | [diff] [blame] | 154 | #endif // CERES_NO_SUITESPARSE |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 155 | } |
| 156 | |
| 157 | SparseSchurComplementSolver::~SparseSchurComplementSolver() { |
Sameer Agarwal | b051873 | 2012-05-29 00:27:57 -0700 | [diff] [blame] | 158 | #ifndef CERES_NO_SUITESPARSE |
Sameer Agarwal | 7a3c43b | 2012-06-05 23:10:59 -0700 | [diff] [blame] | 159 | if (factor_ != NULL) { |
| 160 | ss_.Free(factor_); |
| 161 | factor_ = NULL; |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 162 | } |
Sameer Agarwal | b051873 | 2012-05-29 00:27:57 -0700 | [diff] [blame] | 163 | #endif // CERES_NO_SUITESPARSE |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 164 | } |
| 165 | |
| 166 | // Determine the non-zero blocks in the Schur Complement matrix, and |
| 167 | // initialize a BlockRandomAccessSparseMatrix object. |
| 168 | void SparseSchurComplementSolver::InitStorage( |
| 169 | const CompressedRowBlockStructure* bs) { |
| 170 | const int num_eliminate_blocks = options().num_eliminate_blocks; |
| 171 | const int num_col_blocks = bs->cols.size(); |
| 172 | const int num_row_blocks = bs->rows.size(); |
| 173 | |
Sameer Agarwal | 7a3c43b | 2012-06-05 23:10:59 -0700 | [diff] [blame] | 174 | blocks_.resize(num_col_blocks - num_eliminate_blocks, 0); |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 175 | for (int i = num_eliminate_blocks; i < num_col_blocks; ++i) { |
Sameer Agarwal | 7a3c43b | 2012-06-05 23:10:59 -0700 | [diff] [blame] | 176 | blocks_[i - num_eliminate_blocks] = bs->cols[i].size; |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 177 | } |
| 178 | |
| 179 | set<pair<int, int> > block_pairs; |
Sameer Agarwal | 7a3c43b | 2012-06-05 23:10:59 -0700 | [diff] [blame] | 180 | for (int i = 0; i < blocks_.size(); ++i) { |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 181 | block_pairs.insert(make_pair(i, i)); |
| 182 | } |
| 183 | |
| 184 | int r = 0; |
| 185 | while (r < num_row_blocks) { |
| 186 | int e_block_id = bs->rows[r].cells.front().block_id; |
| 187 | if (e_block_id >= num_eliminate_blocks) { |
| 188 | break; |
| 189 | } |
| 190 | vector<int> f_blocks; |
| 191 | |
| 192 | // Add to the chunk until the first block in the row is |
| 193 | // different than the one in the first row for the chunk. |
| 194 | for (; r < num_row_blocks; ++r) { |
| 195 | const CompressedRow& row = bs->rows[r]; |
| 196 | if (row.cells.front().block_id != e_block_id) { |
| 197 | break; |
| 198 | } |
| 199 | |
| 200 | // Iterate over the blocks in the row, ignoring the first |
| 201 | // block since it is the one to be eliminated. |
| 202 | for (int c = 1; c < row.cells.size(); ++c) { |
| 203 | const Cell& cell = row.cells[c]; |
| 204 | f_blocks.push_back(cell.block_id - num_eliminate_blocks); |
| 205 | } |
| 206 | } |
| 207 | |
| 208 | sort(f_blocks.begin(), f_blocks.end()); |
| 209 | f_blocks.erase(unique(f_blocks.begin(), f_blocks.end()), f_blocks.end()); |
| 210 | for (int i = 0; i < f_blocks.size(); ++i) { |
| 211 | for (int j = i + 1; j < f_blocks.size(); ++j) { |
| 212 | block_pairs.insert(make_pair(f_blocks[i], f_blocks[j])); |
| 213 | } |
| 214 | } |
| 215 | } |
| 216 | |
| 217 | // Remaing rows do not contribute to the chunks and directly go |
| 218 | // into the schur complement via an outer product. |
| 219 | for (; r < num_row_blocks; ++r) { |
| 220 | const CompressedRow& row = bs->rows[r]; |
| 221 | CHECK_GE(row.cells.front().block_id, num_eliminate_blocks); |
| 222 | for (int i = 0; i < row.cells.size(); ++i) { |
| 223 | int r_block1_id = row.cells[i].block_id - num_eliminate_blocks; |
| 224 | for (int j = 0; j < row.cells.size(); ++j) { |
| 225 | int r_block2_id = row.cells[j].block_id - num_eliminate_blocks; |
| 226 | if (r_block1_id <= r_block2_id) { |
| 227 | block_pairs.insert(make_pair(r_block1_id, r_block2_id)); |
| 228 | } |
| 229 | } |
| 230 | } |
| 231 | } |
| 232 | |
Sameer Agarwal | 7a3c43b | 2012-06-05 23:10:59 -0700 | [diff] [blame] | 233 | set_lhs(new BlockRandomAccessSparseMatrix(blocks_, block_pairs)); |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 234 | set_rhs(new double[lhs()->num_rows()]); |
| 235 | } |
| 236 | |
Sameer Agarwal | b051873 | 2012-05-29 00:27:57 -0700 | [diff] [blame] | 237 | bool SparseSchurComplementSolver::SolveReducedLinearSystem(double* solution) { |
| 238 | switch (options().sparse_linear_algebra_library) { |
| 239 | case SUITE_SPARSE: |
| 240 | return SolveReducedLinearSystemUsingSuiteSparse(solution); |
| 241 | case CX_SPARSE: |
| 242 | return SolveReducedLinearSystemUsingCXSparse(solution); |
| 243 | default: |
| 244 | LOG(FATAL) << "Unknown sparse linear algebra library : " |
| 245 | << options().sparse_linear_algebra_library; |
| 246 | } |
| 247 | |
| 248 | LOG(FATAL) << "Unknown sparse linear algebra library : " |
| 249 | << options().sparse_linear_algebra_library; |
| 250 | return false; |
| 251 | } |
| 252 | |
| 253 | #ifndef CERES_NO_SUITESPARSE |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 254 | // Solve the system Sx = r, assuming that the matrix S is stored in a |
| 255 | // BlockRandomAccessSparseMatrix. The linear system is solved using |
| 256 | // CHOLMOD's sparse cholesky factorization routines. |
Sameer Agarwal | b051873 | 2012-05-29 00:27:57 -0700 | [diff] [blame] | 257 | bool SparseSchurComplementSolver::SolveReducedLinearSystemUsingSuiteSparse( |
| 258 | double* solution) { |
Sameer Agarwal | 7a3c43b | 2012-06-05 23:10:59 -0700 | [diff] [blame] | 259 | const time_t start_time = time(NULL); |
| 260 | |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 261 | TripletSparseMatrix* tsm = |
| 262 | const_cast<TripletSparseMatrix*>( |
| 263 | down_cast<const BlockRandomAccessSparseMatrix*>(lhs())->matrix()); |
| 264 | |
| 265 | const int num_rows = tsm->num_rows(); |
| 266 | |
| 267 | // The case where there are no f blocks, and the system is block |
| 268 | // diagonal. |
| 269 | if (num_rows == 0) { |
| 270 | return true; |
| 271 | } |
| 272 | |
| 273 | cholmod_sparse* cholmod_lhs = ss_.CreateSparseMatrix(tsm); |
| 274 | // The matrix is symmetric, and the upper triangular part of the |
| 275 | // matrix contains the values. |
| 276 | cholmod_lhs->stype = 1; |
Sameer Agarwal | 7a3c43b | 2012-06-05 23:10:59 -0700 | [diff] [blame] | 277 | const time_t lhs_time = time(NULL); |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 278 | |
| 279 | cholmod_dense* cholmod_rhs = |
| 280 | ss_.CreateDenseVector(const_cast<double*>(rhs()), num_rows, num_rows); |
Sameer Agarwal | 7a3c43b | 2012-06-05 23:10:59 -0700 | [diff] [blame] | 281 | const time_t rhs_time = time(NULL); |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 282 | |
| 283 | // Symbolic factorization is computed if we don't already have one handy. |
Sameer Agarwal | 7a3c43b | 2012-06-05 23:10:59 -0700 | [diff] [blame] | 284 | if (factor_ == NULL) { |
| 285 | if (options().use_block_amd) { |
| 286 | factor_ = ss_.BlockAnalyzeCholesky(cholmod_lhs, blocks_, blocks_); |
| 287 | } else { |
| 288 | factor_ = ss_.AnalyzeCholesky(cholmod_lhs); |
| 289 | } |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 290 | |
Sameer Agarwal | cb83b28 | 2012-06-06 22:26:09 -0700 | [diff] [blame] | 291 | if (VLOG_IS_ON(2)) { |
| 292 | cholmod_print_common("Symbolic Analysis", ss_.mutable_cc()); |
| 293 | } |
Sameer Agarwal | 7a3c43b | 2012-06-05 23:10:59 -0700 | [diff] [blame] | 294 | } |
| 295 | |
| 296 | CHECK_NOTNULL(factor_); |
| 297 | |
| 298 | const time_t symbolic_time = time(NULL); |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 299 | cholmod_dense* cholmod_solution = |
Sameer Agarwal | 7a3c43b | 2012-06-05 23:10:59 -0700 | [diff] [blame] | 300 | ss_.SolveCholesky(cholmod_lhs, factor_, cholmod_rhs); |
| 301 | |
| 302 | const time_t solve_time = time(NULL); |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 303 | |
| 304 | ss_.Free(cholmod_lhs); |
| 305 | cholmod_lhs = NULL; |
| 306 | ss_.Free(cholmod_rhs); |
| 307 | cholmod_rhs = NULL; |
| 308 | |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 309 | if (cholmod_solution == NULL) { |
| 310 | LOG(ERROR) << "CHOLMOD solve failed."; |
| 311 | return false; |
| 312 | } |
| 313 | |
| 314 | VectorRef(solution, num_rows) |
| 315 | = VectorRef(static_cast<double*>(cholmod_solution->x), num_rows); |
| 316 | ss_.Free(cholmod_solution); |
Sameer Agarwal | 7a3c43b | 2012-06-05 23:10:59 -0700 | [diff] [blame] | 317 | const time_t final_time = time(NULL); |
| 318 | VLOG(2) << "time: " << (final_time - start_time) |
| 319 | << " lhs : " << (lhs_time - start_time) |
| 320 | << " rhs: " << (rhs_time - lhs_time) |
| 321 | << " analyze: " << (symbolic_time - rhs_time) |
| 322 | << " factor_and_solve: " << (solve_time - symbolic_time) |
| 323 | << " cleanup: " << (final_time - solve_time); |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 324 | return true; |
| 325 | } |
Sameer Agarwal | b051873 | 2012-05-29 00:27:57 -0700 | [diff] [blame] | 326 | #else |
| 327 | bool SparseSchurComplementSolver::SolveReducedLinearSystemUsingSuiteSparse( |
| 328 | double* solution) { |
| 329 | LOG(FATAL) << "No SuiteSparse support in Ceres."; |
| 330 | return false; |
| 331 | } |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 332 | #endif // CERES_NO_SUITESPARSE |
| 333 | |
Sameer Agarwal | b051873 | 2012-05-29 00:27:57 -0700 | [diff] [blame] | 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_sparse tsm_wrapper; |
| 354 | tsm_wrapper.nzmax = tsm->num_nonzeros(); |
| 355 | tsm_wrapper.m = num_rows; |
| 356 | tsm_wrapper.n = num_rows; |
| 357 | tsm_wrapper.p = tsm->mutable_cols(); |
| 358 | tsm_wrapper.i = tsm->mutable_rows(); |
| 359 | tsm_wrapper.x = tsm->mutable_values(); |
| 360 | tsm_wrapper.nz = tsm->num_nonzeros(); |
| 361 | |
| 362 | cs_di_sparse* lhs = cs_compress(&tsm_wrapper); |
| 363 | VectorRef(solution, num_rows) = ConstVectorRef(rhs(), num_rows); |
| 364 | |
| 365 | // It maybe worth caching the ordering here, but for now we are |
| 366 | // going to go with the simple cholsol based implementation. |
| 367 | int ok = cs_di_cholsol(1, lhs, solution); |
| 368 | cs_free(lhs); |
| 369 | return ok; |
| 370 | } |
| 371 | #else |
| 372 | bool SparseSchurComplementSolver::SolveReducedLinearSystemUsingCXSparse( |
| 373 | double* solution) { |
| 374 | LOG(FATAL) << "No CXSparse support in Ceres."; |
| 375 | return false; |
| 376 | } |
| 377 | #endif // CERES_NO_CXPARSE |
| 378 | |
Keir Mierle | 8ebb073 | 2012-04-30 23:09:08 -0700 | [diff] [blame] | 379 | } // namespace internal |
| 380 | } // namespace ceres |