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