blob: 42dc6cc54ac4ae5103747592c7ba9c2b82423850 [file] [log] [blame]
Austin Schuh70cc9552019-01-21 19:46:48 -08001// Ceres Solver - A fast non-linear least squares minimizer
Austin Schuh3de38b02024-06-25 18:25:10 -07002// Copyright 2023 Google Inc. All rights reserved.
Austin Schuh70cc9552019-01-21 19:46:48 -08003// 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 Schuh3de38b02024-06-25 18:25:10 -070039#include "ceres/internal/disable_warnings.h"
40#include "ceres/internal/export.h"
Austin Schuh70cc9552019-01-21 19:46:48 -080041#include "ceres/linear_operator.h"
Austin Schuh3de38b02024-06-25 18:25:10 -070042#include "ceres/linear_solver.h"
Austin Schuh70cc9552019-01-21 19:46:48 -080043#include "ceres/sparse_matrix.h"
44#include "ceres/types.h"
45
Austin Schuh3de38b02024-06-25 18:25:10 -070046namespace ceres::internal {
Austin Schuh70cc9552019-01-21 19:46:48 -080047
48class BlockSparseMatrix;
49class SparseMatrix;
50
Austin Schuh3de38b02024-06-25 18:25:10 -070051class CERES_NO_EXPORT Preconditioner : public LinearOperator {
Austin Schuh70cc9552019-01-21 19:46:48 -080052 public:
53 struct Options {
Austin Schuh3de38b02024-06-25 18:25:10 -070054 Options() = default;
55 Options(const LinearSolver::Options& linear_solver_options)
56 : type(linear_solver_options.preconditioner_type),
57 visibility_clustering_type(
58 linear_solver_options.visibility_clustering_type),
59 sparse_linear_algebra_library_type(
60 linear_solver_options.sparse_linear_algebra_library_type),
61 num_threads(linear_solver_options.num_threads),
62 elimination_groups(linear_solver_options.elimination_groups),
63 row_block_size(linear_solver_options.row_block_size),
64 e_block_size(linear_solver_options.e_block_size),
65 f_block_size(linear_solver_options.f_block_size),
66 context(linear_solver_options.context) {}
67
Austin Schuh70cc9552019-01-21 19:46:48 -080068 PreconditionerType type = JACOBI;
69 VisibilityClusteringType visibility_clustering_type = CANONICAL_VIEWS;
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080070 SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type =
71 SUITE_SPARSE;
Austin Schuh3de38b02024-06-25 18:25:10 -070072 OrderingType ordering_type = OrderingType::NATURAL;
Austin Schuh70cc9552019-01-21 19:46:48 -080073
74 // When using the subset preconditioner, all row blocks starting
75 // from this row block are used to construct the preconditioner.
76 //
77 // i.e., the Jacobian matrix A is horizontally partitioned as
78 //
79 // A = [P]
80 // [Q]
81 //
82 // where P has subset_preconditioner_start_row_block row blocks,
83 // and the preconditioner is the inverse of the matrix Q'Q.
84 int subset_preconditioner_start_row_block = -1;
85
Austin Schuh70cc9552019-01-21 19:46:48 -080086 // If possible, how many threads the preconditioner can use.
87 int num_threads = 1;
88
89 // Hints about the order in which the parameter blocks should be
90 // eliminated by the linear solver.
91 //
92 // For example if elimination_groups is a vector of size k, then
93 // the linear solver is informed that it should eliminate the
94 // parameter blocks 0 ... elimination_groups[0] - 1 first, and
95 // then elimination_groups[0] ... elimination_groups[1] - 1 and so
96 // on. Within each elimination group, the linear solver is free to
97 // choose how the parameter blocks are ordered. Different linear
98 // solvers have differing requirements on elimination_groups.
99 //
100 // The most common use is for Schur type solvers, where there
101 // should be at least two elimination groups and the first
102 // elimination group must form an independent set in the normal
103 // equations. The first elimination group corresponds to the
104 // num_eliminate_blocks in the Schur type solvers.
105 std::vector<int> elimination_groups;
106
107 // If the block sizes in a BlockSparseMatrix are fixed, then in
108 // some cases the Schur complement based solvers can detect and
109 // specialize on them.
110 //
111 // It is expected that these parameters are set programmatically
112 // rather than manually.
113 //
114 // Please see schur_complement_solver.h and schur_eliminator.h for
115 // more details.
116 int row_block_size = Eigen::Dynamic;
117 int e_block_size = Eigen::Dynamic;
118 int f_block_size = Eigen::Dynamic;
119
120 ContextImpl* context = nullptr;
121 };
122
123 // If the optimization problem is such that there are no remaining
124 // e-blocks, ITERATIVE_SCHUR with a Schur type preconditioner cannot
125 // be used. This function returns JACOBI if a preconditioner for
126 // ITERATIVE_SCHUR is used. The input preconditioner_type is
127 // returned otherwise.
128 static PreconditionerType PreconditionerForZeroEBlocks(
129 PreconditionerType preconditioner_type);
130
Austin Schuh3de38b02024-06-25 18:25:10 -0700131 ~Preconditioner() override;
Austin Schuh70cc9552019-01-21 19:46:48 -0800132
133 // Update the numerical value of the preconditioner for the linear
134 // system:
135 //
136 // | A | x = |b|
137 // |diag(D)| |0|
138 //
139 // for some vector b. It is important that the matrix A have the
140 // same block structure as the one used to construct this object.
141 //
Austin Schuh3de38b02024-06-25 18:25:10 -0700142 // D can be nullptr, in which case its interpreted as a diagonal matrix
Austin Schuh70cc9552019-01-21 19:46:48 -0800143 // of size zero.
144 virtual bool Update(const LinearOperator& A, const double* D) = 0;
145
146 // LinearOperator interface. Since the operator is symmetric,
Austin Schuh3de38b02024-06-25 18:25:10 -0700147 // LeftMultiplyAndAccumulate and num_cols are just calls to
148 // RightMultiplyAndAccumulate and num_rows respectively. Update() must be
149 // called before RightMultiplyAndAccumulate can be called.
150 void RightMultiplyAndAccumulate(const double* x,
151 double* y) const override = 0;
152 void LeftMultiplyAndAccumulate(const double* x, double* y) const override {
153 return RightMultiplyAndAccumulate(x, y);
Austin Schuh70cc9552019-01-21 19:46:48 -0800154 }
155
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800156 int num_rows() const override = 0;
157 int num_cols() const override { return num_rows(); }
Austin Schuh70cc9552019-01-21 19:46:48 -0800158};
159
Austin Schuh3de38b02024-06-25 18:25:10 -0700160class CERES_NO_EXPORT IdentityPreconditioner : public Preconditioner {
161 public:
162 IdentityPreconditioner(int num_rows) : num_rows_(num_rows) {}
163
164 bool Update(const LinearOperator& /*A*/, const double* /*D*/) final {
165 return true;
166 }
167
168 void RightMultiplyAndAccumulate(const double* x, double* y) const final {
169 VectorRef(y, num_rows_) += ConstVectorRef(x, num_rows_);
170 }
171
172 int num_rows() const final { return num_rows_; }
173
174 private:
175 int num_rows_ = -1;
176};
177
Austin Schuh70cc9552019-01-21 19:46:48 -0800178// This templated subclass of Preconditioner serves as a base class for
179// other preconditioners that depend on the particular matrix layout of
180// the underlying linear operator.
181template <typename MatrixType>
Austin Schuh3de38b02024-06-25 18:25:10 -0700182class CERES_NO_EXPORT TypedPreconditioner : public Preconditioner {
Austin Schuh70cc9552019-01-21 19:46:48 -0800183 public:
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800184 bool Update(const LinearOperator& A, const double* D) final {
Austin Schuh70cc9552019-01-21 19:46:48 -0800185 return UpdateImpl(*down_cast<const MatrixType*>(&A), D);
186 }
187
188 private:
189 virtual bool UpdateImpl(const MatrixType& A, const double* D) = 0;
190};
191
192// Preconditioners that depend on access to the low level structure
193// of a SparseMatrix.
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800194// clang-format off
Austin Schuh3de38b02024-06-25 18:25:10 -0700195using SparseMatrixPreconditioner = TypedPreconditioner<SparseMatrix>;
196using BlockSparseMatrixPreconditioner = TypedPreconditioner<BlockSparseMatrix>;
197using CompressedRowSparseMatrixPreconditioner = TypedPreconditioner<CompressedRowSparseMatrix>;
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800198// clang-format on
Austin Schuh70cc9552019-01-21 19:46:48 -0800199
200// Wrap a SparseMatrix object as a preconditioner.
Austin Schuh3de38b02024-06-25 18:25:10 -0700201class CERES_NO_EXPORT SparseMatrixPreconditionerWrapper final
202 : public SparseMatrixPreconditioner {
Austin Schuh70cc9552019-01-21 19:46:48 -0800203 public:
204 // Wrapper does NOT take ownership of the matrix pointer.
Austin Schuh3de38b02024-06-25 18:25:10 -0700205 explicit SparseMatrixPreconditionerWrapper(
206 const SparseMatrix* matrix, const Preconditioner::Options& options);
207 ~SparseMatrixPreconditionerWrapper() override;
Austin Schuh70cc9552019-01-21 19:46:48 -0800208
209 // Preconditioner interface
Austin Schuh3de38b02024-06-25 18:25:10 -0700210 void RightMultiplyAndAccumulate(const double* x, double* y) const override;
211 int num_rows() const override;
Austin Schuh70cc9552019-01-21 19:46:48 -0800212
213 private:
Austin Schuh3de38b02024-06-25 18:25:10 -0700214 bool UpdateImpl(const SparseMatrix& A, const double* D) override;
Austin Schuh70cc9552019-01-21 19:46:48 -0800215 const SparseMatrix* matrix_;
Austin Schuh3de38b02024-06-25 18:25:10 -0700216 const Preconditioner::Options options_;
Austin Schuh70cc9552019-01-21 19:46:48 -0800217};
218
Austin Schuh3de38b02024-06-25 18:25:10 -0700219} // namespace ceres::internal
220
221#include "ceres/internal/reenable_warnings.h"
Austin Schuh70cc9552019-01-21 19:46:48 -0800222
223#endif // CERES_INTERNAL_PRECONDITIONER_H_