blob: 5a082de1932ef8c687d26a0c11f2afc72de9e71f [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
Keir Mierle8ebb0732012-04-30 23:09:08 -070031#include "ceres/sparse_normal_cholesky_solver.h"
32
33#include <algorithm>
34#include <cstring>
35#include <ctime>
Sameer Agarwalb0518732012-05-29 00:27:57 -070036
37#ifndef CERES_NO_CXSPARSE
38#include "cs.h"
39#endif
40
Keir Mierle8ebb0732012-04-30 23:09:08 -070041#include "ceres/compressed_row_sparse_matrix.h"
42#include "ceres/linear_solver.h"
43#include "ceres/suitesparse.h"
44#include "ceres/triplet_sparse_matrix.h"
45#include "ceres/internal/eigen.h"
46#include "ceres/internal/scoped_ptr.h"
47#include "ceres/types.h"
48
Sameer Agarwalb0518732012-05-29 00:27:57 -070049
50
Keir Mierle8ebb0732012-04-30 23:09:08 -070051namespace ceres {
52namespace internal {
53
54SparseNormalCholeskySolver::SparseNormalCholeskySolver(
55 const LinearSolver::Options& options)
Sameer Agarwalb0518732012-05-29 00:27:57 -070056 : options_(options) {
57#ifndef CERES_NO_SUITESPARSE
58 symbolic_factor_ = NULL;
59#endif
60}
Keir Mierle8ebb0732012-04-30 23:09:08 -070061
62SparseNormalCholeskySolver::~SparseNormalCholeskySolver() {
Sameer Agarwalb0518732012-05-29 00:27:57 -070063#ifndef CERES_NO_SUITESPARSE
Keir Mierle8ebb0732012-04-30 23:09:08 -070064 if (symbolic_factor_ != NULL) {
65 ss_.Free(symbolic_factor_);
66 symbolic_factor_ = NULL;
67 }
Sameer Agarwalb0518732012-05-29 00:27:57 -070068#endif
Keir Mierle8ebb0732012-04-30 23:09:08 -070069}
70
71LinearSolver::Summary SparseNormalCholeskySolver::SolveImpl(
72 CompressedRowSparseMatrix* A,
73 const double* b,
74 const LinearSolver::PerSolveOptions& per_solve_options,
75 double * x) {
Sameer Agarwalb0518732012-05-29 00:27:57 -070076 switch (options_.sparse_linear_algebra_library) {
77 case SUITE_SPARSE:
78 return SolveImplUsingSuiteSparse(A, b, per_solve_options, x);
79 case CX_SPARSE:
80 return SolveImplUsingCXSparse(A, b, per_solve_options, x);
81 default:
82 LOG(FATAL) << "Unknown sparse linear algebra library : "
83 << options_.sparse_linear_algebra_library;
84 }
85
86 LOG(FATAL) << "Unknown sparse linear algebra library : "
87 << options_.sparse_linear_algebra_library;
88 return LinearSolver::Summary();
89}
90
91#ifndef CERES_NO_CXSPARSE
92LinearSolver::Summary SparseNormalCholeskySolver::SolveImplUsingCXSparse(
93 CompressedRowSparseMatrix* A,
94 const double* b,
95 const LinearSolver::PerSolveOptions& per_solve_options,
96 double * x) {
97 LinearSolver::Summary summary;
98 summary.num_iterations = 1;
99 const int num_cols = A->num_cols();
100 Vector Atb = Vector::Zero(num_cols);
101 A->LeftMultiply(b, Atb.data());
102
103 if (per_solve_options.D != NULL) {
104 // Temporarily append a diagonal block to the A matrix, but undo
105 // it before returning the matrix to the user.
106 CompressedRowSparseMatrix D(per_solve_options.D, num_cols);
107 A->AppendRows(D);
108 }
109
110 VectorRef(x, num_cols).setZero();
111
112 // Wrap the augmented Jacobian in a compressed sparse column matrix.
113 cs_di At;
114 At.m = A->num_cols();
115 At.n = A->num_rows();
116 At.nz = -1;
117 At.nzmax = A->num_nonzeros();
118 At.p = A->mutable_rows();
119 At.i = A->mutable_cols();
120 At.x = A->mutable_values();
121
122 // Compute the normal equations. J'J delta = J'f and solve them
123 // using a sparse Cholesky factorization. Notice that when compared
124 // to SuiteSparse we have to explicitly compute the transpose of Jt,
125 // and then the normal equations before they can be
126 // factorized. CHOLMOD/SuiteSparse on the other hand can just work
127 // off of Jt to compute the Cholesky factorization of the normal
128 // equations.
129 cs_di* A2 = cs_transpose(&At, 1);
130 cs_di* AtA = cs_multiply(&At,A2);
131
132 cs_free(A2);
133 if (per_solve_options.D != NULL) {
134 A->DeleteRows(num_cols);
135 }
136
137 // This recomputes the symbolic factorization every time it is
138 // invoked. It will perhaps be worth it to cache the symbolic
139 // factorization the way we do for SuiteSparse.
140 if (cs_cholsol(1, AtA, Atb.data())) {
141 VectorRef(x, Atb.rows()) = Atb;
142 summary.termination_type = TOLERANCE;
143 }
144
145 cs_free(AtA);
146 return summary;
147}
148#else
149LinearSolver::Summary SparseNormalCholeskySolver::SolveImplUsingCXSparse(
150 CompressedRowSparseMatrix* A,
151 const double* b,
152 const LinearSolver::PerSolveOptions& per_solve_options,
153 double * x) {
154 LOG(FATAL) << "No CXSparse support in Ceres.";
155}
156#endif
157
158#ifndef CERES_NO_SUITESPARSE
159LinearSolver::Summary SparseNormalCholeskySolver::SolveImplUsingSuiteSparse(
160 CompressedRowSparseMatrix* A,
161 const double* b,
162 const LinearSolver::PerSolveOptions& per_solve_options,
163 double * x) {
Keir Mierle8ebb0732012-04-30 23:09:08 -0700164 const time_t start_time = time(NULL);
165 const int num_cols = A->num_cols();
166
167 LinearSolver::Summary summary;
168 Vector Atb = Vector::Zero(num_cols);
169 A->LeftMultiply(b, Atb.data());
170
171 if (per_solve_options.D != NULL) {
172 // Temporarily append a diagonal block to the A matrix, but undo it before
173 // returning the matrix to the user.
174 CompressedRowSparseMatrix D(per_solve_options.D, num_cols);
175 A->AppendRows(D);
176 }
177
178 VectorRef(x, num_cols).setZero();
179
180 scoped_ptr<cholmod_sparse> lhs(ss_.CreateSparseMatrixTransposeView(A));
181 CHECK_NOTNULL(lhs.get());
182
183 cholmod_dense* rhs = ss_.CreateDenseVector(Atb.data(), num_cols, num_cols);
184 const time_t init_time = time(NULL);
185
186 if (symbolic_factor_ == NULL) {
187 symbolic_factor_ = CHECK_NOTNULL(ss_.AnalyzeCholesky(lhs.get()));
188 }
189
190 const time_t symbolic_time = time(NULL);
191
192 cholmod_dense* sol = ss_.SolveCholesky(lhs.get(), symbolic_factor_, rhs);
193 const time_t solve_time = time(NULL);
194
195 ss_.Free(rhs);
196 rhs = NULL;
197
198 if (per_solve_options.D != NULL) {
199 A->DeleteRows(num_cols);
200 }
201
Keir Mierle8ebb0732012-04-30 23:09:08 -0700202 summary.num_iterations = 1;
203 if (sol != NULL) {
204 memcpy(x, sol->x, num_cols * sizeof(*x));
205
206 ss_.Free(sol);
207 sol = NULL;
208 summary.termination_type = TOLERANCE;
209 }
210
211 const time_t cleanup_time = time(NULL);
212 VLOG(2) << "time (sec) total: " << cleanup_time - start_time
213 << " init: " << init_time - start_time
214 << " symbolic: " << symbolic_time - init_time
215 << " solve: " << solve_time - symbolic_time
216 << " cleanup: " << cleanup_time - solve_time;
217 return summary;
218}
Sameer Agarwalb0518732012-05-29 00:27:57 -0700219#else
220LinearSolver::Summary SparseNormalCholeskySolver::SolveImplUsingSuiteSparse(
221 CompressedRowSparseMatrix* A,
222 const double* b,
223 const LinearSolver::PerSolveOptions& per_solve_options,
224 double * x) {
225 LOG(FATAL) << "No SuiteSparse support in Ceres.";
226}
227#endif
Keir Mierle8ebb0732012-04-30 23:09:08 -0700228
229} // namespace internal
230} // namespace ceres