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/iterative_schur_complement_solver.h"
     32 
     33 #include <algorithm>
     34 #include <cstring>
     35 #include <vector>
     36 
     37 #include "Eigen/Dense"
     38 #include "ceres/block_sparse_matrix.h"
     39 #include "ceres/block_structure.h"
     40 #include "ceres/conjugate_gradients_solver.h"
     41 #include "ceres/implicit_schur_complement.h"
     42 #include "ceres/internal/eigen.h"
     43 #include "ceres/internal/scoped_ptr.h"
     44 #include "ceres/linear_solver.h"
     45 #include "ceres/preconditioner.h"
     46 #include "ceres/schur_jacobi_preconditioner.h"
     47 #include "ceres/triplet_sparse_matrix.h"
     48 #include "ceres/types.h"
     49 #include "ceres/visibility_based_preconditioner.h"
     50 #include "ceres/wall_time.h"
     51 #include "glog/logging.h"
     52 
     53 namespace ceres {
     54 namespace internal {
     55 
     56 IterativeSchurComplementSolver::IterativeSchurComplementSolver(
     57     const LinearSolver::Options& options)
     58     : options_(options) {
     59 }
     60 
     61 IterativeSchurComplementSolver::~IterativeSchurComplementSolver() {
     62 }
     63 
     64 LinearSolver::Summary IterativeSchurComplementSolver::SolveImpl(
     65     BlockSparseMatrix* A,
     66     const double* b,
     67     const LinearSolver::PerSolveOptions& per_solve_options,
     68     double* x) {
     69   EventLogger event_logger("IterativeSchurComplementSolver::Solve");
     70 
     71   CHECK_NOTNULL(A->block_structure());
     72 
     73   // Initialize a ImplicitSchurComplement object.
     74   if (schur_complement_ == NULL) {
     75     schur_complement_.reset(
     76         new ImplicitSchurComplement(options_.elimination_groups[0],
     77                                     options_.preconditioner_type == JACOBI));
     78   }
     79   schur_complement_->Init(*A, per_solve_options.D, b);
     80 
     81   const int num_schur_complement_blocks =
     82       A->block_structure()->cols.size() - options_.elimination_groups[0];
     83   if (num_schur_complement_blocks == 0) {
     84     VLOG(2) << "No parameter blocks left in the schur complement.";
     85     LinearSolver::Summary cg_summary;
     86     cg_summary.num_iterations = 0;
     87     cg_summary.termination_type = TOLERANCE;
     88     schur_complement_->BackSubstitute(NULL, x);
     89     return cg_summary;
     90   }
     91 
     92   // Initialize the solution to the Schur complement system to zero.
     93   //
     94   // TODO(sameeragarwal): There maybe a better initialization than an
     95   // all zeros solution. Explore other cheap starting points.
     96   reduced_linear_system_solution_.resize(schur_complement_->num_rows());
     97   reduced_linear_system_solution_.setZero();
     98 
     99   // Instantiate a conjugate gradient solver that runs on the Schur complement
    100   // matrix with the block diagonal of the matrix F'F as the preconditioner.
    101   LinearSolver::Options cg_options;
    102   cg_options.max_num_iterations = options_.max_num_iterations;
    103   ConjugateGradientsSolver cg_solver(cg_options);
    104   LinearSolver::PerSolveOptions cg_per_solve_options;
    105 
    106   cg_per_solve_options.r_tolerance = per_solve_options.r_tolerance;
    107   cg_per_solve_options.q_tolerance = per_solve_options.q_tolerance;
    108 
    109   Preconditioner::Options preconditioner_options;
    110   preconditioner_options.type = options_.preconditioner_type;
    111   preconditioner_options.sparse_linear_algebra_library_type =
    112       options_.sparse_linear_algebra_library_type;
    113   preconditioner_options.num_threads = options_.num_threads;
    114   preconditioner_options.row_block_size = options_.row_block_size;
    115   preconditioner_options.e_block_size = options_.e_block_size;
    116   preconditioner_options.f_block_size = options_.f_block_size;
    117   preconditioner_options.elimination_groups = options_.elimination_groups;
    118 
    119   switch (options_.preconditioner_type) {
    120     case IDENTITY:
    121       break;
    122     case JACOBI:
    123       preconditioner_.reset(
    124           new SparseMatrixPreconditionerWrapper(
    125               schur_complement_->block_diagonal_FtF_inverse()));
    126       break;
    127     case SCHUR_JACOBI:
    128       if (preconditioner_.get() == NULL) {
    129         preconditioner_.reset(
    130             new SchurJacobiPreconditioner(*A->block_structure(),
    131                                           preconditioner_options));
    132       }
    133       break;
    134     case CLUSTER_JACOBI:
    135     case CLUSTER_TRIDIAGONAL:
    136       if (preconditioner_.get() == NULL) {
    137         preconditioner_.reset(
    138             new VisibilityBasedPreconditioner(*A->block_structure(),
    139                                               preconditioner_options));
    140       }
    141       break;
    142     default:
    143       LOG(FATAL) << "Unknown Preconditioner Type";
    144   }
    145 
    146   bool preconditioner_update_was_successful = true;
    147   if (preconditioner_.get() != NULL) {
    148     preconditioner_update_was_successful =
    149         preconditioner_->Update(*A, per_solve_options.D);
    150     cg_per_solve_options.preconditioner = preconditioner_.get();
    151   }
    152 
    153   event_logger.AddEvent("Setup");
    154 
    155   LinearSolver::Summary cg_summary;
    156   cg_summary.num_iterations = 0;
    157   cg_summary.termination_type = FAILURE;
    158 
    159   if (preconditioner_update_was_successful) {
    160     cg_summary = cg_solver.Solve(schur_complement_.get(),
    161                                  schur_complement_->rhs().data(),
    162                                  cg_per_solve_options,
    163                                  reduced_linear_system_solution_.data());
    164     if (cg_summary.termination_type != FAILURE) {
    165       schur_complement_->BackSubstitute(
    166           reduced_linear_system_solution_.data(), x);
    167     }
    168   }
    169 
    170   VLOG(2) << "CG Iterations : " << cg_summary.num_iterations;
    171 
    172   event_logger.AddEvent("Solve");
    173   return cg_summary;
    174 }
    175 
    176 }  // namespace internal
    177 }  // namespace ceres
    178