Home | History | Annotate | Download | only in ceres
      1 // 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 (at) google.com (Sameer Agarwal)
     30 
     31 #include "ceres/dense_qr_solver.h"
     32 
     33 
     34 #include <cstddef>
     35 #include "Eigen/Dense"
     36 #include "ceres/dense_sparse_matrix.h"
     37 #include "ceres/internal/eigen.h"
     38 #include "ceres/internal/scoped_ptr.h"
     39 #include "ceres/lapack.h"
     40 #include "ceres/linear_solver.h"
     41 #include "ceres/types.h"
     42 #include "ceres/wall_time.h"
     43 
     44 namespace ceres {
     45 namespace internal {
     46 
     47 DenseQRSolver::DenseQRSolver(const LinearSolver::Options& options)
     48     : options_(options) {
     49   work_.resize(1);
     50 }
     51 
     52 LinearSolver::Summary DenseQRSolver::SolveImpl(
     53     DenseSparseMatrix* A,
     54     const double* b,
     55     const LinearSolver::PerSolveOptions& per_solve_options,
     56     double* x) {
     57   if (options_.dense_linear_algebra_library_type == EIGEN) {
     58     return SolveUsingEigen(A, b, per_solve_options, x);
     59   } else {
     60     return SolveUsingLAPACK(A, b, per_solve_options, x);
     61   }
     62 }
     63 
     64 LinearSolver::Summary DenseQRSolver::SolveUsingLAPACK(
     65     DenseSparseMatrix* A,
     66     const double* b,
     67     const LinearSolver::PerSolveOptions& per_solve_options,
     68     double* x) {
     69   EventLogger event_logger("DenseQRSolver::Solve");
     70 
     71   const int num_rows = A->num_rows();
     72   const int num_cols = A->num_cols();
     73 
     74   if (per_solve_options.D != NULL) {
     75     // Temporarily append a diagonal block to the A matrix, but undo
     76     // it before returning the matrix to the user.
     77     A->AppendDiagonal(per_solve_options.D);
     78   }
     79 
     80   // TODO(sameeragarwal): Since we are copying anyways, the diagonal
     81   // can be appended to the matrix instead of doing it on A.
     82   lhs_ =  A->matrix();
     83 
     84   if (per_solve_options.D != NULL) {
     85     // Undo the modifications to the matrix A.
     86     A->RemoveDiagonal();
     87   }
     88 
     89   // rhs = [b;0] to account for the additional rows in the lhs.
     90   if (rhs_.rows() != lhs_.rows()) {
     91     rhs_.resize(lhs_.rows());
     92   }
     93   rhs_.setZero();
     94   rhs_.head(num_rows) = ConstVectorRef(b, num_rows);
     95 
     96   if (work_.rows() == 1) {
     97     const int work_size =
     98         LAPACK::EstimateWorkSizeForQR(lhs_.rows(), lhs_.cols());
     99     VLOG(3) << "Working memory for Dense QR factorization: "
    100             << work_size * sizeof(double);
    101     work_.resize(work_size);
    102   }
    103 
    104   LinearSolver::Summary summary;
    105   summary.num_iterations = 1;
    106   summary.termination_type = LAPACK::SolveInPlaceUsingQR(lhs_.rows(),
    107                                                          lhs_.cols(),
    108                                                          lhs_.data(),
    109                                                          work_.rows(),
    110                                                          work_.data(),
    111                                                          rhs_.data(),
    112                                                          &summary.message);
    113   event_logger.AddEvent("Solve");
    114   if (summary.termination_type == LINEAR_SOLVER_SUCCESS) {
    115     VectorRef(x, num_cols) = rhs_.head(num_cols);
    116   }
    117 
    118   event_logger.AddEvent("TearDown");
    119   return summary;
    120 }
    121 
    122 LinearSolver::Summary DenseQRSolver::SolveUsingEigen(
    123     DenseSparseMatrix* A,
    124     const double* b,
    125     const LinearSolver::PerSolveOptions& per_solve_options,
    126     double* x) {
    127   EventLogger event_logger("DenseQRSolver::Solve");
    128 
    129   const int num_rows = A->num_rows();
    130   const int num_cols = A->num_cols();
    131 
    132   if (per_solve_options.D != NULL) {
    133     // Temporarily append a diagonal block to the A matrix, but undo
    134     // it before returning the matrix to the user.
    135     A->AppendDiagonal(per_solve_options.D);
    136   }
    137 
    138   // rhs = [b;0] to account for the additional rows in the lhs.
    139   const int augmented_num_rows =
    140       num_rows + ((per_solve_options.D != NULL) ? num_cols : 0);
    141   if (rhs_.rows() != augmented_num_rows) {
    142     rhs_.resize(augmented_num_rows);
    143     rhs_.setZero();
    144   }
    145   rhs_.head(num_rows) = ConstVectorRef(b, num_rows);
    146   event_logger.AddEvent("Setup");
    147 
    148   // Solve the system.
    149   VectorRef(x, num_cols) = A->matrix().householderQr().solve(rhs_);
    150   event_logger.AddEvent("Solve");
    151 
    152   if (per_solve_options.D != NULL) {
    153     // Undo the modifications to the matrix A.
    154     A->RemoveDiagonal();
    155   }
    156 
    157   // We always succeed, since the QR solver returns the best solution
    158   // it can. It is the job of the caller to determine if the solution
    159   // is good enough or not.
    160   LinearSolver::Summary summary;
    161   summary.num_iterations = 1;
    162   summary.termination_type = LINEAR_SOLVER_SUCCESS;
    163   summary.message = "Success.";
    164 
    165   event_logger.AddEvent("TearDown");
    166   return summary;
    167 }
    168 
    169 }   // namespace internal
    170 }   // namespace ceres
    171