Add a dense Cholesky factorization based linear solver.

For problems with a small number of variables, but a large
number of residuals, it is sometimes beneficial to use the
Cholesky factorization on the normal equations, instead of
the dense QR factorization of the Jacobian, even though it
is numerically the better thing to do.

Change-Id: I3506b006195754018deec964e6e190b7e8c9ac8f
diff --git a/internal/ceres/CMakeLists.txt b/internal/ceres/CMakeLists.txt
index ef0bd68..a914005 100644
--- a/internal/ceres/CMakeLists.txt
+++ b/internal/ceres/CMakeLists.txt
@@ -45,6 +45,7 @@
     conditioned_cost_function.cc
     conjugate_gradients_solver.cc
     corrector.cc
+    dense_normal_cholesky_solver.cc
     dense_qr_solver.cc
     dense_sparse_matrix.cc
     detect_structure.cc
diff --git a/internal/ceres/dense_normal_cholesky_solver.cc b/internal/ceres/dense_normal_cholesky_solver.cc
new file mode 100644
index 0000000..f6bb99a
--- /dev/null
+++ b/internal/ceres/dense_normal_cholesky_solver.cc
@@ -0,0 +1,86 @@
+// Ceres Solver - A fast non-linear least squares minimizer
+// Copyright 2012 Google Inc. All rights reserved.
+// http://code.google.com/p/ceres-solver/
+//
+// Redistribution and use in source and binary forms, with or without
+// modification, are permitted provided that the following conditions are met:
+//
+// * Redistributions of source code must retain the above copyright notice,
+//   this list of conditions and the following disclaimer.
+// * Redistributions in binary form must reproduce the above copyright notice,
+//   this list of conditions and the following disclaimer in the documentation
+//   and/or other materials provided with the distribution.
+// * Neither the name of Google Inc. nor the names of its contributors may be
+//   used to endorse or promote products derived from this software without
+//   specific prior written permission.
+//
+// THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
+// AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
+// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
+// ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE
+// LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
+// CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
+// SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
+// INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
+// CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
+// ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
+// POSSIBILITY OF SUCH DAMAGE.
+//
+// Author: sameeragarwal@google.com (Sameer Agarwal)
+
+#include "ceres/dense_normal_cholesky_solver.h"
+
+#include <cstddef>
+
+#include "Eigen/Dense"
+#include "ceres/dense_sparse_matrix.h"
+#include "ceres/linear_solver.h"
+#include "ceres/internal/eigen.h"
+#include "ceres/internal/scoped_ptr.h"
+#include "ceres/types.h"
+
+namespace ceres {
+namespace internal {
+
+DenseNormalCholeskySolver::DenseNormalCholeskySolver(
+    const LinearSolver::Options& options)
+    : options_(options) {}
+
+LinearSolver::Summary DenseNormalCholeskySolver::SolveImpl(
+    DenseSparseMatrix* A,
+    const double* b,
+    const LinearSolver::PerSolveOptions& per_solve_options,
+    double* x) {
+  const int num_rows = A->num_rows();
+  const int num_cols = A->num_cols();
+
+  ConstAlignedMatrixRef Aref = A->matrix();
+  Matrix lhs(num_cols, num_cols);
+  lhs.setZero();
+
+  //   lhs += A'A
+  //
+  // Using rankUpdate instead of GEMM, exposes the fact that its the
+  // same matrix being multiplied with itself and that the product is
+  // symmetric.
+  lhs.selfadjointView<Eigen::Upper>().rankUpdate(Aref.transpose());
+
+  //   rhs = A'b
+  Vector rhs = Aref.transpose() * ConstVectorRef(b, num_rows);
+
+  if (per_solve_options.D != NULL) {
+    ConstVectorRef D(per_solve_options.D, num_cols);
+    lhs += D.array().square().matrix().asDiagonal();
+  }
+
+  VectorRef(x, num_cols) =
+      lhs.selfadjointView<Eigen::Upper>().ldlt().solve(rhs);
+
+  LinearSolver::Summary summary;
+  summary.num_iterations = 1;
+  summary.termination_type = TOLERANCE;
+  return summary;
+}
+
+}   // namespace internal
+}   // namespace ceres
diff --git a/internal/ceres/dense_normal_cholesky_solver.h b/internal/ceres/dense_normal_cholesky_solver.h
new file mode 100644
index 0000000..de47740
--- /dev/null
+++ b/internal/ceres/dense_normal_cholesky_solver.h
@@ -0,0 +1,95 @@
+// Ceres Solver - A fast non-linear least squares minimizer
+// Copyright 2012 Google Inc. All rights reserved.
+// http://code.google.com/p/ceres-solver/
+//
+// Redistribution and use in source and binary forms, with or without
+// modification, are permitted provided that the following conditions are met:
+//
+// * Redistributions of source code must retain the above copyright notice,
+//   this list of conditions and the following disclaimer.
+// * Redistributions in binary form must reproduce the above copyright notice,
+//   this list of conditions and the following disclaimer in the documentation
+//   and/or other materials provided with the distribution.
+// * Neither the name of Google Inc. nor the names of its contributors may be
+//   used to endorse or promote products derived from this software without
+//   specific prior written permission.
+//
+// THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
+// AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
+// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
+// ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE
+// LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
+// CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
+// SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
+// INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
+// CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
+// ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
+// POSSIBILITY OF SUCH DAMAGE.
+//
+// Author: sameeragarwal@google.com (Sameer Agarwal)
+//
+// Solve dense rectangular systems Ax = b by forming the normal
+// equations and solving them using the Cholesky factorization.
+
+#ifndef CERES_INTERNAL_DENSE_NORMAL_CHOLESKY_SOLVER_H_
+#define CERES_INTERNAL_DENSE_NORMAL_CHOLESKY_SOLVER_H_
+
+#include "ceres/linear_solver.h"
+#include "ceres/internal/macros.h"
+
+namespace ceres {
+namespace internal {
+
+class DenseSparseMatrix;
+
+// This class implements the LinearSolver interface for solving
+// rectangular/unsymmetric (well constrained) linear systems of the
+// form
+//
+//   Ax = b
+//
+// Since there does not usually exist a solution that satisfies these
+// equations, the solver instead solves the linear least squares
+// problem
+//
+//   min_x |Ax - b|^2
+//
+// Setting the gradient of the above optimization problem to zero
+// gives us the normal equations
+//
+//   A'Ax = A'b
+//
+// A'A is a positive definite matrix (hopefully), and the resulting
+// linear system can be solved using Cholesky factorization.
+//
+// If the PerSolveOptions struct has a non-null array D, then the
+// augmented/regularized linear system
+//
+//   [    A    ]x = [b]
+//   [ diag(D) ]    [0]
+//
+// is solved.
+//
+// This class uses the LDLT factorization routines from the Eigen
+// library. This solver always returns a solution, it is the user's
+// responsibility to judge if the solution is good enough for their
+// purposes.
+class DenseNormalCholeskySolver: public DenseSparseMatrixSolver {
+ public:
+  explicit DenseNormalCholeskySolver(const LinearSolver::Options& options);
+
+ private:
+  virtual LinearSolver::Summary SolveImpl(
+      DenseSparseMatrix* A,
+      const double* b,
+      const LinearSolver::PerSolveOptions& per_solve_options,
+      double* x);
+
+  const LinearSolver::Options options_;
+  CERES_DISALLOW_COPY_AND_ASSIGN(DenseNormalCholeskySolver);
+};
+
+}  // namespace internal
+}  // namespace ceres
+
+#endif  // CERES_INTERNAL_DENSE_NORMAL_CHOLESKY_SOLVER_H_
diff --git a/internal/ceres/dense_qr_solver.cc b/internal/ceres/dense_qr_solver.cc
index 21c5fb5..2b329ee 100644
--- a/internal/ceres/dense_qr_solver.cc
+++ b/internal/ceres/dense_qr_solver.cc
@@ -33,8 +33,8 @@
 #include <cstddef>
 
 #include "Eigen/Dense"
+#include "ceres/dense_sparse_matrix.h"
 #include "ceres/linear_solver.h"
-#include "ceres/triplet_sparse_matrix.h"
 #include "ceres/internal/eigen.h"
 #include "ceres/internal/scoped_ptr.h"
 #include "ceres/types.h"
diff --git a/internal/ceres/dense_qr_solver.h b/internal/ceres/dense_qr_solver.h
index 2e0b283..dd683a8 100644
--- a/internal/ceres/dense_qr_solver.h
+++ b/internal/ceres/dense_qr_solver.h
@@ -28,7 +28,7 @@
 //
 // Author: sameeragarwal@google.com (Sameer Agarwal)
 //
-// Solve dense rectangular systems Ax = b using the QR factoriztion.
+// Solve dense rectangular systems Ax = b using the QR factorization.
 #ifndef CERES_INTERNAL_DENSE_QR_SOLVER_H_
 #define CERES_INTERNAL_DENSE_QR_SOLVER_H_
 
diff --git a/internal/ceres/evaluator.cc b/internal/ceres/evaluator.cc
index 204198c..a3ce6f0 100644
--- a/internal/ceres/evaluator.cc
+++ b/internal/ceres/evaluator.cc
@@ -51,6 +51,7 @@
                              string* error) {
   switch (options.linear_solver_type) {
     case DENSE_QR:
+    case DENSE_NORMAL_CHOLESKY:
       return new ProgramEvaluator<ScratchEvaluatePreparer,
                                   DenseJacobianWriter>(options,
                                                        program);
diff --git a/internal/ceres/linear_solver.cc b/internal/ceres/linear_solver.cc
index 97dbba4..08c3ba1 100644
--- a/internal/ceres/linear_solver.cc
+++ b/internal/ceres/linear_solver.cc
@@ -31,6 +31,7 @@
 #include "ceres/linear_solver.h"
 
 #include "ceres/cgnr_solver.h"
+#include "ceres/dense_normal_cholesky_solver.h"
 #include "ceres/dense_qr_solver.h"
 #include "ceres/iterative_schur_complement_solver.h"
 #include "ceres/schur_complement_solver.h"
@@ -78,6 +79,9 @@
     case DENSE_QR:
       return new DenseQRSolver(options);
 
+    case DENSE_NORMAL_CHOLESKY:
+      return new DenseNormalCholeskySolver(options);
+
     default:
       LOG(FATAL) << "Unknown linear solver type :"
                  << options.type;
diff --git a/internal/ceres/solver_impl_test.cc b/internal/ceres/solver_impl_test.cc
index a6d6aac..76ece8d 100644
--- a/internal/ceres/solver_impl_test.cc
+++ b/internal/ceres/solver_impl_test.cc
@@ -577,6 +577,11 @@
   EXPECT_EQ(options.linear_solver_type, DENSE_QR);
   EXPECT_TRUE(solver.get() != NULL);
 
+  options.linear_solver_type = DENSE_NORMAL_CHOLESKY;
+  solver.reset(SolverImpl::CreateLinearSolver(&options, &error));
+  EXPECT_EQ(options.linear_solver_type, DENSE_NORMAL_CHOLESKY);
+  EXPECT_TRUE(solver.get() != NULL);
+
 #ifndef CERES_NO_SUITESPARSE
   options.linear_solver_type = SPARSE_NORMAL_CHOLESKY;
   options.sparse_linear_algebra_library = SUITE_SPARSE;
diff --git a/internal/ceres/system_test.cc b/internal/ceres/system_test.cc
index 7cdff21..3dfbc00 100644
--- a/internal/ceres/system_test.cc
+++ b/internal/ceres/system_test.cc
@@ -279,7 +279,8 @@
                                  sparse_linear_algebra_library,           \
                                  ordering))
 
-  CONFIGURE(DENSE_QR,    SUITE_SPARSE, NATURAL);
+  CONFIGURE(DENSE_QR, SUITE_SPARSE, NATURAL);
+  CONFIGURE(DENSE_NORMAL_CHOLESKY, SUITE_SPARSE, NATURAL);
   CONFIGURE(DENSE_SCHUR, SUITE_SPARSE, SCHUR);
 
 #ifndef CERES_NO_SUITESPARSE
diff --git a/internal/ceres/types.cc b/internal/ceres/types.cc
index 18fefad..2e950c5 100644
--- a/internal/ceres/types.cc
+++ b/internal/ceres/types.cc
@@ -37,8 +37,9 @@
 
 const char* LinearSolverTypeToString(LinearSolverType solver_type) {
   switch (solver_type) {
-    CASESTR(SPARSE_NORMAL_CHOLESKY);
+    CASESTR(DENSE_NORMAL_CHOLESKY);
     CASESTR(DENSE_QR);
+    CASESTR(SPARSE_NORMAL_CHOLESKY);
     CASESTR(DENSE_SCHUR);
     CASESTR(SPARSE_SCHUR);
     CASESTR(ITERATIVE_SCHUR);
diff --git a/internal/ceres/unsymmetric_linear_solver_test.cc b/internal/ceres/unsymmetric_linear_solver_test.cc
index f13c5a4..0b0d593 100644
--- a/internal/ceres/unsymmetric_linear_solver_test.cc
+++ b/internal/ceres/unsymmetric_linear_solver_test.cc
@@ -73,7 +73,8 @@
 
     scoped_ptr<SparseMatrix> transformed_A;
 
-    if (linear_solver_type == DENSE_QR) {
+    if (linear_solver_type == DENSE_QR ||
+        linear_solver_type == DENSE_NORMAL_CHOLESKY) {
       transformed_A.reset(new DenseSparseMatrix(*A_));
     } else if (linear_solver_type == SPARSE_NORMAL_CHOLESKY) {
       transformed_A.reset(new   CompressedRowSparseMatrix(*A_));
@@ -118,6 +119,10 @@
   TestSolver(DENSE_QR, SUITE_SPARSE);
 }
 
+TEST_F(UnsymmetricLinearSolverTest, DenseNormalCholesky) {
+  TestSolver(DENSE_NORMAL_CHOLESKY, SUITE_SPARSE);
+}
+
 #ifndef CERES_NO_SUITESPARSE
 TEST_F(UnsymmetricLinearSolverTest, SparseNormalCholeskyUsingSuiteSparse) {
   TestSolver(SPARSE_NORMAL_CHOLESKY, SUITE_SPARSE);