blob: 6cb8e9e4ed1e7c72501134249fa70685f2fbfba8 [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/program.h"
32
33#include <cmath>
34#include <limits>
35#include <memory>
36#include <vector>
37
38#include "ceres/internal/integer_sequence_algorithm.h"
39#include "ceres/problem_impl.h"
40#include "ceres/residual_block.h"
41#include "ceres/sized_cost_function.h"
42#include "ceres/triplet_sparse_matrix.h"
43#include "gtest/gtest.h"
44
45namespace ceres {
46namespace internal {
47
48using std::string;
49using std::vector;
50
51// A cost function that simply returns its argument.
52class UnaryIdentityCostFunction : public SizedCostFunction<1, 1> {
53 public:
54 virtual bool Evaluate(double const* const* parameters,
55 double* residuals,
56 double** jacobians) const {
57 residuals[0] = parameters[0][0];
58 if (jacobians != nullptr && jacobians[0] != nullptr) {
59 jacobians[0][0] = 1.0;
60 }
61 return true;
62 }
63};
64
65// Templated base class for the CostFunction signatures.
66template <int kNumResiduals, int... Ns>
67class MockCostFunctionBase : public SizedCostFunction<kNumResiduals, Ns...> {
68 public:
69 virtual bool Evaluate(double const* const* parameters,
70 double* residuals,
71 double** jacobians) const {
72 const int kNumParameters = Sum<integer_sequence<int, Ns...>>::Value;
73
74 for (int i = 0; i < kNumResiduals; ++i) {
75 residuals[i] = kNumResiduals + kNumParameters;
76 }
77 return true;
78 }
79};
80
81class UnaryCostFunction : public MockCostFunctionBase<2, 1> {};
82class BinaryCostFunction : public MockCostFunctionBase<2, 1, 1> {};
83class TernaryCostFunction : public MockCostFunctionBase<2, 1, 1, 1> {};
84
85TEST(Program, RemoveFixedBlocksNothingConstant) {
86 ProblemImpl problem;
87 double x;
88 double y;
89 double z;
90
91 problem.AddParameterBlock(&x, 1);
92 problem.AddParameterBlock(&y, 1);
93 problem.AddParameterBlock(&z, 1);
94 problem.AddResidualBlock(new UnaryCostFunction(), nullptr, &x);
95 problem.AddResidualBlock(new BinaryCostFunction(), nullptr, &x, &y);
96 problem.AddResidualBlock(new TernaryCostFunction(), nullptr, &x, &y, &z);
97
98 vector<double*> removed_parameter_blocks;
99 double fixed_cost = 0.0;
100 string message;
101 std::unique_ptr<Program> reduced_program(problem.program().CreateReducedProgram(
102 &removed_parameter_blocks, &fixed_cost, &message));
103
104 EXPECT_EQ(reduced_program->NumParameterBlocks(), 3);
105 EXPECT_EQ(reduced_program->NumResidualBlocks(), 3);
106 EXPECT_EQ(removed_parameter_blocks.size(), 0);
107 EXPECT_EQ(fixed_cost, 0.0);
108}
109
110TEST(Program, RemoveFixedBlocksAllParameterBlocksConstant) {
111 ProblemImpl problem;
112 double x = 1.0;
113
114 problem.AddParameterBlock(&x, 1);
115 problem.AddResidualBlock(new UnaryCostFunction(), nullptr, &x);
116 problem.SetParameterBlockConstant(&x);
117
118 vector<double*> removed_parameter_blocks;
119 double fixed_cost = 0.0;
120 string message;
121 std::unique_ptr<Program> reduced_program(
122 problem.program().CreateReducedProgram(
123 &removed_parameter_blocks, &fixed_cost, &message));
124
125 EXPECT_EQ(reduced_program->NumParameterBlocks(), 0);
126 EXPECT_EQ(reduced_program->NumResidualBlocks(), 0);
127 EXPECT_EQ(removed_parameter_blocks.size(), 1);
128 EXPECT_EQ(removed_parameter_blocks[0], &x);
129 EXPECT_EQ(fixed_cost, 9.0);
130}
131
132
133TEST(Program, RemoveFixedBlocksNoResidualBlocks) {
134 ProblemImpl problem;
135 double x;
136 double y;
137 double z;
138
139 problem.AddParameterBlock(&x, 1);
140 problem.AddParameterBlock(&y, 1);
141 problem.AddParameterBlock(&z, 1);
142
143 vector<double*> removed_parameter_blocks;
144 double fixed_cost = 0.0;
145 string message;
146 std::unique_ptr<Program> reduced_program(
147 problem.program().CreateReducedProgram(
148 &removed_parameter_blocks, &fixed_cost, &message));
149 EXPECT_EQ(reduced_program->NumParameterBlocks(), 0);
150 EXPECT_EQ(reduced_program->NumResidualBlocks(), 0);
151 EXPECT_EQ(removed_parameter_blocks.size(), 3);
152 EXPECT_EQ(fixed_cost, 0.0);
153}
154
155TEST(Program, RemoveFixedBlocksOneParameterBlockConstant) {
156 ProblemImpl problem;
157 double x;
158 double y;
159 double z;
160
161 problem.AddParameterBlock(&x, 1);
162 problem.AddParameterBlock(&y, 1);
163 problem.AddParameterBlock(&z, 1);
164
165 problem.AddResidualBlock(new UnaryCostFunction(), nullptr, &x);
166 problem.AddResidualBlock(new BinaryCostFunction(), nullptr, &x, &y);
167 problem.SetParameterBlockConstant(&x);
168
169 vector<double*> removed_parameter_blocks;
170 double fixed_cost = 0.0;
171 string message;
172 std::unique_ptr<Program> reduced_program(
173 problem.program().CreateReducedProgram(
174 &removed_parameter_blocks, &fixed_cost, &message));
175 EXPECT_EQ(reduced_program->NumParameterBlocks(), 1);
176 EXPECT_EQ(reduced_program->NumResidualBlocks(), 1);
177}
178
179TEST(Program, RemoveFixedBlocksNumEliminateBlocks) {
180 ProblemImpl problem;
181 double x;
182 double y;
183 double z;
184
185 problem.AddParameterBlock(&x, 1);
186 problem.AddParameterBlock(&y, 1);
187 problem.AddParameterBlock(&z, 1);
188 problem.AddResidualBlock(new UnaryCostFunction(), nullptr, &x);
189 problem.AddResidualBlock(new TernaryCostFunction(), nullptr, &x, &y, &z);
190 problem.AddResidualBlock(new BinaryCostFunction(), nullptr, &x, &y);
191 problem.SetParameterBlockConstant(&x);
192
193 vector<double*> removed_parameter_blocks;
194 double fixed_cost = 0.0;
195 string message;
196 std::unique_ptr<Program> reduced_program(
197 problem.program().CreateReducedProgram(
198 &removed_parameter_blocks, &fixed_cost, &message));
199 EXPECT_EQ(reduced_program->NumParameterBlocks(), 2);
200 EXPECT_EQ(reduced_program->NumResidualBlocks(), 2);
201}
202
203TEST(Program, RemoveFixedBlocksFixedCost) {
204 ProblemImpl problem;
205 double x = 1.23;
206 double y = 4.56;
207 double z = 7.89;
208
209 problem.AddParameterBlock(&x, 1);
210 problem.AddParameterBlock(&y, 1);
211 problem.AddParameterBlock(&z, 1);
212 problem.AddResidualBlock(new UnaryIdentityCostFunction(), nullptr, &x);
213 problem.AddResidualBlock(new TernaryCostFunction(), nullptr, &x, &y, &z);
214 problem.AddResidualBlock(new BinaryCostFunction(), nullptr, &x, &y);
215 problem.SetParameterBlockConstant(&x);
216
217 ResidualBlock *expected_removed_block =
218 problem.program().residual_blocks()[0];
219 std::unique_ptr<double[]> scratch(
220 new double[expected_removed_block->NumScratchDoublesForEvaluate()]);
221 double expected_fixed_cost;
222 expected_removed_block->Evaluate(true,
223 &expected_fixed_cost,
224 nullptr,
225 nullptr,
226 scratch.get());
227
228
229 vector<double*> removed_parameter_blocks;
230 double fixed_cost = 0.0;
231 string message;
232 std::unique_ptr<Program> reduced_program(
233 problem.program().CreateReducedProgram(
234 &removed_parameter_blocks, &fixed_cost, &message));
235
236 EXPECT_EQ(reduced_program->NumParameterBlocks(), 2);
237 EXPECT_EQ(reduced_program->NumResidualBlocks(), 2);
238 EXPECT_DOUBLE_EQ(fixed_cost, expected_fixed_cost);
239}
240
241TEST(Program, CreateJacobianBlockSparsityTranspose) {
242 ProblemImpl problem;
243 double x[2];
244 double y[3];
245 double z;
246
247 problem.AddParameterBlock(x, 2);
248 problem.AddParameterBlock(y, 3);
249 problem.AddParameterBlock(&z, 1);
250
251 problem.AddResidualBlock(new MockCostFunctionBase<2, 2>(), nullptr, x);
252 problem.AddResidualBlock(new MockCostFunctionBase<3, 1, 2>(), nullptr, &z, x);
253 problem.AddResidualBlock(new MockCostFunctionBase<4, 1, 3>(), nullptr, &z, y);
254 problem.AddResidualBlock(new MockCostFunctionBase<5, 1, 3>(), nullptr, &z, y);
255 problem.AddResidualBlock(new MockCostFunctionBase<1, 2, 1>(), nullptr, x, &z);
256 problem.AddResidualBlock(new MockCostFunctionBase<2, 1, 3>(), nullptr, &z, y);
257 problem.AddResidualBlock(new MockCostFunctionBase<2, 2, 1>(), nullptr, x, &z);
258 problem.AddResidualBlock(new MockCostFunctionBase<1, 3>(), nullptr, y);
259
260 TripletSparseMatrix expected_block_sparse_jacobian(3, 8, 14);
261 {
262 int* rows = expected_block_sparse_jacobian.mutable_rows();
263 int* cols = expected_block_sparse_jacobian.mutable_cols();
264 double* values = expected_block_sparse_jacobian.mutable_values();
265 rows[0] = 0;
266 cols[0] = 0;
267
268 rows[1] = 2;
269 cols[1] = 1;
270 rows[2] = 0;
271 cols[2] = 1;
272
273 rows[3] = 2;
274 cols[3] = 2;
275 rows[4] = 1;
276 cols[4] = 2;
277
278 rows[5] = 2;
279 cols[5] = 3;
280 rows[6] = 1;
281 cols[6] = 3;
282
283 rows[7] = 0;
284 cols[7] = 4;
285 rows[8] = 2;
286 cols[8] = 4;
287
288 rows[9] = 2;
289 cols[9] = 5;
290 rows[10] = 1;
291 cols[10] = 5;
292
293 rows[11] = 0;
294 cols[11] = 6;
295 rows[12] = 2;
296 cols[12] = 6;
297
298 rows[13] = 1;
299 cols[13] = 7;
300 std::fill(values, values + 14, 1.0);
301 expected_block_sparse_jacobian.set_num_nonzeros(14);
302 }
303
304 Program* program = problem.mutable_program();
305 program->SetParameterOffsetsAndIndex();
306
307 std::unique_ptr<TripletSparseMatrix> actual_block_sparse_jacobian(
308 program->CreateJacobianBlockSparsityTranspose());
309
310 Matrix expected_dense_jacobian;
311 expected_block_sparse_jacobian.ToDenseMatrix(&expected_dense_jacobian);
312
313 Matrix actual_dense_jacobian;
314 actual_block_sparse_jacobian->ToDenseMatrix(&actual_dense_jacobian);
315 EXPECT_EQ((expected_dense_jacobian - actual_dense_jacobian).norm(), 0.0);
316}
317
318template <int kNumResiduals, int kNumParameterBlocks>
319class NumParameterBlocksCostFunction : public CostFunction {
320 public:
321 NumParameterBlocksCostFunction() {
322 set_num_residuals(kNumResiduals);
323 for (int i = 0; i < kNumParameterBlocks; ++i) {
324 mutable_parameter_block_sizes()->push_back(1);
325 }
326 }
327
328 virtual ~NumParameterBlocksCostFunction() {
329 }
330
331 virtual bool Evaluate(double const* const* parameters,
332 double* residuals,
333 double** jacobians) const {
334 return true;
335 }
336};
337
338TEST(Program, ReallocationInCreateJacobianBlockSparsityTranspose) {
339 // CreateJacobianBlockSparsityTranspose starts with a conservative
340 // estimate of the size of the sparsity pattern. This test ensures
341 // that when those estimates are violated, the reallocation/resizing
342 // logic works correctly.
343
344 ProblemImpl problem;
345 double x[20];
346
347 vector<double*> parameter_blocks;
348 for (int i = 0; i < 20; ++i) {
349 problem.AddParameterBlock(x + i, 1);
350 parameter_blocks.push_back(x + i);
351 }
352
353 problem.AddResidualBlock(new NumParameterBlocksCostFunction<1, 20>(),
354 nullptr,
355 parameter_blocks.data(),
356 static_cast<int>(parameter_blocks.size()));
357
358 TripletSparseMatrix expected_block_sparse_jacobian(20, 1, 20);
359 {
360 int* rows = expected_block_sparse_jacobian.mutable_rows();
361 int* cols = expected_block_sparse_jacobian.mutable_cols();
362 for (int i = 0; i < 20; ++i) {
363 rows[i] = i;
364 cols[i] = 0;
365 }
366
367 double* values = expected_block_sparse_jacobian.mutable_values();
368 std::fill(values, values + 20, 1.0);
369 expected_block_sparse_jacobian.set_num_nonzeros(20);
370 }
371
372 Program* program = problem.mutable_program();
373 program->SetParameterOffsetsAndIndex();
374
375 std::unique_ptr<TripletSparseMatrix> actual_block_sparse_jacobian(
376 program->CreateJacobianBlockSparsityTranspose());
377
378 Matrix expected_dense_jacobian;
379 expected_block_sparse_jacobian.ToDenseMatrix(&expected_dense_jacobian);
380
381 Matrix actual_dense_jacobian;
382 actual_block_sparse_jacobian->ToDenseMatrix(&actual_dense_jacobian);
383 EXPECT_EQ((expected_dense_jacobian - actual_dense_jacobian).norm(), 0.0);
384}
385
386TEST(Program, ProblemHasNanParameterBlocks) {
387 ProblemImpl problem;
388 double x[2];
389 x[0] = 1.0;
390 x[1] = std::numeric_limits<double>::quiet_NaN();
391 problem.AddResidualBlock(new MockCostFunctionBase<1, 2>(), nullptr, x);
392 string error;
393 EXPECT_FALSE(problem.program().ParameterBlocksAreFinite(&error));
394 EXPECT_NE(error.find("has at least one invalid value"),
395 string::npos) << error;
396}
397
398TEST(Program, InfeasibleParameterBlock) {
399 ProblemImpl problem;
400 double x[] = {0.0, 0.0};
401 problem.AddResidualBlock(new MockCostFunctionBase<1, 2>(), nullptr, x);
402 problem.SetParameterLowerBound(x, 0, 2.0);
403 problem.SetParameterUpperBound(x, 0, 1.0);
404 string error;
405 EXPECT_FALSE(problem.program().IsFeasible(&error));
406 EXPECT_NE(error.find("infeasible bound"), string::npos) << error;
407}
408
409TEST(Program, InfeasibleConstantParameterBlock) {
410 ProblemImpl problem;
411 double x[] = {0.0, 0.0};
412 problem.AddResidualBlock(new MockCostFunctionBase<1, 2>(), nullptr, x);
413 problem.SetParameterLowerBound(x, 0, 1.0);
414 problem.SetParameterUpperBound(x, 0, 2.0);
415 problem.SetParameterBlockConstant(x);
416 string error;
417 EXPECT_FALSE(problem.program().IsFeasible(&error));
418 EXPECT_NE(error.find("infeasible value"), string::npos) << error;
419}
420
421} // namespace internal
422} // namespace ceres