aboutsummaryrefslogtreecommitdiff
path: root/internal/ceres/sparse_normal_cholesky_solver.cc
diff options
context:
space:
mode:
authorSascha Haeberling <haeberling@google.com>2013-07-23 19:00:21 -0700
committerSascha Haeberling <haeberling@google.com>2013-07-24 12:00:09 -0700
commit1d2624a10e2c559f8ba9ef89eaa30832c0a83a96 (patch)
treef43667ef858dd0f377b15a58a9d5c9a126762c55 /internal/ceres/sparse_normal_cholesky_solver.cc
parent0ae28bd5885b5daa526898fcf7c323dc2c3e1963 (diff)
downloadceres-solver-1d2624a10e2c559f8ba9ef89eaa30832c0a83a96.tar.gz
Update ceres to the latest version in google3.
Change-Id: I0165fffa55f60714f23e0096eac89fa68df75a05
Diffstat (limited to 'internal/ceres/sparse_normal_cholesky_solver.cc')
-rw-r--r--internal/ceres/sparse_normal_cholesky_solver.cc70
1 files changed, 39 insertions, 31 deletions
diff --git a/internal/ceres/sparse_normal_cholesky_solver.cc b/internal/ceres/sparse_normal_cholesky_solver.cc
index 9e00b44..9601142 100644
--- a/internal/ceres/sparse_normal_cholesky_solver.cc
+++ b/internal/ceres/sparse_normal_cholesky_solver.cc
@@ -28,6 +28,8 @@
//
// Author: sameeragarwal@google.com (Sameer Agarwal)
+#if !defined(CERES_NO_SUITESPARSE) || !defined(CERES_NO_CXSPARSE)
+
#include "ceres/sparse_normal_cholesky_solver.h"
#include <algorithm>
@@ -39,12 +41,13 @@
#endif
#include "ceres/compressed_row_sparse_matrix.h"
+#include "ceres/internal/eigen.h"
+#include "ceres/internal/scoped_ptr.h"
#include "ceres/linear_solver.h"
#include "ceres/suitesparse.h"
#include "ceres/triplet_sparse_matrix.h"
-#include "ceres/internal/eigen.h"
-#include "ceres/internal/scoped_ptr.h"
#include "ceres/types.h"
+#include "ceres/wall_time.h"
namespace ceres {
namespace internal {
@@ -103,6 +106,8 @@ LinearSolver::Summary SparseNormalCholeskySolver::SolveImplUsingCXSparse(
const double* b,
const LinearSolver::PerSolveOptions& per_solve_options,
double * x) {
+ EventLogger event_logger("SparseNormalCholeskySolver::CXSparse::Solve");
+
LinearSolver::Summary summary;
summary.num_iterations = 1;
const int num_cols = A->num_cols();
@@ -128,26 +133,38 @@ LinearSolver::Summary SparseNormalCholeskySolver::SolveImplUsingCXSparse(
// factorized. CHOLMOD/SuiteSparse on the other hand can just work
// off of Jt to compute the Cholesky factorization of the normal
// equations.
- cs_di* A2 = cs_transpose(&At, 1);
- cs_di* AtA = cs_multiply(&At,A2);
+ cs_di* A2 = cxsparse_.TransposeMatrix(&At);
+ cs_di* AtA = cxsparse_.MatrixMatrixMultiply(&At, A2);
cxsparse_.Free(A2);
if (per_solve_options.D != NULL) {
A->DeleteRows(num_cols);
}
+ event_logger.AddEvent("Setup");
// Compute symbolic factorization if not available.
if (cxsparse_factor_ == NULL) {
- cxsparse_factor_ = CHECK_NOTNULL(cxsparse_.AnalyzeCholesky(AtA));
+ if (options_.use_postordering) {
+ cxsparse_factor_ =
+ CHECK_NOTNULL(cxsparse_.BlockAnalyzeCholesky(AtA,
+ A->col_blocks(),
+ A->col_blocks()));
+ } else {
+ cxsparse_factor_ =
+ CHECK_NOTNULL(cxsparse_.AnalyzeCholeskyWithNaturalOrdering(AtA));
+ }
}
+ event_logger.AddEvent("Analysis");
// Solve the linear system.
if (cxsparse_.SolveCholesky(AtA, cxsparse_factor_, Atb.data())) {
VectorRef(x, Atb.rows()) = Atb;
summary.termination_type = TOLERANCE;
}
+ event_logger.AddEvent("Solve");
cxsparse_.Free(AtA);
+ event_logger.AddEvent("Teardown");
return summary;
}
#else
@@ -169,9 +186,9 @@ LinearSolver::Summary SparseNormalCholeskySolver::SolveImplUsingSuiteSparse(
const double* b,
const LinearSolver::PerSolveOptions& per_solve_options,
double * x) {
- const time_t start_time = time(NULL);
- const int num_cols = A->num_cols();
+ EventLogger event_logger("SparseNormalCholeskySolver::SuiteSparse::Solve");
+ const int num_cols = A->num_cols();
LinearSolver::Summary summary;
Vector Atb = Vector::Zero(num_cols);
A->LeftMultiply(b, Atb.data());
@@ -185,32 +202,26 @@ LinearSolver::Summary SparseNormalCholeskySolver::SolveImplUsingSuiteSparse(
VectorRef(x, num_cols).setZero();
- scoped_ptr<cholmod_sparse> lhs(ss_.CreateSparseMatrixTransposeView(A));
- CHECK_NOTNULL(lhs.get());
-
+ cholmod_sparse lhs = ss_.CreateSparseMatrixTransposeView(A);
cholmod_dense* rhs = ss_.CreateDenseVector(Atb.data(), num_cols, num_cols);
- const time_t init_time = time(NULL);
+ event_logger.AddEvent("Setup");
if (factor_ == NULL) {
- if (options_.use_block_amd) {
- factor_ = ss_.BlockAnalyzeCholesky(lhs.get(),
- A->col_blocks(),
- A->row_blocks());
+ if (options_.use_postordering) {
+ factor_ =
+ CHECK_NOTNULL(ss_.BlockAnalyzeCholesky(&lhs,
+ A->col_blocks(),
+ A->row_blocks()));
} else {
- factor_ = ss_.AnalyzeCholesky(lhs.get());
- }
-
- if (VLOG_IS_ON(2)) {
- cholmod_print_common("Symbolic Analysis", ss_.mutable_cc());
+ factor_ =
+ CHECK_NOTNULL(ss_.AnalyzeCholeskyWithNaturalOrdering(&lhs));
}
}
- CHECK_NOTNULL(factor_);
-
- const time_t symbolic_time = time(NULL);
+ event_logger.AddEvent("Analysis");
- cholmod_dense* sol = ss_.SolveCholesky(lhs.get(), factor_, rhs);
- const time_t solve_time = time(NULL);
+ cholmod_dense* sol = ss_.SolveCholesky(&lhs, factor_, rhs);
+ event_logger.AddEvent("Solve");
ss_.Free(rhs);
rhs = NULL;
@@ -228,12 +239,7 @@ LinearSolver::Summary SparseNormalCholeskySolver::SolveImplUsingSuiteSparse(
summary.termination_type = TOLERANCE;
}
- const time_t cleanup_time = time(NULL);
- VLOG(2) << "time (sec) total: " << (cleanup_time - start_time)
- << " init: " << (init_time - start_time)
- << " symbolic: " << (symbolic_time - init_time)
- << " solve: " << (solve_time - symbolic_time)
- << " cleanup: " << (cleanup_time - solve_time);
+ event_logger.AddEvent("Teardown");
return summary;
}
#else
@@ -251,3 +257,5 @@ LinearSolver::Summary SparseNormalCholeskySolver::SolveImplUsingSuiteSparse(
} // namespace internal
} // namespace ceres
+
+#endif // !defined(CERES_NO_SUITESPARSE) || !defined(CERES_NO_CXSPARSE)