blob: 1a5001f9c71754bee572973162a6365cc0106cc2 [file] [log] [blame]
Keir Mierlef7898fb2012-05-05 20:55:08 -07001// Ceres Solver - A fast non-linear least squares minimizer
2// Copyright 2012 Google Inc. All rights reserved.
3// http://code.google.com/p/ceres-solver/
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: keir@google.com (Keir Mierle)
30
Keir Mierlee2a6cdc2012-05-07 06:39:56 -070031#include "ceres/block_jacobi_preconditioner.h"
Keir Mierlef7898fb2012-05-05 20:55:08 -070032
33#include "Eigen/Cholesky"
34#include "ceres/block_sparse_matrix.h"
35#include "ceres/block_structure.h"
36#include "ceres/casts.h"
Keir Mierlee2a6cdc2012-05-07 06:39:56 -070037#include "ceres/integral_types.h"
Keir Mierlef7898fb2012-05-05 20:55:08 -070038#include "ceres/internal/eigen.h"
39
40namespace ceres {
41namespace internal {
42
Keir Mierlee2a6cdc2012-05-07 06:39:56 -070043BlockJacobiPreconditioner::BlockJacobiPreconditioner(
Keir Mierlef7898fb2012-05-05 20:55:08 -070044 const LinearOperator& A)
45 : block_structure_(
Keir Mierlee2a6cdc2012-05-07 06:39:56 -070046 *(down_cast<const BlockSparseMatrix*>(&A)->block_structure())),
47 num_rows_(A.num_rows()) {
Keir Mierlef7898fb2012-05-05 20:55:08 -070048 // Calculate the amount of storage needed.
49 int storage_needed = 0;
50 for (int c = 0; c < block_structure_.cols.size(); ++c) {
51 int size = block_structure_.cols[c].size;
52 storage_needed += size * size;
53 }
54
55 // Size the offsets and storage.
56 blocks_.resize(block_structure_.cols.size());
57 block_storage_.resize(storage_needed);
58
59 // Put pointers to the storage in the offsets.
Keir Mierlee2a6cdc2012-05-07 06:39:56 -070060 double* block_cursor = &block_storage_[0];
Keir Mierlef7898fb2012-05-05 20:55:08 -070061 for (int c = 0; c < block_structure_.cols.size(); ++c) {
62 int size = block_structure_.cols[c].size;
63 blocks_[c] = block_cursor;
64 block_cursor += size * size;
65 }
66}
67
Keir Mierlee2a6cdc2012-05-07 06:39:56 -070068BlockJacobiPreconditioner::~BlockJacobiPreconditioner() {
Keir Mierlef7898fb2012-05-05 20:55:08 -070069}
70
Keir Mierlee2a6cdc2012-05-07 06:39:56 -070071void BlockJacobiPreconditioner::Update(const LinearOperator& matrix, const double* D) {
Keir Mierlef7898fb2012-05-05 20:55:08 -070072 const BlockSparseMatrix& A = *(down_cast<const BlockSparseMatrix*>(&matrix));
73 const CompressedRowBlockStructure* bs = A.block_structure();
74
75 // Compute the diagonal blocks by block inner products.
76 std::fill(block_storage_.begin(), block_storage_.end(), 0.0);
77 for (int r = 0; r < bs->rows.size(); ++r) {
78 const int row_block_size = bs->rows[r].block.size;
79 const vector<Cell>& cells = bs->rows[r].cells;
80 const double* row_values = A.RowBlockValues(r);
81 for (int c = 0; c < cells.size(); ++c) {
82 const int col_block_size = bs->cols[cells[c].block_id].size;
83 ConstMatrixRef m(row_values + cells[c].position,
84 row_block_size,
85 col_block_size);
86
87 MatrixRef(blocks_[cells[c].block_id],
88 col_block_size,
89 col_block_size).noalias() += m.transpose() * m;
Keir Mierlee2a6cdc2012-05-07 06:39:56 -070090
91 // TODO(keir): Figure out when the below expression is actually faster
92 // than doing the full rank update. The issue is that for smaller sizes,
93 // the rankUpdate() function is slower than the full product done above.
94 //
95 // On the typical bundling problems, the above product is ~5% faster.
96 //
97 // MatrixRef(blocks_[cells[c].block_id],
98 // col_block_size,
99 // col_block_size).selfadjointView<Eigen::Upper>().rankUpdate(m);
100 //
Keir Mierlef7898fb2012-05-05 20:55:08 -0700101 }
102 }
103
Keir Mierlee2a6cdc2012-05-07 06:39:56 -0700104 // Add the diagonal and invert each block.
Keir Mierlef7898fb2012-05-05 20:55:08 -0700105 for (int c = 0; c < bs->cols.size(); ++c) {
106 const int size = block_structure_.cols[c].size;
Keir Mierlee2a6cdc2012-05-07 06:39:56 -0700107 const int position = block_structure_.cols[c].position;
Keir Mierlef8bd7fa2012-05-08 05:39:32 -0700108 MatrixRef block(blocks_[c], size, size);
Keir Mierlee2a6cdc2012-05-07 06:39:56 -0700109
Keir Mierlef8bd7fa2012-05-08 05:39:32 -0700110 if (D != NULL) {
111 block.diagonal() += ConstVectorRef(D + position, size).array().square().matrix();
112 }
Keir Mierlee2a6cdc2012-05-07 06:39:56 -0700113
Keir Mierlef8bd7fa2012-05-08 05:39:32 -0700114 block = block.selfadjointView<Eigen::Upper>()
115 .ldlt()
116 .solve(Matrix::Identity(size, size));
Keir Mierlef7898fb2012-05-05 20:55:08 -0700117 }
118}
119
Keir Mierlee2a6cdc2012-05-07 06:39:56 -0700120void BlockJacobiPreconditioner::RightMultiply(const double* x, double* y) const {
Keir Mierlef7898fb2012-05-05 20:55:08 -0700121 for (int c = 0; c < block_structure_.cols.size(); ++c) {
122 const int size = block_structure_.cols[c].size;
123 const int position = block_structure_.cols[c].position;
124 ConstMatrixRef D(blocks_[c], size, size);
125 ConstVectorRef x_block(x + position, size);
126 VectorRef y_block(y + position, size);
127 y_block += D * x_block;
128 }
129}
130
Keir Mierlee2a6cdc2012-05-07 06:39:56 -0700131void BlockJacobiPreconditioner::LeftMultiply(const double* x, double* y) const {
Keir Mierlef7898fb2012-05-05 20:55:08 -0700132 RightMultiply(x, y);
133}
134
135} // namespace internal
136} // namespace ceres