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