blob: 17537596c75be6dc3a0382f29efe955bbaf0e330 [file] [log] [blame]
Keir Mierle8ebb0732012-04-30 23:09:08 -07001// 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 Agarwalb0518732012-05-29 00:27:57 -070035
36#ifndef CERES_NO_CXSPARSE
37#include "cs.h"
38#endif // CERES_NO_CXSPARSE
39
Keir Mierle8ebb0732012-04-30 23:09:08 -070040#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"
Sameer Agarwal42a84b82013-02-01 12:22:53 -080047#include "ceres/internal/eigen.h"
48#include "ceres/internal/port.h"
49#include "ceres/internal/scoped_ptr.h"
Keir Mierle8ebb0732012-04-30 23:09:08 -070050#include "ceres/linear_solver.h"
51#include "ceres/schur_complement_solver.h"
52#include "ceres/suitesparse.h"
53#include "ceres/triplet_sparse_matrix.h"
Keir Mierle8ebb0732012-04-30 23:09:08 -070054#include "ceres/types.h"
Sameer Agarwal42a84b82013-02-01 12:22:53 -080055#include "ceres/wall_time.h"
Sameer Agarwalb0518732012-05-29 00:27:57 -070056
Keir Mierle8ebb0732012-04-30 23:09:08 -070057namespace ceres {
58namespace internal {
59
60LinearSolver::Summary SchurComplementSolver::SolveImpl(
61 BlockSparseMatrixBase* A,
62 const double* b,
63 const LinearSolver::PerSolveOptions& per_solve_options,
64 double* x) {
Sameer Agarwal42a84b82013-02-01 12:22:53 -080065 EventLogger event_logger("SchurComplementSolver::Solve");
66
Sameer Agarwala9d8ef82012-05-14 02:28:05 -070067 if (eliminator_.get() == NULL) {
Keir Mierle8ebb0732012-04-30 23:09:08 -070068 InitStorage(A->block_structure());
69 DetectStructure(*A->block_structure(),
Sameer Agarwal0c52f1e2012-09-17 11:30:14 -070070 options_.elimination_groups[0],
Keir Mierle8ebb0732012-04-30 23:09:08 -070071 &options_.row_block_size,
72 &options_.e_block_size,
73 &options_.f_block_size);
74 eliminator_.reset(CHECK_NOTNULL(SchurEliminatorBase::Create(options_)));
Sameer Agarwal0c52f1e2012-09-17 11:30:14 -070075 eliminator_->Init(options_.elimination_groups[0], A->block_structure());
Keir Mierle8ebb0732012-04-30 23:09:08 -070076 };
Keir Mierle8ebb0732012-04-30 23:09:08 -070077 fill(x, x + A->num_cols(), 0.0);
Sameer Agarwal42a84b82013-02-01 12:22:53 -080078 event_logger.AddEvent("Setup");
Keir Mierle8ebb0732012-04-30 23:09:08 -070079
80 LinearSolver::Summary summary;
81 summary.num_iterations = 1;
82 summary.termination_type = FAILURE;
83 eliminator_->Eliminate(A, b, per_solve_options.D, lhs_.get(), rhs_.get());
Sameer Agarwal42a84b82013-02-01 12:22:53 -080084 event_logger.AddEvent("Eliminate");
Keir Mierle8ebb0732012-04-30 23:09:08 -070085
86 double* reduced_solution = x + A->num_cols() - lhs_->num_cols();
87 const bool status = SolveReducedLinearSystem(reduced_solution);
Sameer Agarwal42a84b82013-02-01 12:22:53 -080088 event_logger.AddEvent("ReducedSolve");
Keir Mierle8ebb0732012-04-30 23:09:08 -070089
90 if (!status) {
91 return summary;
92 }
93
94 eliminator_->BackSubstitute(A, b, per_solve_options.D, reduced_solution, x);
Keir Mierle8ebb0732012-04-30 23:09:08 -070095 summary.termination_type = TOLERANCE;
96
Sameer Agarwal42a84b82013-02-01 12:22:53 -080097 event_logger.AddEvent("BackSubstitute");
Keir Mierle8ebb0732012-04-30 23:09:08 -070098 return summary;
99}
100
101// Initialize a BlockRandomAccessDenseMatrix to store the Schur
102// complement.
103void DenseSchurComplementSolver::InitStorage(
104 const CompressedRowBlockStructure* bs) {
Sameer Agarwal0c52f1e2012-09-17 11:30:14 -0700105 const int num_eliminate_blocks = options().elimination_groups[0];
Keir Mierle8ebb0732012-04-30 23:09:08 -0700106 const int num_col_blocks = bs->cols.size();
107
108 vector<int> blocks(num_col_blocks - num_eliminate_blocks, 0);
109 for (int i = num_eliminate_blocks, j = 0;
110 i < num_col_blocks;
111 ++i, ++j) {
112 blocks[j] = bs->cols[i].size;
113 }
114
115 set_lhs(new BlockRandomAccessDenseMatrix(blocks));
116 set_rhs(new double[lhs()->num_rows()]);
117}
118
119// Solve the system Sx = r, assuming that the matrix S is stored in a
120// BlockRandomAccessDenseMatrix. The linear system is solved using
121// Eigen's Cholesky factorization.
122bool DenseSchurComplementSolver::SolveReducedLinearSystem(double* solution) {
123 const BlockRandomAccessDenseMatrix* m =
124 down_cast<const BlockRandomAccessDenseMatrix*>(lhs());
125 const int num_rows = m->num_rows();
126
127 // The case where there are no f blocks, and the system is block
128 // diagonal.
129 if (num_rows == 0) {
130 return true;
131 }
132
133 // TODO(sameeragarwal): Add proper error handling; this completely ignores
134 // the quality of the solution to the solve.
135 VectorRef(solution, num_rows) =
136 ConstMatrixRef(m->values(), num_rows, num_rows)
137 .selfadjointView<Eigen::Upper>()
138 .ldlt()
139 .solve(ConstVectorRef(rhs(), num_rows));
140
141 return true;
142}
143
Sameer Agarwalb0518732012-05-29 00:27:57 -0700144
Keir Mierle8ebb0732012-04-30 23:09:08 -0700145SparseSchurComplementSolver::SparseSchurComplementSolver(
146 const LinearSolver::Options& options)
Sameer Agarwalb0518732012-05-29 00:27:57 -0700147 : SchurComplementSolver(options) {
148#ifndef CERES_NO_SUITESPARSE
Sameer Agarwal7a3c43b2012-06-05 23:10:59 -0700149 factor_ = NULL;
Sameer Agarwalb0518732012-05-29 00:27:57 -0700150#endif // CERES_NO_SUITESPARSE
Petter Strandmark1e3cbd92012-08-29 09:39:56 -0700151
152#ifndef CERES_NO_CXSPARSE
153 cxsparse_factor_ = NULL;
154#endif // CERES_NO_CXSPARSE
Keir Mierle8ebb0732012-04-30 23:09:08 -0700155}
156
157SparseSchurComplementSolver::~SparseSchurComplementSolver() {
Sameer Agarwalb0518732012-05-29 00:27:57 -0700158#ifndef CERES_NO_SUITESPARSE
Sameer Agarwal7a3c43b2012-06-05 23:10:59 -0700159 if (factor_ != NULL) {
160 ss_.Free(factor_);
161 factor_ = NULL;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700162 }
Sameer Agarwalb0518732012-05-29 00:27:57 -0700163#endif // CERES_NO_SUITESPARSE
Petter Strandmark1e3cbd92012-08-29 09:39:56 -0700164
165#ifndef CERES_NO_CXSPARSE
166 if (cxsparse_factor_ != NULL) {
167 cxsparse_.Free(cxsparse_factor_);
168 cxsparse_factor_ = NULL;
169 }
170#endif // CERES_NO_CXSPARSE
Keir Mierle8ebb0732012-04-30 23:09:08 -0700171}
172
173// Determine the non-zero blocks in the Schur Complement matrix, and
174// initialize a BlockRandomAccessSparseMatrix object.
175void SparseSchurComplementSolver::InitStorage(
176 const CompressedRowBlockStructure* bs) {
Sameer Agarwal0c52f1e2012-09-17 11:30:14 -0700177 const int num_eliminate_blocks = options().elimination_groups[0];
Keir Mierle8ebb0732012-04-30 23:09:08 -0700178 const int num_col_blocks = bs->cols.size();
179 const int num_row_blocks = bs->rows.size();
180
Sameer Agarwal7a3c43b2012-06-05 23:10:59 -0700181 blocks_.resize(num_col_blocks - num_eliminate_blocks, 0);
Keir Mierle8ebb0732012-04-30 23:09:08 -0700182 for (int i = num_eliminate_blocks; i < num_col_blocks; ++i) {
Sameer Agarwal7a3c43b2012-06-05 23:10:59 -0700183 blocks_[i - num_eliminate_blocks] = bs->cols[i].size;
Keir Mierle8ebb0732012-04-30 23:09:08 -0700184 }
185
186 set<pair<int, int> > block_pairs;
Sameer Agarwal7a3c43b2012-06-05 23:10:59 -0700187 for (int i = 0; i < blocks_.size(); ++i) {
Keir Mierle8ebb0732012-04-30 23:09:08 -0700188 block_pairs.insert(make_pair(i, i));
189 }
190
191 int r = 0;
192 while (r < num_row_blocks) {
193 int e_block_id = bs->rows[r].cells.front().block_id;
194 if (e_block_id >= num_eliminate_blocks) {
195 break;
196 }
197 vector<int> f_blocks;
198
199 // Add to the chunk until the first block in the row is
200 // different than the one in the first row for the chunk.
201 for (; r < num_row_blocks; ++r) {
202 const CompressedRow& row = bs->rows[r];
203 if (row.cells.front().block_id != e_block_id) {
204 break;
205 }
206
207 // Iterate over the blocks in the row, ignoring the first
208 // block since it is the one to be eliminated.
209 for (int c = 1; c < row.cells.size(); ++c) {
210 const Cell& cell = row.cells[c];
211 f_blocks.push_back(cell.block_id - num_eliminate_blocks);
212 }
213 }
214
215 sort(f_blocks.begin(), f_blocks.end());
216 f_blocks.erase(unique(f_blocks.begin(), f_blocks.end()), f_blocks.end());
217 for (int i = 0; i < f_blocks.size(); ++i) {
218 for (int j = i + 1; j < f_blocks.size(); ++j) {
219 block_pairs.insert(make_pair(f_blocks[i], f_blocks[j]));
220 }
221 }
222 }
223
224 // Remaing rows do not contribute to the chunks and directly go
225 // into the schur complement via an outer product.
226 for (; r < num_row_blocks; ++r) {
227 const CompressedRow& row = bs->rows[r];
228 CHECK_GE(row.cells.front().block_id, num_eliminate_blocks);
229 for (int i = 0; i < row.cells.size(); ++i) {
230 int r_block1_id = row.cells[i].block_id - num_eliminate_blocks;
231 for (int j = 0; j < row.cells.size(); ++j) {
232 int r_block2_id = row.cells[j].block_id - num_eliminate_blocks;
233 if (r_block1_id <= r_block2_id) {
234 block_pairs.insert(make_pair(r_block1_id, r_block2_id));
235 }
236 }
237 }
238 }
239
Sameer Agarwal7a3c43b2012-06-05 23:10:59 -0700240 set_lhs(new BlockRandomAccessSparseMatrix(blocks_, block_pairs));
Keir Mierle8ebb0732012-04-30 23:09:08 -0700241 set_rhs(new double[lhs()->num_rows()]);
242}
243
Sameer Agarwalb0518732012-05-29 00:27:57 -0700244bool SparseSchurComplementSolver::SolveReducedLinearSystem(double* solution) {
245 switch (options().sparse_linear_algebra_library) {
246 case SUITE_SPARSE:
247 return SolveReducedLinearSystemUsingSuiteSparse(solution);
248 case CX_SPARSE:
249 return SolveReducedLinearSystemUsingCXSparse(solution);
250 default:
251 LOG(FATAL) << "Unknown sparse linear algebra library : "
252 << options().sparse_linear_algebra_library;
253 }
254
255 LOG(FATAL) << "Unknown sparse linear algebra library : "
256 << options().sparse_linear_algebra_library;
257 return false;
258}
259
260#ifndef CERES_NO_SUITESPARSE
Keir Mierle8ebb0732012-04-30 23:09:08 -0700261// Solve the system Sx = r, assuming that the matrix S is stored in a
262// BlockRandomAccessSparseMatrix. The linear system is solved using
263// CHOLMOD's sparse cholesky factorization routines.
Sameer Agarwalb0518732012-05-29 00:27:57 -0700264bool SparseSchurComplementSolver::SolveReducedLinearSystemUsingSuiteSparse(
265 double* solution) {
Keir Mierle8ebb0732012-04-30 23:09:08 -0700266 TripletSparseMatrix* tsm =
267 const_cast<TripletSparseMatrix*>(
268 down_cast<const BlockRandomAccessSparseMatrix*>(lhs())->matrix());
269
270 const int num_rows = tsm->num_rows();
271
272 // The case where there are no f blocks, and the system is block
273 // diagonal.
274 if (num_rows == 0) {
275 return true;
276 }
277
278 cholmod_sparse* cholmod_lhs = ss_.CreateSparseMatrix(tsm);
279 // The matrix is symmetric, and the upper triangular part of the
280 // matrix contains the values.
281 cholmod_lhs->stype = 1;
282
283 cholmod_dense* cholmod_rhs =
284 ss_.CreateDenseVector(const_cast<double*>(rhs()), num_rows, num_rows);
285
286 // Symbolic factorization is computed if we don't already have one handy.
Sameer Agarwal7a3c43b2012-06-05 23:10:59 -0700287 if (factor_ == NULL) {
288 if (options().use_block_amd) {
289 factor_ = ss_.BlockAnalyzeCholesky(cholmod_lhs, blocks_, blocks_);
290 } else {
291 factor_ = ss_.AnalyzeCholesky(cholmod_lhs);
292 }
Keir Mierle8ebb0732012-04-30 23:09:08 -0700293
Sameer Agarwalcb83b282012-06-06 22:26:09 -0700294 if (VLOG_IS_ON(2)) {
295 cholmod_print_common("Symbolic Analysis", ss_.mutable_cc());
296 }
Sameer Agarwal7a3c43b2012-06-05 23:10:59 -0700297 }
298
299 CHECK_NOTNULL(factor_);
300
Keir Mierle8ebb0732012-04-30 23:09:08 -0700301 cholmod_dense* cholmod_solution =
Sameer Agarwal7a3c43b2012-06-05 23:10:59 -0700302 ss_.SolveCholesky(cholmod_lhs, factor_, cholmod_rhs);
303
Keir Mierle8ebb0732012-04-30 23:09:08 -0700304 ss_.Free(cholmod_lhs);
305 cholmod_lhs = NULL;
306 ss_.Free(cholmod_rhs);
307 cholmod_rhs = NULL;
308
Keir Mierle8ebb0732012-04-30 23:09:08 -0700309 if (cholmod_solution == NULL) {
Sameer Agarwalbdabc292012-11-07 11:26:32 -0800310 LOG(WARNING) << "CHOLMOD solve failed.";
Keir Mierle8ebb0732012-04-30 23:09:08 -0700311 return false;
312 }
313
314 VectorRef(solution, num_rows)
315 = VectorRef(static_cast<double*>(cholmod_solution->x), num_rows);
316 ss_.Free(cholmod_solution);
317 return true;
318}
Sameer Agarwalb0518732012-05-29 00:27:57 -0700319#else
320bool SparseSchurComplementSolver::SolveReducedLinearSystemUsingSuiteSparse(
321 double* solution) {
322 LOG(FATAL) << "No SuiteSparse support in Ceres.";
323 return false;
324}
Keir Mierle8ebb0732012-04-30 23:09:08 -0700325#endif // CERES_NO_SUITESPARSE
326
Sameer Agarwalb0518732012-05-29 00:27:57 -0700327#ifndef CERES_NO_CXSPARSE
328// Solve the system Sx = r, assuming that the matrix S is stored in a
329// BlockRandomAccessSparseMatrix. The linear system is solved using
330// CXSparse's sparse cholesky factorization routines.
331bool SparseSchurComplementSolver::SolveReducedLinearSystemUsingCXSparse(
332 double* solution) {
333 // Extract the TripletSparseMatrix that is used for actually storing S.
334 TripletSparseMatrix* tsm =
335 const_cast<TripletSparseMatrix*>(
336 down_cast<const BlockRandomAccessSparseMatrix*>(lhs())->matrix());
337
338 const int num_rows = tsm->num_rows();
339
340 // The case where there are no f blocks, and the system is block
341 // diagonal.
342 if (num_rows == 0) {
343 return true;
344 }
345
Petter Strandmark1e3cbd92012-08-29 09:39:56 -0700346 cs_di* lhs = CHECK_NOTNULL(cxsparse_.CreateSparseMatrix(tsm));
Sameer Agarwalb0518732012-05-29 00:27:57 -0700347 VectorRef(solution, num_rows) = ConstVectorRef(rhs(), num_rows);
348
Petter Strandmark1e3cbd92012-08-29 09:39:56 -0700349 // Compute symbolic factorization if not available.
350 if (cxsparse_factor_ == NULL) {
351 cxsparse_factor_ = CHECK_NOTNULL(cxsparse_.AnalyzeCholesky(lhs));
352 }
353
354 // Solve the linear system.
355 bool ok = cxsparse_.SolveCholesky(lhs, cxsparse_factor_, solution);
356
357 cxsparse_.Free(lhs);
Sameer Agarwalb0518732012-05-29 00:27:57 -0700358 return ok;
359}
360#else
361bool SparseSchurComplementSolver::SolveReducedLinearSystemUsingCXSparse(
362 double* solution) {
363 LOG(FATAL) << "No CXSparse support in Ceres.";
364 return false;
365}
366#endif // CERES_NO_CXPARSE
367
Keir Mierle8ebb0732012-04-30 23:09:08 -0700368} // namespace internal
369} // namespace ceres