blob: b7381a447b5925faf83fc072c8cbb2285de2fd6f [file] [log] [blame]
Keir Mierle8ebb0732012-04-30 23:09:08 -07001// 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/
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//
31// Tests shared across evaluators. The tests try all combinations of linear
32// solver and num_eliminate_blocks (for schur-based solvers).
33
34#include "ceres/evaluator.h"
35
36#include "gtest/gtest.h"
37#include "ceres/casts.h"
38#include "ceres/problem_impl.h"
39#include "ceres/program.h"
40#include "ceres/sparse_matrix.h"
41#include "ceres/internal/scoped_ptr.h"
42#include "ceres/local_parameterization.h"
43#include "ceres/types.h"
44#include "ceres/sized_cost_function.h"
45#include "ceres/internal/eigen.h"
46
47namespace ceres {
48namespace internal {
49
50// TODO(keir): Consider pushing this into a common test utils file.
51template<int kFactor, int kNumResiduals,
52 int N0 = 0, int N1 = 0, int N2 = 0, bool kSucceeds = true>
53class ParameterIgnoringCostFunction
54 : public SizedCostFunction<kNumResiduals, N0, N1, N2> {
55 typedef SizedCostFunction<kNumResiduals, N0, N1, N2> Base;
56 public:
57 virtual bool Evaluate(double const* const* parameters,
58 double* residuals,
59 double** jacobians) const {
60 for (int i = 0; i < Base::num_residuals(); ++i) {
61 residuals[i] = i + 1;
62 }
63 if (jacobians) {
64 for (int k = 0; k < Base::parameter_block_sizes().size(); ++k) {
65 // The jacobians here are full sized, but they are transformed in the
66 // evaluator into the "local" jacobian. In the tests, the "subset
67 // constant" parameterization is used, which should pick out columns
68 // from these jacobians. Put values in the jacobian that make this
69 // obvious; in particular, make the jacobians like this:
70 //
71 // 1 2 3 4 ...
72 // 1 2 3 4 ... .* kFactor
73 // 1 2 3 4 ...
74 //
75 // where the multiplication by kFactor makes it easier to distinguish
76 // between Jacobians of different residuals for the same parameter.
77 if (jacobians[k] != NULL) {
78 MatrixRef jacobian(jacobians[k],
79 Base::num_residuals(),
80 Base::parameter_block_sizes()[k]);
81 for (int j = 0; j < Base::parameter_block_sizes()[k]; ++j) {
82 jacobian.col(j).setConstant(kFactor * (j + 1));
83 }
84 }
85 }
86 }
87 return kSucceeds;
88 }
89};
90
91struct EvaluatorTest
92 : public ::testing::TestWithParam<pair<LinearSolverType, int> > {
93 Evaluator* CreateEvaluator(Program* program) {
94 // This program is straight from the ProblemImpl, and so has no index/offset
95 // yet; compute it here as required by the evalutor implementations.
96 program->SetParameterOffsetsAndIndex();
97
98 VLOG(1) << "Creating evaluator with type: " << GetParam().first
99 << " and num_eliminate_blocks: " << GetParam().second;
100 Evaluator::Options options;
101 options.linear_solver_type = GetParam().first;
102 options.num_eliminate_blocks = GetParam().second;
103 string error;
104 return Evaluator::Create(options, program, &error);
105 }
106};
107
108void SetSparseMatrixConstant(SparseMatrix* sparse_matrix, double value) {
109 VectorRef(sparse_matrix->mutable_values(),
110 sparse_matrix->num_nonzeros()).setConstant(value);
111}
112
113TEST_P(EvaluatorTest, SingleResidualProblem) {
114 ProblemImpl problem;
115
116 // The values are ignored completely by the cost function.
117 double x[2];
118 double y[3];
119 double z[4];
120 double state[9];
121
122 problem.AddResidualBlock(new ParameterIgnoringCostFunction<1, 3, 2, 3, 4>,
123 NULL,
124 x, y, z);
125
126 scoped_ptr<Evaluator> evaluator(CreateEvaluator(problem.mutable_program()));
127 scoped_ptr<SparseMatrix> jacobian(evaluator->CreateJacobian());
128 ASSERT_EQ(3, jacobian->num_rows());
129 ASSERT_EQ(9, jacobian->num_cols());
130
131 // Cost only; no residuals and no jacobian.
132 {
133 double cost = -1;
134 ASSERT_TRUE(evaluator->Evaluate(state, &cost, NULL, NULL));
135 EXPECT_EQ(7.0, cost);
136 }
137
138 // Cost and residuals, no jacobian.
139 {
140 double cost = -1;
141 double residuals[3] = { -2, -2, -2 };
142 ASSERT_TRUE(evaluator->Evaluate(state, &cost, residuals, NULL));
143 EXPECT_EQ(7.0, cost);
144 EXPECT_EQ(1.0, residuals[0]);
145 EXPECT_EQ(2.0, residuals[1]);
146 EXPECT_EQ(3.0, residuals[2]);
147 }
148
149 // Cost, residuals, and jacobian.
150 {
151 double cost = -1;
152 double residuals[3] = { -2, -2, -2 };
153 SetSparseMatrixConstant(jacobian.get(), -1);
154 ASSERT_TRUE(evaluator->Evaluate(state, &cost, residuals, jacobian.get()));
155 EXPECT_EQ(7.0, cost);
156 EXPECT_EQ(1.0, residuals[0]);
157 EXPECT_EQ(2.0, residuals[1]);
158 EXPECT_EQ(3.0, residuals[2]);
159
160 Matrix actual_jacobian;
161 jacobian->ToDenseMatrix(&actual_jacobian);
162
163 Matrix expected_jacobian(3, 9);
164 expected_jacobian
165 // x y z
166 << 1, 2, 1, 2, 3, 1, 2, 3, 4,
167 1, 2, 1, 2, 3, 1, 2, 3, 4,
168 1, 2, 1, 2, 3, 1, 2, 3, 4;
169
170 EXPECT_TRUE((actual_jacobian.array() == expected_jacobian.array()).all())
171 << "Actual:\n" << actual_jacobian
172 << "\nExpected:\n" << expected_jacobian;
173 }
174}
175
176TEST_P(EvaluatorTest, SingleResidualProblemWithPermutedParameters) {
177 ProblemImpl problem;
178
179 // The values are ignored completely by the cost function.
180 double x[2];
181 double y[3];
182 double z[4];
183 double state[9];
184
185 // Add the parameters in explicit order to force the ordering in the program.
186 problem.AddParameterBlock(x, 2);
187 problem.AddParameterBlock(y, 3);
188 problem.AddParameterBlock(z, 4);
189
190 // Then use a cost function which is similar to the others, but swap around
191 // the ordering of the parameters to the cost function. This shouldn't affect
192 // the jacobian evaluation, but requires explicit handling in the evaluators.
193 // At one point the compressed row evaluator had a bug that went undetected
194 // for a long time, since by chance most users added parameters to the problem
195 // in the same order that they occured as parameters to a cost function.
196 problem.AddResidualBlock(new ParameterIgnoringCostFunction<1, 3, 4, 3, 2>,
197 NULL,
198 z, y, x);
199
200 scoped_ptr<Evaluator> evaluator(CreateEvaluator(problem.mutable_program()));
201 scoped_ptr<SparseMatrix> jacobian(evaluator->CreateJacobian());
202 ASSERT_EQ(3, jacobian->num_rows());
203 ASSERT_EQ(9, jacobian->num_cols());
204
205 // Cost only; no residuals and no jacobian.
206 {
207 double cost = -1;
208 ASSERT_TRUE(evaluator->Evaluate(state, &cost, NULL, NULL));
209 EXPECT_EQ(7.0, cost);
210 }
211
212 // Cost and residuals, no jacobian.
213 {
214 double cost = -1;
215 double residuals[3] = { -2, -2, -2 };
216 ASSERT_TRUE(evaluator->Evaluate(state, &cost, residuals, NULL));
217 EXPECT_EQ(7.0, cost);
218 EXPECT_EQ(1.0, residuals[0]);
219 EXPECT_EQ(2.0, residuals[1]);
220 EXPECT_EQ(3.0, residuals[2]);
221 }
222
223 // Cost, residuals, and jacobian.
224 {
225 double cost = -1;
226 double residuals[3] = { -2, -2, -2 };
227 SetSparseMatrixConstant(jacobian.get(), -1);
228 ASSERT_TRUE(evaluator->Evaluate(state, &cost, residuals, jacobian.get()));
229 EXPECT_EQ(7.0, cost);
230 EXPECT_EQ(1.0, residuals[0]);
231 EXPECT_EQ(2.0, residuals[1]);
232 EXPECT_EQ(3.0, residuals[2]);
233
234 Matrix actual_jacobian;
235 jacobian->ToDenseMatrix(&actual_jacobian);
236
237 Matrix expected_jacobian(3, 9);
238 expected_jacobian
239 // x y z
240 << 1, 2, 1, 2, 3, 1, 2, 3, 4,
241 1, 2, 1, 2, 3, 1, 2, 3, 4,
242 1, 2, 1, 2, 3, 1, 2, 3, 4;
243
244 EXPECT_TRUE((actual_jacobian.array() == expected_jacobian.array()).all())
245 << "Actual:\n" << actual_jacobian
246 << "\nExpected:\n" << expected_jacobian;
247 }
248}
249TEST_P(EvaluatorTest, SingleResidualProblemWithNuisanceParameters) {
250 ProblemImpl problem;
251
252 // The values are ignored completely by the cost function.
253 double x[2];
254 double y[3];
255 double z[4];
256 double state[9];
257
258 // These parameters are not used.
259 double w1[2];
260 double w2[1];
261 double w3[1];
262 double w4[3];
263
264 // Add the parameters in a mixed order so the Jacobian is "checkered" with the
265 // values from the other parameters.
266 problem.AddParameterBlock(w1, 2);
267 problem.AddParameterBlock(x, 2);
268 problem.AddParameterBlock(w2, 1);
269 problem.AddParameterBlock(y, 3);
270 problem.AddParameterBlock(w3, 1);
271 problem.AddParameterBlock(z, 4);
272 problem.AddParameterBlock(w4, 3);
273
274 problem.AddResidualBlock(new ParameterIgnoringCostFunction<1, 3, 2, 3, 4>,
275 NULL,
276 x, y, z);
277
278 scoped_ptr<Evaluator> evaluator(CreateEvaluator(problem.mutable_program()));
279 scoped_ptr<SparseMatrix> jacobian(evaluator->CreateJacobian());
280 ASSERT_EQ(3, jacobian->num_rows());
281 ASSERT_EQ(16, jacobian->num_cols());
282
283 // Cost only; no residuals and no jacobian.
284 {
285 double cost = -1;
286 ASSERT_TRUE(evaluator->Evaluate(state, &cost, NULL, NULL));
287 EXPECT_EQ(7.0, cost);
288 }
289
290 // Cost and residuals, no jacobian.
291 {
292 double cost = -1;
293 double residuals[3] = { -2, -2, -2 };
294 ASSERT_TRUE(evaluator->Evaluate(state, &cost, residuals, NULL));
295 EXPECT_EQ(7.0, cost);
296 EXPECT_EQ(1.0, residuals[0]);
297 EXPECT_EQ(2.0, residuals[1]);
298 EXPECT_EQ(3.0, residuals[2]);
299 }
300
301 // Cost, residuals, and jacobian.
302 {
303 double cost = -1;
304 double residuals[3] = { -2, -2, -2 };
305 SetSparseMatrixConstant(jacobian.get(), -1);
306 ASSERT_TRUE(evaluator->Evaluate(state, &cost, residuals, jacobian.get()));
307 EXPECT_EQ(7.0, cost);
308 EXPECT_EQ(1.0, residuals[0]);
309 EXPECT_EQ(2.0, residuals[1]);
310 EXPECT_EQ(3.0, residuals[2]);
311
312 Matrix actual_jacobian;
313 jacobian->ToDenseMatrix(&actual_jacobian);
314
315 Matrix expected_jacobian(3, 16);
316 expected_jacobian
317 // w1 x w2 y w2 z w3
318 << 0, 0, 1, 2, 0, 1, 2, 3, 0, 1, 2, 3, 4, 0, 0, 0,
319 0, 0, 1, 2, 0, 1, 2, 3, 0, 1, 2, 3, 4, 0, 0, 0,
320 0, 0, 1, 2, 0, 1, 2, 3, 0, 1, 2, 3, 4, 0, 0, 0;
321
322 EXPECT_TRUE((actual_jacobian.array() == expected_jacobian.array()).all())
323 << "Actual:\n" << actual_jacobian
324 << "\nExpected:\n" << expected_jacobian;
325 }
326}
327
328TEST_P(EvaluatorTest, MultipleResidualProblem) {
329 ProblemImpl problem;
330
331 // The values are ignored completely by the cost function.
332 double x[2];
333 double y[3];
334 double z[4];
335 double state[9];
336
337 // Add the parameters in explicit order to force the ordering in the program.
338 problem.AddParameterBlock(x, 2);
339 problem.AddParameterBlock(y, 3);
340 problem.AddParameterBlock(z, 4);
341
342 // f(x, y) in R^2
343 problem.AddResidualBlock(new ParameterIgnoringCostFunction<1, 2, 2, 3>,
344 NULL,
345 x, y);
346
347 // g(x, z) in R^3
348 problem.AddResidualBlock(new ParameterIgnoringCostFunction<2, 3, 2, 4>,
349 NULL,
350 x, z);
351
352 // h(y, z) in R^4
353 problem.AddResidualBlock(new ParameterIgnoringCostFunction<3, 4, 3, 4>,
354 NULL,
355 y, z);
356
357 scoped_ptr<Evaluator> evaluator(CreateEvaluator(problem.mutable_program()));
358 scoped_ptr<SparseMatrix> jacobian(evaluator->CreateJacobian());
359 ASSERT_EQ(9, jacobian->num_rows());
360 ASSERT_EQ(9, jacobian->num_cols());
361
362 // f g h
363 double expected_cost = (1 + 4 + 1 + 4 + 9 + 1 + 4 + 9 + 16) / 2.0;
364
365
366 // Cost only; no residuals and no jacobian.
367 {
368 double cost = -1;
369 ASSERT_TRUE(evaluator->Evaluate(state, &cost, NULL, NULL));
370 EXPECT_EQ(expected_cost, cost);
371 }
372
373 // Cost and residuals, no jacobian.
374 {
375 double cost = -1;
376 double residuals[9] = { -2, -2, -2, -2, -2, -2, -2, -2, -2 };
377 ASSERT_TRUE(evaluator->Evaluate(state, &cost, residuals, NULL));
378 EXPECT_EQ(expected_cost, cost);
379 EXPECT_EQ(1.0, residuals[0]);
380 EXPECT_EQ(2.0, residuals[1]);
381 EXPECT_EQ(1.0, residuals[2]);
382 EXPECT_EQ(2.0, residuals[3]);
383 EXPECT_EQ(3.0, residuals[4]);
384 EXPECT_EQ(1.0, residuals[5]);
385 EXPECT_EQ(2.0, residuals[6]);
386 EXPECT_EQ(3.0, residuals[7]);
387 EXPECT_EQ(4.0, residuals[8]);
388 }
389
390 // Cost, residuals, and jacobian.
391 {
392 double cost = -1;
393 double residuals[9] = { -2, -2, -2, -2, -2, -2, -2, -2, -2 };
394 SetSparseMatrixConstant(jacobian.get(), -1);
395 ASSERT_TRUE(evaluator->Evaluate(state, &cost, residuals, jacobian.get()));
396 EXPECT_EQ(expected_cost, cost);
397 EXPECT_EQ(1.0, residuals[0]);
398 EXPECT_EQ(2.0, residuals[1]);
399 EXPECT_EQ(1.0, residuals[2]);
400 EXPECT_EQ(2.0, residuals[3]);
401 EXPECT_EQ(3.0, residuals[4]);
402 EXPECT_EQ(1.0, residuals[5]);
403 EXPECT_EQ(2.0, residuals[6]);
404 EXPECT_EQ(3.0, residuals[7]);
405 EXPECT_EQ(4.0, residuals[8]);
406
407 Matrix actual_jacobian;
408 jacobian->ToDenseMatrix(&actual_jacobian);
409
410 Matrix expected_jacobian(9, 9);
411 expected_jacobian <<
412 // x y z
413 /* f(x, y) */ 1, 2, 1, 2, 3, 0, 0, 0, 0,
414 1, 2, 1, 2, 3, 0, 0, 0, 0,
415
416 /* g(x, z) */ 2, 4, 0, 0, 0, 2, 4, 6, 8,
417 2, 4, 0, 0, 0, 2, 4, 6, 8,
418 2, 4, 0, 0, 0, 2, 4, 6, 8,
419
420 /* h(y, z) */ 0, 0, 3, 6, 9, 3, 6, 9, 12,
421 0, 0, 3, 6, 9, 3, 6, 9, 12,
422 0, 0, 3, 6, 9, 3, 6, 9, 12,
423 0, 0, 3, 6, 9, 3, 6, 9, 12;
424
425 EXPECT_TRUE((actual_jacobian.array() == expected_jacobian.array()).all())
426 << "Actual:\n" << actual_jacobian
427 << "\nExpected:\n" << expected_jacobian;
428 }
429}
430
431TEST_P(EvaluatorTest, MultipleResidualsWithLocalParameterizations) {
432 ProblemImpl problem;
433
434 // The values are ignored completely by the cost function.
435 double x[2];
436 double y[3];
437 double z[4];
438 double state[9];
439
440 // Add the parameters in explicit order to force the ordering in the program.
441 problem.AddParameterBlock(x, 2);
442
443 // Fix y's first dimension.
444 vector<int> y_fixed;
445 y_fixed.push_back(0);
446 problem.AddParameterBlock(y, 3, new SubsetParameterization(3, y_fixed));
447
448 // Fix z's second dimension.
449 vector<int> z_fixed;
450 z_fixed.push_back(1);
451 problem.AddParameterBlock(z, 4, new SubsetParameterization(4, z_fixed));
452
453 // f(x, y) in R^2
454 problem.AddResidualBlock(new ParameterIgnoringCostFunction<1, 2, 2, 3>,
455 NULL,
456 x, y);
457
458 // g(x, z) in R^3
459 problem.AddResidualBlock(new ParameterIgnoringCostFunction<2, 3, 2, 4>,
460 NULL,
461 x, z);
462
463 // h(y, z) in R^4
464 problem.AddResidualBlock(new ParameterIgnoringCostFunction<3, 4, 3, 4>,
465 NULL,
466 y, z);
467
468 scoped_ptr<Evaluator> evaluator(CreateEvaluator(problem.mutable_program()));
469 scoped_ptr<SparseMatrix> jacobian(evaluator->CreateJacobian());
470 ASSERT_EQ(9, jacobian->num_rows());
471 ASSERT_EQ(7, jacobian->num_cols());
472
473 // f g h
474 double expected_cost = (1 + 4 + 1 + 4 + 9 + 1 + 4 + 9 + 16) / 2.0;
475
476
477 // Cost only; no residuals and no jacobian.
478 {
479 double cost = -1;
480 ASSERT_TRUE(evaluator->Evaluate(state, &cost, NULL, NULL));
481 EXPECT_EQ(expected_cost, cost);
482 }
483
484 // Cost and residuals, no jacobian.
485 {
486 double cost = -1;
487 double residuals[9] = { -2, -2, -2, -2, -2, -2, -2, -2, -2 };
488 ASSERT_TRUE(evaluator->Evaluate(state, &cost, residuals, NULL));
489 EXPECT_EQ(expected_cost, cost);
490 EXPECT_EQ(1.0, residuals[0]);
491 EXPECT_EQ(2.0, residuals[1]);
492 EXPECT_EQ(1.0, residuals[2]);
493 EXPECT_EQ(2.0, residuals[3]);
494 EXPECT_EQ(3.0, residuals[4]);
495 EXPECT_EQ(1.0, residuals[5]);
496 EXPECT_EQ(2.0, residuals[6]);
497 EXPECT_EQ(3.0, residuals[7]);
498 EXPECT_EQ(4.0, residuals[8]);
499 }
500
501 // Cost, residuals, and jacobian.
502 {
503 double cost = -1;
504 double residuals[9] = { -2, -2, -2, -2, -2, -2, -2, -2, -2 };
505 SetSparseMatrixConstant(jacobian.get(), -1);
506 ASSERT_TRUE(evaluator->Evaluate(state, &cost, residuals, jacobian.get()));
507 EXPECT_EQ(expected_cost, cost);
508 EXPECT_EQ(1.0, residuals[0]);
509 EXPECT_EQ(2.0, residuals[1]);
510 EXPECT_EQ(1.0, residuals[2]);
511 EXPECT_EQ(2.0, residuals[3]);
512 EXPECT_EQ(3.0, residuals[4]);
513 EXPECT_EQ(1.0, residuals[5]);
514 EXPECT_EQ(2.0, residuals[6]);
515 EXPECT_EQ(3.0, residuals[7]);
516 EXPECT_EQ(4.0, residuals[8]);
517
518 Matrix actual_jacobian;
519 jacobian->ToDenseMatrix(&actual_jacobian);
520
521 // Note y and z are missing columns due to the subset parameterization.
522 Matrix expected_jacobian(9, 7);
523 expected_jacobian <<
524 // x y z
525 /* f(x, y) */ 1, 2, 2, 3, 0, 0, 0,
526 1, 2, 2, 3, 0, 0, 0,
527
528 /* g(x, z) */ 2, 4, 0, 0, 2, 6, 8,
529 2, 4, 0, 0, 2, 6, 8,
530 2, 4, 0, 0, 2, 6, 8,
531
532 /* h(y, z) */ 0, 0, 6, 9, 3, 9, 12,
533 0, 0, 6, 9, 3, 9, 12,
534 0, 0, 6, 9, 3, 9, 12,
535 0, 0, 6, 9, 3, 9, 12;
536
537 EXPECT_TRUE((actual_jacobian.array() == expected_jacobian.array()).all())
538 << "Actual:\n" << actual_jacobian
539 << "\nExpected:\n" << expected_jacobian;
540 }
541}
542
543TEST_P(EvaluatorTest, MultipleResidualProblemWithSomeConstantParameters) {
544 ProblemImpl problem;
545
546 // The values are ignored completely by the cost function.
547 double x[2];
548 double y[3];
549 double z[4];
550 double state[9];
551
552 // Add the parameters in explicit order to force the ordering in the program.
553 problem.AddParameterBlock(x, 2);
554 problem.AddParameterBlock(y, 3);
555 problem.AddParameterBlock(z, 4);
556
557 // f(x, y) in R^2
558 problem.AddResidualBlock(new ParameterIgnoringCostFunction<1, 2, 2, 3>,
559 NULL,
560 x, y);
561
562 // g(x, z) in R^3
563 problem.AddResidualBlock(new ParameterIgnoringCostFunction<2, 3, 2, 4>,
564 NULL,
565 x, z);
566
567 // h(y, z) in R^4
568 problem.AddResidualBlock(new ParameterIgnoringCostFunction<3, 4, 3, 4>,
569 NULL,
570 y, z);
571
572 // For this test, "z" is constant.
573 problem.SetParameterBlockConstant(z);
574
575 // Create the reduced program which is missing the fixed "z" variable.
576 // Normally, the preprocessing of the program that happens in solver_impl
577 // takes care of this, but we don't want to invoke the solver here.
578 Program reduced_program;
579 *reduced_program.mutable_residual_blocks() =
580 problem.program().residual_blocks();
581 *reduced_program.mutable_parameter_blocks() =
582 problem.program().parameter_blocks();
583
584 // "z" is the last parameter; pop it off.
585 reduced_program.mutable_parameter_blocks()->pop_back();
586
587 scoped_ptr<Evaluator> evaluator(CreateEvaluator(&reduced_program));
588 scoped_ptr<SparseMatrix> jacobian(evaluator->CreateJacobian());
589 ASSERT_EQ(9, jacobian->num_rows());
590 ASSERT_EQ(5, jacobian->num_cols());
591
592 // f g h
593 double expected_cost = (1 + 4 + 1 + 4 + 9 + 1 + 4 + 9 + 16) / 2.0;
594
595
596 // Cost only; no residuals and no jacobian.
597 {
598 double cost = -1;
599 ASSERT_TRUE(evaluator->Evaluate(state, &cost, NULL, NULL));
600 EXPECT_EQ(expected_cost, cost);
601 }
602
603 // Cost and residuals, no jacobian.
604 {
605 double cost = -1;
606 double residuals[9] = { -2, -2, -2, -2, -2, -2, -2, -2, -2 };
607 ASSERT_TRUE(evaluator->Evaluate(state, &cost, residuals, NULL));
608 EXPECT_EQ(expected_cost, cost);
609 EXPECT_EQ(1.0, residuals[0]);
610 EXPECT_EQ(2.0, residuals[1]);
611 EXPECT_EQ(1.0, residuals[2]);
612 EXPECT_EQ(2.0, residuals[3]);
613 EXPECT_EQ(3.0, residuals[4]);
614 EXPECT_EQ(1.0, residuals[5]);
615 EXPECT_EQ(2.0, residuals[6]);
616 EXPECT_EQ(3.0, residuals[7]);
617 EXPECT_EQ(4.0, residuals[8]);
618 }
619
620 // Cost, residuals, and jacobian.
621 {
622 double cost = -1;
623 double residuals[9] = { -2, -2, -2, -2, -2, -2, -2, -2, -2 };
624 SetSparseMatrixConstant(jacobian.get(), -1);
625 ASSERT_TRUE(evaluator->Evaluate(state, &cost, residuals, jacobian.get()));
626 EXPECT_EQ(expected_cost, cost);
627 EXPECT_EQ(1.0, residuals[0]);
628 EXPECT_EQ(2.0, residuals[1]);
629 EXPECT_EQ(1.0, residuals[2]);
630 EXPECT_EQ(2.0, residuals[3]);
631 EXPECT_EQ(3.0, residuals[4]);
632 EXPECT_EQ(1.0, residuals[5]);
633 EXPECT_EQ(2.0, residuals[6]);
634 EXPECT_EQ(3.0, residuals[7]);
635 EXPECT_EQ(4.0, residuals[8]);
636
637 Matrix actual_jacobian;
638 jacobian->ToDenseMatrix(&actual_jacobian);
639
640 Matrix expected_jacobian(9, 5);
641 expected_jacobian <<
642 // x y
643 /* f(x, y) */ 1, 2, 1, 2, 3,
644 1, 2, 1, 2, 3,
645
646 /* g(x, z) */ 2, 4, 0, 0, 0,
647 2, 4, 0, 0, 0,
648 2, 4, 0, 0, 0,
649
650 /* h(y, z) */ 0, 0, 3, 6, 9,
651 0, 0, 3, 6, 9,
652 0, 0, 3, 6, 9,
653 0, 0, 3, 6, 9;
654
655 EXPECT_TRUE((actual_jacobian.array() == expected_jacobian.array()).all())
656 << "Actual:\n" << actual_jacobian
657 << "\nExpected:\n" << expected_jacobian;
658 }
659}
660
661TEST_P(EvaluatorTest, EvaluatorAbortsForResidualsThatFailToEvaluate) {
662 ProblemImpl problem;
663
664 // The values are ignored completely by the cost function.
665 double x[2];
666 double y[3];
667 double z[4];
668 double state[9];
669
670 // Switch the return value to failure.
671 problem.AddResidualBlock(
672 new ParameterIgnoringCostFunction<20, 3, 2, 3, 4, false>, NULL, x, y, z);
673
674 scoped_ptr<Evaluator> evaluator(CreateEvaluator(problem.mutable_program()));
675 scoped_ptr<SparseMatrix> jacobian(evaluator->CreateJacobian());
676 double cost;
677 EXPECT_FALSE(evaluator->Evaluate(state, &cost, NULL, NULL));
678}
679
680// In the pairs, the first argument is the linear solver type, and the second
681// argument is num_eliminate_blocks. Changing the num_eliminate_blocks only
682// makes sense for the schur-based solvers.
683//
684// Try all values of num_eliminate_blocks that make sense given that in the
685// tests a maximum of 4 parameter blocks are present.
686INSTANTIATE_TEST_CASE_P(
687 LinearSolvers,
688 EvaluatorTest,
689 ::testing::Values(make_pair(DENSE_QR, 0),
690 make_pair(DENSE_SCHUR, 0),
691 make_pair(DENSE_SCHUR, 1),
692 make_pair(DENSE_SCHUR, 2),
693 make_pair(DENSE_SCHUR, 3),
694 make_pair(DENSE_SCHUR, 4),
695 make_pair(SPARSE_SCHUR, 0),
696 make_pair(SPARSE_SCHUR, 1),
697 make_pair(SPARSE_SCHUR, 2),
698 make_pair(SPARSE_SCHUR, 3),
699 make_pair(SPARSE_SCHUR, 4),
700 make_pair(ITERATIVE_SCHUR, 0),
701 make_pair(ITERATIVE_SCHUR, 1),
702 make_pair(ITERATIVE_SCHUR, 2),
703 make_pair(ITERATIVE_SCHUR, 3),
704 make_pair(ITERATIVE_SCHUR, 4),
705 make_pair(SPARSE_NORMAL_CHOLESKY, 0)));
706
707// Simple cost function used to check if the evaluator is sensitive to
708// state changes.
709class ParameterSensitiveCostFunction : public SizedCostFunction<2, 2> {
710 public:
711 virtual bool Evaluate(double const* const* parameters,
712 double* residuals,
713 double** jacobians) const {
714 double x1 = parameters[0][0];
715 double x2 = parameters[0][1];
716 residuals[0] = x1 * x1;
717 residuals[1] = x2 * x2;
718
719 if (jacobians != NULL) {
720 double* jacobian = jacobians[0];
721 if (jacobian != NULL) {
722 jacobian[0] = 2.0 * x1;
723 jacobian[1] = 0.0;
724 jacobian[2] = 0.0;
725 jacobian[3] = 2.0 * x2;
726 }
727 }
728 return true;
729 }
730};
731
732TEST(Evaluator, EvaluatorRespectsParameterChanges) {
733 ProblemImpl problem;
734
735 double x[2];
736 x[0] = 1.0;
737 x[1] = 1.0;
738
739 problem.AddResidualBlock(new ParameterSensitiveCostFunction(), NULL, x);
740 Program* program = problem.mutable_program();
741 program->SetParameterOffsetsAndIndex();
742
743 Evaluator::Options options;
744 options.linear_solver_type = DENSE_QR;
745 options.num_eliminate_blocks = 0;
746 string error;
747 scoped_ptr<Evaluator> evaluator(Evaluator::Create(options, program, &error));
748 scoped_ptr<SparseMatrix> jacobian(evaluator->CreateJacobian());
749
750 ASSERT_EQ(2, jacobian->num_rows());
751 ASSERT_EQ(2, jacobian->num_cols());
752
753 double state[2];
754 state[0] = 2.0;
755 state[1] = 3.0;
756
757 // The original state of a residual block comes from the user's
758 // state. So the original state is 1.0, 1.0, and the only way we get
759 // the 2.0, 3.0 results in the following tests is if it respects the
760 // values in the state vector.
761
762 // Cost only; no residuals and no jacobian.
763 {
764 double cost = -1;
765 ASSERT_TRUE(evaluator->Evaluate(state, &cost, NULL, NULL));
766 EXPECT_EQ(48.5, cost);
767 }
768
769 // Cost and residuals, no jacobian.
770 {
771 double cost = -1;
772 double residuals[2] = { -2, -2 };
773 ASSERT_TRUE(evaluator->Evaluate(state, &cost, residuals, NULL));
774 EXPECT_EQ(48.5, cost);
775 EXPECT_EQ(4, residuals[0]);
776 EXPECT_EQ(9, residuals[1]);
777 }
778
779 // Cost, residuals, and jacobian.
780 {
781 double cost = -1;
782 double residuals[2] = { -2, -2};
783 SetSparseMatrixConstant(jacobian.get(), -1);
784 ASSERT_TRUE(evaluator->Evaluate(state, &cost, residuals, jacobian.get()));
785 EXPECT_EQ(48.5, cost);
786 EXPECT_EQ(4, residuals[0]);
787 EXPECT_EQ(9, residuals[1]);
788 Matrix actual_jacobian;
789 jacobian->ToDenseMatrix(&actual_jacobian);
790
791 Matrix expected_jacobian(2, 2);
792 expected_jacobian
793 << 2 * state[0], 0,
794 0, 2 * state[1];
795
796 EXPECT_TRUE((actual_jacobian.array() == expected_jacobian.array()).all())
797 << "Actual:\n" << actual_jacobian
798 << "\nExpected:\n" << expected_jacobian;
799 }
800}
801
802} // namespace internal
803} // namespace ceres