1
// Ceres Solver - A fast non-linear least squares minimizer
2
// Copyright 2010, 2011, 2012 Google Inc. All rights reserved.
3
// http://code.google.com/p/ceres-solver/
5
// Redistribution and use in source and binary forms, with or without
6
// modification, are permitted provided that the following conditions are met:
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.
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.
29
// Author: sameeragarwal@google.com (Sameer Agarwal)
31
#include "ceres/schur_eliminator.h"
33
#include <glog/logging.h>
34
#include "ceres/file.h"
35
#include "gtest/gtest.h"
36
#include "Eigen/Dense"
37
#include "ceres/block_random_access_dense_matrix.h"
38
#include "ceres/block_sparse_matrix.h"
39
#include "ceres/casts.h"
40
#include "ceres/detect_structure.h"
41
#include "ceres/linear_least_squares_problems.h"
42
#include "ceres/triplet_sparse_matrix.h"
43
#include "ceres/internal/eigen.h"
44
#include "ceres/internal/scoped_ptr.h"
45
#include "ceres/types.h"
47
// TODO(sameeragarwal): Reduce the size of these tests and redo the
48
// parameterization to be more efficient.
50
DECLARE_string(test_srcdir);
55
class SchurEliminatorTest : public ::testing::Test {
57
void SetUpFromId(int id) {
58
scoped_ptr<LinearLeastSquaresProblem>
59
problem(CreateLinearLeastSquaresProblemFromId(id));
60
CHECK_NOTNULL(problem.get());
61
SetupHelper(problem.get());
64
void SetUpFromFilename(const string& filename) {
65
scoped_ptr<LinearLeastSquaresProblem>
66
problem(CreateLinearLeastSquaresProblemFromFile(filename));
67
CHECK_NOTNULL(problem.get());
68
SetupHelper(problem.get());
71
void SetupHelper(LinearLeastSquaresProblem* problem) {
72
A.reset(down_cast<BlockSparseMatrix*>(problem->A.release()));
73
b.reset(problem->b.release());
74
D.reset(problem->D.release());
76
num_eliminate_blocks = problem->num_eliminate_blocks;
77
num_eliminate_cols = 0;
78
const CompressedRowBlockStructure* bs = A->block_structure();
80
for (int i = 0; i < num_eliminate_blocks; ++i) {
81
num_eliminate_cols += bs->cols[i].size;
85
// Compute the golden values for the reduced linear system and the
86
// solution to the linear least squares problem using dense linear
88
void ComputeReferenceSolution(const Vector& D) {
91
VectorRef f(b.get(), J.rows());
93
Matrix H = (D.cwiseProduct(D)).asDiagonal();
94
H.noalias() += J.transpose() * J;
96
const Vector g = J.transpose() * f;
97
const int schur_size = J.cols() - num_eliminate_cols;
99
lhs_expected.resize(schur_size, schur_size);
100
lhs_expected.setZero();
102
rhs_expected.resize(schur_size);
103
rhs_expected.setZero();
105
sol_expected.resize(J.cols());
106
sol_expected.setZero();
108
Matrix P = H.block(0, 0, num_eliminate_cols, num_eliminate_cols);
109
Matrix Q = H.block(0,
113
Matrix R = H.block(num_eliminate_cols,
118
const CompressedRowBlockStructure* bs = A->block_structure();
119
for (int i = 0; i < num_eliminate_blocks; ++i) {
120
const int block_size = bs->cols[i].size;
121
P.block(row, row, block_size, block_size) =
123
.block(row, row, block_size, block_size)
125
.solve(Matrix::Identity(block_size, block_size));
130
.triangularView<Eigen::Upper>() = R - Q.transpose() * P * Q;
132
g.tail(schur_size) - Q.transpose() * P * g.head(num_eliminate_cols);
133
sol_expected = H.ldlt().solve(g);
136
void EliminateSolveAndCompare(const VectorRef& diagonal,
137
bool use_static_structure,
138
const double relative_tolerance) {
139
const CompressedRowBlockStructure* bs = A->block_structure();
140
const int num_col_blocks = bs->cols.size();
141
vector<int> blocks(num_col_blocks - num_eliminate_blocks, 0);
142
for (int i = num_eliminate_blocks; i < num_col_blocks; ++i) {
143
blocks[i - num_eliminate_blocks] = bs->cols[i].size;
146
BlockRandomAccessDenseMatrix lhs(blocks);
148
const int num_cols = A->num_cols();
149
const int schur_size = lhs.num_rows();
151
Vector rhs(schur_size);
153
LinearSolver::Options options;
154
options.num_eliminate_blocks = num_eliminate_blocks;
155
if (use_static_structure) {
157
options.num_eliminate_blocks,
158
&options.row_block_size,
159
&options.e_block_size,
160
&options.f_block_size);
163
scoped_ptr<SchurEliminatorBase> eliminator;
164
eliminator.reset(SchurEliminatorBase::Create(options));
165
eliminator->Init(num_eliminate_blocks, A->block_structure());
166
eliminator->Eliminate(A.get(), b.get(), diagonal.data(), &lhs, rhs.data());
168
MatrixRef lhs_ref(lhs.mutable_values(), lhs.num_rows(), lhs.num_cols());
171
.selfadjointView<Eigen::Upper>()
175
// Solution to the linear least squares problem.
176
Vector sol(num_cols);
178
sol.tail(schur_size) = reduced_sol;
179
eliminator->BackSubstitute(A.get(),
185
Matrix delta = (lhs_ref - lhs_expected).selfadjointView<Eigen::Upper>();
186
double diff = delta.norm();
187
EXPECT_NEAR(diff / lhs_expected.norm(), 0.0, relative_tolerance);
188
EXPECT_NEAR((rhs - rhs_expected).norm() / rhs_expected.norm(), 0.0,
190
EXPECT_NEAR((sol - sol_expected).norm() / sol_expected.norm(), 0.0,
194
scoped_ptr<BlockSparseMatrix> A;
195
scoped_array<double> b;
196
scoped_array<double> D;
197
int num_eliminate_blocks;
198
int num_eliminate_cols;
205
TEST_F(SchurEliminatorTest, ScalarProblem) {
207
Vector zero(A->num_cols());
210
ComputeReferenceSolution(VectorRef(zero.data(), A->num_cols()));
211
EliminateSolveAndCompare(VectorRef(zero.data(), A->num_cols()), true, 1e-14);
212
EliminateSolveAndCompare(VectorRef(zero.data(), A->num_cols()), false, 1e-14);
214
ComputeReferenceSolution(VectorRef(D.get(), A->num_cols()));
215
EliminateSolveAndCompare(VectorRef(D.get(), A->num_cols()), true, 1e-14);
216
EliminateSolveAndCompare(VectorRef(D.get(), A->num_cols()), false, 1e-14);
219
#ifndef CERES_DONT_HAVE_PROTOCOL_BUFFERS
220
TEST_F(SchurEliminatorTest, BlockProblem) {
221
const string input_file =
222
JoinPath(FLAGS_test_srcdir,
223
"problem-6-1384-000.lsqp"); // NOLINT
225
SetUpFromFilename(input_file);
226
ComputeReferenceSolution(VectorRef(D.get(), A->num_cols()));
227
EliminateSolveAndCompare(VectorRef(D.get(), A->num_cols()), true, 1e-10);
228
EliminateSolveAndCompare(VectorRef(D.get(), A->num_cols()), false, 1e-10);
230
#endif // CERES_DONT_HAVE_PROTOCOL_BUFFERS
232
} // namespace internal