blob: 476697d3285bdbc1720fd49769d2dafedf0fbd3f [file] [log] [blame]
Austin Schuh70cc9552019-01-21 19:46:48 -08001// Ceres Solver - A fast non-linear least squares minimizer
2// Copyright 2015 Google Inc. All rights reserved.
3// http://ceres-solver.org/
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#ifndef CERES_INTERNAL_PRECONDITIONER_H_
32#define CERES_INTERNAL_PRECONDITIONER_H_
33
34#include <vector>
35#include "ceres/casts.h"
36#include "ceres/compressed_row_sparse_matrix.h"
37#include "ceres/context_impl.h"
38#include "ceres/linear_operator.h"
39#include "ceres/sparse_matrix.h"
40#include "ceres/types.h"
41
42namespace ceres {
43namespace internal {
44
45class BlockSparseMatrix;
46class SparseMatrix;
47
48class Preconditioner : public LinearOperator {
49 public:
50 struct Options {
51 PreconditionerType type = JACOBI;
52 VisibilityClusteringType visibility_clustering_type = CANONICAL_VIEWS;
53 SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type = SUITE_SPARSE;
54
55 // When using the subset preconditioner, all row blocks starting
56 // from this row block are used to construct the preconditioner.
57 //
58 // i.e., the Jacobian matrix A is horizontally partitioned as
59 //
60 // A = [P]
61 // [Q]
62 //
63 // where P has subset_preconditioner_start_row_block row blocks,
64 // and the preconditioner is the inverse of the matrix Q'Q.
65 int subset_preconditioner_start_row_block = -1;
66
67 // See solver.h for information about these flags.
68 bool use_postordering = false;
69
70 // If possible, how many threads the preconditioner can use.
71 int num_threads = 1;
72
73 // Hints about the order in which the parameter blocks should be
74 // eliminated by the linear solver.
75 //
76 // For example if elimination_groups is a vector of size k, then
77 // the linear solver is informed that it should eliminate the
78 // parameter blocks 0 ... elimination_groups[0] - 1 first, and
79 // then elimination_groups[0] ... elimination_groups[1] - 1 and so
80 // on. Within each elimination group, the linear solver is free to
81 // choose how the parameter blocks are ordered. Different linear
82 // solvers have differing requirements on elimination_groups.
83 //
84 // The most common use is for Schur type solvers, where there
85 // should be at least two elimination groups and the first
86 // elimination group must form an independent set in the normal
87 // equations. The first elimination group corresponds to the
88 // num_eliminate_blocks in the Schur type solvers.
89 std::vector<int> elimination_groups;
90
91 // If the block sizes in a BlockSparseMatrix are fixed, then in
92 // some cases the Schur complement based solvers can detect and
93 // specialize on them.
94 //
95 // It is expected that these parameters are set programmatically
96 // rather than manually.
97 //
98 // Please see schur_complement_solver.h and schur_eliminator.h for
99 // more details.
100 int row_block_size = Eigen::Dynamic;
101 int e_block_size = Eigen::Dynamic;
102 int f_block_size = Eigen::Dynamic;
103
104 ContextImpl* context = nullptr;
105 };
106
107 // If the optimization problem is such that there are no remaining
108 // e-blocks, ITERATIVE_SCHUR with a Schur type preconditioner cannot
109 // be used. This function returns JACOBI if a preconditioner for
110 // ITERATIVE_SCHUR is used. The input preconditioner_type is
111 // returned otherwise.
112 static PreconditionerType PreconditionerForZeroEBlocks(
113 PreconditionerType preconditioner_type);
114
115 virtual ~Preconditioner();
116
117 // Update the numerical value of the preconditioner for the linear
118 // system:
119 //
120 // | A | x = |b|
121 // |diag(D)| |0|
122 //
123 // for some vector b. It is important that the matrix A have the
124 // same block structure as the one used to construct this object.
125 //
126 // D can be NULL, in which case its interpreted as a diagonal matrix
127 // of size zero.
128 virtual bool Update(const LinearOperator& A, const double* D) = 0;
129
130 // LinearOperator interface. Since the operator is symmetric,
131 // LeftMultiply and num_cols are just calls to RightMultiply and
132 // num_rows respectively. Update() must be called before
133 // RightMultiply can be called.
134 virtual void RightMultiply(const double* x, double* y) const = 0;
135 virtual void LeftMultiply(const double* x, double* y) const {
136 return RightMultiply(x, y);
137 }
138
139 virtual int num_rows() const = 0;
140 virtual int num_cols() const {
141 return num_rows();
142 }
143};
144
145// This templated subclass of Preconditioner serves as a base class for
146// other preconditioners that depend on the particular matrix layout of
147// the underlying linear operator.
148template <typename MatrixType>
149class TypedPreconditioner : public Preconditioner {
150 public:
151 virtual ~TypedPreconditioner() {}
152 virtual bool Update(const LinearOperator& A, const double* D) {
153 return UpdateImpl(*down_cast<const MatrixType*>(&A), D);
154 }
155
156 private:
157 virtual bool UpdateImpl(const MatrixType& A, const double* D) = 0;
158};
159
160// Preconditioners that depend on access to the low level structure
161// of a SparseMatrix.
162typedef TypedPreconditioner<SparseMatrix> SparseMatrixPreconditioner; // NOLINT
163typedef TypedPreconditioner<BlockSparseMatrix> BlockSparseMatrixPreconditioner; // NOLINT
164typedef TypedPreconditioner<CompressedRowSparseMatrix> CompressedRowSparseMatrixPreconditioner; // NOLINT
165
166// Wrap a SparseMatrix object as a preconditioner.
167class SparseMatrixPreconditionerWrapper : public SparseMatrixPreconditioner {
168 public:
169 // Wrapper does NOT take ownership of the matrix pointer.
170 explicit SparseMatrixPreconditionerWrapper(const SparseMatrix* matrix);
171 virtual ~SparseMatrixPreconditionerWrapper();
172
173 // Preconditioner interface
174 virtual void RightMultiply(const double* x, double* y) const;
175 virtual int num_rows() const;
176
177 private:
178 virtual bool UpdateImpl(const SparseMatrix& A, const double* D);
179 const SparseMatrix* matrix_;
180};
181
182} // namespace internal
183} // namespace ceres
184
185#endif // CERES_INTERNAL_PRECONDITIONER_H_