blob: f2196d4ef9c11bac128f21f80652384da73ee5ed [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#include "ceres/implicit_schur_complement.h"
32
33#include "Eigen/Dense"
34#include "ceres/block_sparse_matrix.h"
35#include "ceres/block_structure.h"
36#include "ceres/internal/eigen.h"
37#include "ceres/linear_solver.h"
38#include "ceres/types.h"
39#include "glog/logging.h"
40
41namespace ceres {
42namespace internal {
43
44ImplicitSchurComplement::ImplicitSchurComplement(
45 const LinearSolver::Options& options)
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080046 : options_(options), D_(NULL), b_(NULL) {}
Austin Schuh70cc9552019-01-21 19:46:48 -080047
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080048ImplicitSchurComplement::~ImplicitSchurComplement() {}
Austin Schuh70cc9552019-01-21 19:46:48 -080049
50void ImplicitSchurComplement::Init(const BlockSparseMatrix& A,
51 const double* D,
52 const double* b) {
53 // Since initialization is reasonably heavy, perhaps we can save on
54 // constructing a new object everytime.
55 if (A_ == NULL) {
56 A_.reset(PartitionedMatrixViewBase::Create(options_, A));
57 }
58
59 D_ = D;
60 b_ = b;
61
62 // Initialize temporary storage and compute the block diagonals of
63 // E'E and F'E.
64 if (block_diagonal_EtE_inverse_ == NULL) {
65 block_diagonal_EtE_inverse_.reset(A_->CreateBlockDiagonalEtE());
66 if (options_.preconditioner_type == JACOBI) {
67 block_diagonal_FtF_inverse_.reset(A_->CreateBlockDiagonalFtF());
68 }
69 rhs_.resize(A_->num_cols_f());
70 rhs_.setZero();
71 tmp_rows_.resize(A_->num_rows());
72 tmp_e_cols_.resize(A_->num_cols_e());
73 tmp_e_cols_2_.resize(A_->num_cols_e());
74 tmp_f_cols_.resize(A_->num_cols_f());
75 } else {
76 A_->UpdateBlockDiagonalEtE(block_diagonal_EtE_inverse_.get());
77 if (options_.preconditioner_type == JACOBI) {
78 A_->UpdateBlockDiagonalFtF(block_diagonal_FtF_inverse_.get());
79 }
80 }
81
82 // The block diagonals of the augmented linear system contain
83 // contributions from the diagonal D if it is non-null. Add that to
84 // the block diagonals and invert them.
85 AddDiagonalAndInvert(D_, block_diagonal_EtE_inverse_.get());
86 if (options_.preconditioner_type == JACOBI) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080087 AddDiagonalAndInvert((D_ == NULL) ? NULL : D_ + A_->num_cols_e(),
Austin Schuh70cc9552019-01-21 19:46:48 -080088 block_diagonal_FtF_inverse_.get());
89 }
90
91 // Compute the RHS of the Schur complement system.
92 UpdateRhs();
93}
94
95// Evaluate the product
96//
97// Sx = [F'F - F'E (E'E)^-1 E'F]x
98//
99// By breaking it down into individual matrix vector products
100// involving the matrices E and F. This is implemented using a
101// PartitionedMatrixView of the input matrix A.
102void ImplicitSchurComplement::RightMultiply(const double* x, double* y) const {
103 // y1 = F x
104 tmp_rows_.setZero();
105 A_->RightMultiplyF(x, tmp_rows_.data());
106
107 // y2 = E' y1
108 tmp_e_cols_.setZero();
109 A_->LeftMultiplyE(tmp_rows_.data(), tmp_e_cols_.data());
110
111 // y3 = -(E'E)^-1 y2
112 tmp_e_cols_2_.setZero();
113 block_diagonal_EtE_inverse_->RightMultiply(tmp_e_cols_.data(),
114 tmp_e_cols_2_.data());
115 tmp_e_cols_2_ *= -1.0;
116
117 // y1 = y1 + E y3
118 A_->RightMultiplyE(tmp_e_cols_2_.data(), tmp_rows_.data());
119
120 // y5 = D * x
121 if (D_ != NULL) {
122 ConstVectorRef Dref(D_ + A_->num_cols_e(), num_cols());
123 VectorRef(y, num_cols()) =
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800124 (Dref.array().square() * ConstVectorRef(x, num_cols()).array())
125 .matrix();
Austin Schuh70cc9552019-01-21 19:46:48 -0800126 } else {
127 VectorRef(y, num_cols()).setZero();
128 }
129
130 // y = y5 + F' y1
131 A_->LeftMultiplyF(tmp_rows_.data(), y);
132}
133
134// Given a block diagonal matrix and an optional array of diagonal
135// entries D, add them to the diagonal of the matrix and compute the
136// inverse of each diagonal block.
137void ImplicitSchurComplement::AddDiagonalAndInvert(
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800138 const double* D, BlockSparseMatrix* block_diagonal) {
Austin Schuh70cc9552019-01-21 19:46:48 -0800139 const CompressedRowBlockStructure* block_diagonal_structure =
140 block_diagonal->block_structure();
141 for (int r = 0; r < block_diagonal_structure->rows.size(); ++r) {
142 const int row_block_pos = block_diagonal_structure->rows[r].block.position;
143 const int row_block_size = block_diagonal_structure->rows[r].block.size;
144 const Cell& cell = block_diagonal_structure->rows[r].cells[0];
145 MatrixRef m(block_diagonal->mutable_values() + cell.position,
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800146 row_block_size,
147 row_block_size);
Austin Schuh70cc9552019-01-21 19:46:48 -0800148
149 if (D != NULL) {
150 ConstVectorRef d(D + row_block_pos, row_block_size);
151 m += d.array().square().matrix().asDiagonal();
152 }
153
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800154 m = m.selfadjointView<Eigen::Upper>().llt().solve(
155 Matrix::Identity(row_block_size, row_block_size));
Austin Schuh70cc9552019-01-21 19:46:48 -0800156 }
157}
158
159// Similar to RightMultiply, use the block structure of the matrix A
160// to compute y = (E'E)^-1 (E'b - E'F x).
161void ImplicitSchurComplement::BackSubstitute(const double* x, double* y) {
162 const int num_cols_e = A_->num_cols_e();
163 const int num_cols_f = A_->num_cols_f();
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800164 const int num_cols = A_->num_cols();
Austin Schuh70cc9552019-01-21 19:46:48 -0800165 const int num_rows = A_->num_rows();
166
167 // y1 = F x
168 tmp_rows_.setZero();
169 A_->RightMultiplyF(x, tmp_rows_.data());
170
171 // y2 = b - y1
172 tmp_rows_ = ConstVectorRef(b_, num_rows) - tmp_rows_;
173
174 // y3 = E' y2
175 tmp_e_cols_.setZero();
176 A_->LeftMultiplyE(tmp_rows_.data(), tmp_e_cols_.data());
177
178 // y = (E'E)^-1 y3
179 VectorRef(y, num_cols).setZero();
180 block_diagonal_EtE_inverse_->RightMultiply(tmp_e_cols_.data(), y);
181
182 // The full solution vector y has two blocks. The first block of
183 // variables corresponds to the eliminated variables, which we just
184 // computed via back substitution. The second block of variables
185 // corresponds to the Schur complement system, so we just copy those
186 // values from the solution to the Schur complement.
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800187 VectorRef(y + num_cols_e, num_cols_f) = ConstVectorRef(x, num_cols_f);
Austin Schuh70cc9552019-01-21 19:46:48 -0800188}
189
190// Compute the RHS of the Schur complement system.
191//
192// rhs = F'b - F'E (E'E)^-1 E'b
193//
194// Like BackSubstitute, we use the block structure of A to implement
195// this using a series of matrix vector products.
196void ImplicitSchurComplement::UpdateRhs() {
197 // y1 = E'b
198 tmp_e_cols_.setZero();
199 A_->LeftMultiplyE(b_, tmp_e_cols_.data());
200
201 // y2 = (E'E)^-1 y1
202 Vector y2 = Vector::Zero(A_->num_cols_e());
203 block_diagonal_EtE_inverse_->RightMultiply(tmp_e_cols_.data(), y2.data());
204
205 // y3 = E y2
206 tmp_rows_.setZero();
207 A_->RightMultiplyE(y2.data(), tmp_rows_.data());
208
209 // y3 = b - y3
210 tmp_rows_ = ConstVectorRef(b_, A_->num_rows()) - tmp_rows_;
211
212 // rhs = F' y3
213 rhs_.setZero();
214 A_->LeftMultiplyF(tmp_rows_.data(), rhs_.data());
215}
216
217} // namespace internal
218} // namespace ceres