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