blob: 0d6f6bdfb88873f0558e31ce680d3f1d395d03d1 [file] [log] [blame]
Austin Schuh70cc9552019-01-21 19:46:48 -08001// Ceres Solver - A fast non-linear least squares minimizer
2// Copyright 2015 Google Inc. All rights reserved.
3// http://ceres-solver.org/
4//
5// Redistribution and use in source and binary forms, with or without
6// modification, are permitted provided that the following conditions are met:
7//
8// * Redistributions of source code must retain the above copyright notice,
9// this list of conditions and the following disclaimer.
10// * Redistributions in binary form must reproduce the above copyright notice,
11// this list of conditions and the following disclaimer in the documentation
12// and/or other materials provided with the distribution.
13// * Neither the name of Google Inc. nor the names of its contributors may be
14// used to endorse or promote products derived from this software without
15// specific prior written permission.
16//
17// THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
18// AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
19// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
20// ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE
21// LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
22// CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
23// SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
24// INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
25// CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
26// ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
27// POSSIBILITY OF SUCH DAMAGE.
28//
29// Author: sameeragarwal@google.com (Sameer Agarwal)
30
31// This include must come before any #ifndef check on Ceres compile options.
32#include "ceres/internal/port.h"
33
34#ifndef CERES_NO_SUITESPARSE
Austin Schuh70cc9552019-01-21 19:46:48 -080035#include <vector>
36
37#include "ceres/compressed_col_sparse_matrix_utils.h"
38#include "ceres/compressed_row_sparse_matrix.h"
39#include "ceres/linear_solver.h"
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080040#include "ceres/suitesparse.h"
Austin Schuh70cc9552019-01-21 19:46:48 -080041#include "ceres/triplet_sparse_matrix.h"
42#include "cholmod.h"
43
44namespace ceres {
45namespace internal {
46
47using std::string;
48using std::vector;
49
50SuiteSparse::SuiteSparse() { cholmod_start(&cc_); }
51
52SuiteSparse::~SuiteSparse() { cholmod_finish(&cc_); }
53
54cholmod_sparse* SuiteSparse::CreateSparseMatrix(TripletSparseMatrix* A) {
55 cholmod_triplet triplet;
56
57 triplet.nrow = A->num_rows();
58 triplet.ncol = A->num_cols();
59 triplet.nzmax = A->max_num_nonzeros();
60 triplet.nnz = A->num_nonzeros();
61 triplet.i = reinterpret_cast<void*>(A->mutable_rows());
62 triplet.j = reinterpret_cast<void*>(A->mutable_cols());
63 triplet.x = reinterpret_cast<void*>(A->mutable_values());
64 triplet.stype = 0; // Matrix is not symmetric.
65 triplet.itype = CHOLMOD_INT;
66 triplet.xtype = CHOLMOD_REAL;
67 triplet.dtype = CHOLMOD_DOUBLE;
68
69 return cholmod_triplet_to_sparse(&triplet, triplet.nnz, &cc_);
70}
71
72cholmod_sparse* SuiteSparse::CreateSparseMatrixTranspose(
73 TripletSparseMatrix* A) {
74 cholmod_triplet triplet;
75
76 triplet.ncol = A->num_rows(); // swap row and columns
77 triplet.nrow = A->num_cols();
78 triplet.nzmax = A->max_num_nonzeros();
79 triplet.nnz = A->num_nonzeros();
80
81 // swap rows and columns
82 triplet.j = reinterpret_cast<void*>(A->mutable_rows());
83 triplet.i = reinterpret_cast<void*>(A->mutable_cols());
84 triplet.x = reinterpret_cast<void*>(A->mutable_values());
85 triplet.stype = 0; // Matrix is not symmetric.
86 triplet.itype = CHOLMOD_INT;
87 triplet.xtype = CHOLMOD_REAL;
88 triplet.dtype = CHOLMOD_DOUBLE;
89
90 return cholmod_triplet_to_sparse(&triplet, triplet.nnz, &cc_);
91}
92
93cholmod_sparse SuiteSparse::CreateSparseMatrixTransposeView(
94 CompressedRowSparseMatrix* A) {
95 cholmod_sparse m;
96 m.nrow = A->num_cols();
97 m.ncol = A->num_rows();
98 m.nzmax = A->num_nonzeros();
99 m.nz = nullptr;
100 m.p = reinterpret_cast<void*>(A->mutable_rows());
101 m.i = reinterpret_cast<void*>(A->mutable_cols());
102 m.x = reinterpret_cast<void*>(A->mutable_values());
103 m.z = nullptr;
104
105 if (A->storage_type() == CompressedRowSparseMatrix::LOWER_TRIANGULAR) {
106 m.stype = 1;
107 } else if (A->storage_type() == CompressedRowSparseMatrix::UPPER_TRIANGULAR) {
108 m.stype = -1;
109 } else {
110 m.stype = 0;
111 }
112
113 m.itype = CHOLMOD_INT;
114 m.xtype = CHOLMOD_REAL;
115 m.dtype = CHOLMOD_DOUBLE;
116 m.sorted = 1;
117 m.packed = 1;
118
119 return m;
120}
121
122cholmod_dense SuiteSparse::CreateDenseVectorView(const double* x, int size) {
123 cholmod_dense v;
124 v.nrow = size;
125 v.ncol = 1;
126 v.nzmax = size;
127 v.d = size;
128 v.x = const_cast<void*>(reinterpret_cast<const void*>(x));
129 v.xtype = CHOLMOD_REAL;
130 v.dtype = CHOLMOD_DOUBLE;
131 return v;
132}
133
134cholmod_dense* SuiteSparse::CreateDenseVector(const double* x,
135 int in_size,
136 int out_size) {
137 CHECK_LE(in_size, out_size);
138 cholmod_dense* v = cholmod_zeros(out_size, 1, CHOLMOD_REAL, &cc_);
139 if (x != nullptr) {
140 memcpy(v->x, x, in_size * sizeof(*x));
141 }
142 return v;
143}
144
145cholmod_factor* SuiteSparse::AnalyzeCholesky(cholmod_sparse* A,
146 string* message) {
147 // Cholmod can try multiple re-ordering strategies to find a fill
148 // reducing ordering. Here we just tell it use AMD with automatic
149 // matrix dependence choice of supernodal versus simplicial
150 // factorization.
151 cc_.nmethods = 1;
152 cc_.method[0].ordering = CHOLMOD_AMD;
153 cc_.supernodal = CHOLMOD_AUTO;
154
155 cholmod_factor* factor = cholmod_analyze(A, &cc_);
156 if (VLOG_IS_ON(2)) {
157 cholmod_print_common(const_cast<char*>("Symbolic Analysis"), &cc_);
158 }
159
160 if (cc_.status != CHOLMOD_OK) {
161 *message =
162 StringPrintf("cholmod_analyze failed. error code: %d", cc_.status);
163 return nullptr;
164 }
165
166 CHECK(factor != nullptr);
167 return factor;
168}
169
170cholmod_factor* SuiteSparse::BlockAnalyzeCholesky(cholmod_sparse* A,
171 const vector<int>& row_blocks,
172 const vector<int>& col_blocks,
173 string* message) {
174 vector<int> ordering;
175 if (!BlockAMDOrdering(A, row_blocks, col_blocks, &ordering)) {
176 return nullptr;
177 }
178 return AnalyzeCholeskyWithUserOrdering(A, ordering, message);
179}
180
181cholmod_factor* SuiteSparse::AnalyzeCholeskyWithUserOrdering(
182 cholmod_sparse* A, const vector<int>& ordering, string* message) {
183 CHECK_EQ(ordering.size(), A->nrow);
184
185 cc_.nmethods = 1;
186 cc_.method[0].ordering = CHOLMOD_GIVEN;
187
188 cholmod_factor* factor =
189 cholmod_analyze_p(A, const_cast<int*>(&ordering[0]), nullptr, 0, &cc_);
190 if (VLOG_IS_ON(2)) {
191 cholmod_print_common(const_cast<char*>("Symbolic Analysis"), &cc_);
192 }
193 if (cc_.status != CHOLMOD_OK) {
194 *message =
195 StringPrintf("cholmod_analyze failed. error code: %d", cc_.status);
196 return nullptr;
197 }
198
199 CHECK(factor != nullptr);
200 return factor;
201}
202
203cholmod_factor* SuiteSparse::AnalyzeCholeskyWithNaturalOrdering(
204 cholmod_sparse* A, string* message) {
205 cc_.nmethods = 1;
206 cc_.method[0].ordering = CHOLMOD_NATURAL;
207 cc_.postorder = 0;
208
209 cholmod_factor* factor = cholmod_analyze(A, &cc_);
210 if (VLOG_IS_ON(2)) {
211 cholmod_print_common(const_cast<char*>("Symbolic Analysis"), &cc_);
212 }
213 if (cc_.status != CHOLMOD_OK) {
214 *message =
215 StringPrintf("cholmod_analyze failed. error code: %d", cc_.status);
216 return nullptr;
217 }
218
219 CHECK(factor != nullptr);
220 return factor;
221}
222
223bool SuiteSparse::BlockAMDOrdering(const cholmod_sparse* A,
224 const vector<int>& row_blocks,
225 const vector<int>& col_blocks,
226 vector<int>* ordering) {
227 const int num_row_blocks = row_blocks.size();
228 const int num_col_blocks = col_blocks.size();
229
230 // Arrays storing the compressed column structure of the matrix
231 // incoding the block sparsity of A.
232 vector<int> block_cols;
233 vector<int> block_rows;
234
235 CompressedColumnScalarMatrixToBlockMatrix(reinterpret_cast<const int*>(A->i),
236 reinterpret_cast<const int*>(A->p),
237 row_blocks,
238 col_blocks,
239 &block_rows,
240 &block_cols);
241 cholmod_sparse_struct block_matrix;
242 block_matrix.nrow = num_row_blocks;
243 block_matrix.ncol = num_col_blocks;
244 block_matrix.nzmax = block_rows.size();
245 block_matrix.p = reinterpret_cast<void*>(&block_cols[0]);
246 block_matrix.i = reinterpret_cast<void*>(&block_rows[0]);
247 block_matrix.x = nullptr;
248 block_matrix.stype = A->stype;
249 block_matrix.itype = CHOLMOD_INT;
250 block_matrix.xtype = CHOLMOD_PATTERN;
251 block_matrix.dtype = CHOLMOD_DOUBLE;
252 block_matrix.sorted = 1;
253 block_matrix.packed = 1;
254
255 vector<int> block_ordering(num_row_blocks);
256 if (!cholmod_amd(&block_matrix, nullptr, 0, &block_ordering[0], &cc_)) {
257 return false;
258 }
259
260 BlockOrderingToScalarOrdering(row_blocks, block_ordering, ordering);
261 return true;
262}
263
264LinearSolverTerminationType SuiteSparse::Cholesky(cholmod_sparse* A,
265 cholmod_factor* L,
266 string* message) {
267 CHECK(A != nullptr);
268 CHECK(L != nullptr);
269
270 // Save the current print level and silence CHOLMOD, otherwise
271 // CHOLMOD is prone to dumping stuff to stderr, which can be
272 // distracting when the error (matrix is indefinite) is not a fatal
273 // failure.
274 const int old_print_level = cc_.print;
275 cc_.print = 0;
276
277 cc_.quick_return_if_not_posdef = 1;
278 int cholmod_status = cholmod_factorize(A, L, &cc_);
279 cc_.print = old_print_level;
280
281 switch (cc_.status) {
282 case CHOLMOD_NOT_INSTALLED:
283 *message = "CHOLMOD failure: Method not installed.";
284 return LINEAR_SOLVER_FATAL_ERROR;
285 case CHOLMOD_OUT_OF_MEMORY:
286 *message = "CHOLMOD failure: Out of memory.";
287 return LINEAR_SOLVER_FATAL_ERROR;
288 case CHOLMOD_TOO_LARGE:
289 *message = "CHOLMOD failure: Integer overflow occurred.";
290 return LINEAR_SOLVER_FATAL_ERROR;
291 case CHOLMOD_INVALID:
292 *message = "CHOLMOD failure: Invalid input.";
293 return LINEAR_SOLVER_FATAL_ERROR;
294 case CHOLMOD_NOT_POSDEF:
295 *message = "CHOLMOD warning: Matrix not positive definite.";
296 return LINEAR_SOLVER_FAILURE;
297 case CHOLMOD_DSMALL:
298 *message =
299 "CHOLMOD warning: D for LDL' or diag(L) or "
300 "LL' has tiny absolute value.";
301 return LINEAR_SOLVER_FAILURE;
302 case CHOLMOD_OK:
303 if (cholmod_status != 0) {
304 return LINEAR_SOLVER_SUCCESS;
305 }
306
307 *message =
308 "CHOLMOD failure: cholmod_factorize returned false "
309 "but cholmod_common::status is CHOLMOD_OK."
310 "Please report this to ceres-solver@googlegroups.com.";
311 return LINEAR_SOLVER_FATAL_ERROR;
312 default:
313 *message = StringPrintf(
314 "Unknown cholmod return code: %d. "
315 "Please report this to ceres-solver@googlegroups.com.",
316 cc_.status);
317 return LINEAR_SOLVER_FATAL_ERROR;
318 }
319
320 return LINEAR_SOLVER_FATAL_ERROR;
321}
322
323cholmod_dense* SuiteSparse::Solve(cholmod_factor* L,
324 cholmod_dense* b,
325 string* message) {
326 if (cc_.status != CHOLMOD_OK) {
327 *message = "cholmod_solve failed. CHOLMOD status is not CHOLMOD_OK";
328 return nullptr;
329 }
330
331 return cholmod_solve(CHOLMOD_A, L, b, &cc_);
332}
333
334bool SuiteSparse::ApproximateMinimumDegreeOrdering(cholmod_sparse* matrix,
335 int* ordering) {
336 return cholmod_amd(matrix, nullptr, 0, ordering, &cc_);
337}
338
339bool SuiteSparse::ConstrainedApproximateMinimumDegreeOrdering(
340 cholmod_sparse* matrix, int* constraints, int* ordering) {
341#ifndef CERES_NO_CAMD
342 return cholmod_camd(matrix, nullptr, 0, constraints, ordering, &cc_);
343#else
344 LOG(FATAL) << "Congratulations you have found a bug in Ceres."
345 << "Ceres Solver was compiled with SuiteSparse "
346 << "version 4.1.0 or less. Calling this function "
347 << "in that case is a bug. Please contact the"
348 << "the Ceres Solver developers.";
349 return false;
350#endif
351}
352
353std::unique_ptr<SparseCholesky> SuiteSparseCholesky::Create(
354 const OrderingType ordering_type) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800355 return std::unique_ptr<SparseCholesky>(
356 new SuiteSparseCholesky(ordering_type));
Austin Schuh70cc9552019-01-21 19:46:48 -0800357}
358
359SuiteSparseCholesky::SuiteSparseCholesky(const OrderingType ordering_type)
360 : ordering_type_(ordering_type), factor_(nullptr) {}
361
362SuiteSparseCholesky::~SuiteSparseCholesky() {
363 if (factor_ != nullptr) {
364 ss_.Free(factor_);
365 }
366}
367
368LinearSolverTerminationType SuiteSparseCholesky::Factorize(
369 CompressedRowSparseMatrix* lhs, string* message) {
370 if (lhs == nullptr) {
371 *message = "Failure: Input lhs is NULL.";
372 return LINEAR_SOLVER_FATAL_ERROR;
373 }
374
375 cholmod_sparse cholmod_lhs = ss_.CreateSparseMatrixTransposeView(lhs);
376
377 if (factor_ == nullptr) {
378 if (ordering_type_ == NATURAL) {
379 factor_ = ss_.AnalyzeCholeskyWithNaturalOrdering(&cholmod_lhs, message);
380 } else {
381 if (!lhs->col_blocks().empty() && !(lhs->row_blocks().empty())) {
382 factor_ = ss_.BlockAnalyzeCholesky(
383 &cholmod_lhs, lhs->col_blocks(), lhs->row_blocks(), message);
384 } else {
385 factor_ = ss_.AnalyzeCholesky(&cholmod_lhs, message);
386 }
387 }
388
389 if (factor_ == nullptr) {
390 return LINEAR_SOLVER_FATAL_ERROR;
391 }
392 }
393
394 return ss_.Cholesky(&cholmod_lhs, factor_, message);
395}
396
397CompressedRowSparseMatrix::StorageType SuiteSparseCholesky::StorageType()
398 const {
399 return ((ordering_type_ == NATURAL)
400 ? CompressedRowSparseMatrix::UPPER_TRIANGULAR
401 : CompressedRowSparseMatrix::LOWER_TRIANGULAR);
402}
403
404LinearSolverTerminationType SuiteSparseCholesky::Solve(const double* rhs,
405 double* solution,
406 string* message) {
407 // Error checking
408 if (factor_ == nullptr) {
409 *message = "Solve called without a call to Factorize first.";
410 return LINEAR_SOLVER_FATAL_ERROR;
411 }
412
413 const int num_cols = factor_->n;
414 cholmod_dense cholmod_rhs = ss_.CreateDenseVectorView(rhs, num_cols);
415 cholmod_dense* cholmod_dense_solution =
416 ss_.Solve(factor_, &cholmod_rhs, message);
417
418 if (cholmod_dense_solution == nullptr) {
419 return LINEAR_SOLVER_FAILURE;
420 }
421
422 memcpy(solution, cholmod_dense_solution->x, num_cols * sizeof(*solution));
423 ss_.Free(cholmod_dense_solution);
424 return LINEAR_SOLVER_SUCCESS;
425}
426
427} // namespace internal
428} // namespace ceres
429
430#endif // CERES_NO_SUITESPARSE