blob: 5ddb73375690bb830493ab9f565a96490499a68e [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: 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 <memory>
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080037
Austin Schuh70cc9552019-01-21 19:46:48 -080038#include "ceres/casts.h"
39#include "ceres/cost_function.h"
40#include "ceres/crs_matrix.h"
41#include "ceres/evaluator_test_utils.h"
42#include "ceres/internal/eigen.h"
43#include "ceres/local_parameterization.h"
44#include "ceres/problem_impl.h"
45#include "ceres/program.h"
46#include "ceres/sized_cost_function.h"
47#include "ceres/sparse_matrix.h"
48#include "ceres/stringprintf.h"
49#include "ceres/types.h"
50#include "gtest/gtest.h"
51
52namespace ceres {
53namespace internal {
54
55using std::string;
56using std::vector;
57
58// TODO(keir): Consider pushing this into a common test utils file.
59template <int kFactor, int kNumResiduals, int... Ns>
60class ParameterIgnoringCostFunction
61 : public SizedCostFunction<kNumResiduals, Ns...> {
62 typedef SizedCostFunction<kNumResiduals, Ns...> Base;
63
64 public:
65 explicit ParameterIgnoringCostFunction(bool succeeds = true)
66 : succeeds_(succeeds) {}
67
Austin Schuh1d1e6ea2020-12-23 21:56:30 -080068 bool Evaluate(double const* const* parameters,
69 double* residuals,
70 double** jacobians) const final {
Austin Schuh70cc9552019-01-21 19:46:48 -080071 for (int i = 0; i < Base::num_residuals(); ++i) {
72 residuals[i] = i + 1;
73 }
74 if (jacobians) {
75 for (int k = 0; k < Base::parameter_block_sizes().size(); ++k) {
76 // The jacobians here are full sized, but they are transformed in the
77 // evaluator into the "local" jacobian. In the tests, the "subset
78 // constant" parameterization is used, which should pick out columns
79 // from these jacobians. Put values in the jacobian that make this
80 // obvious; in particular, make the jacobians like this:
81 //
82 // 1 2 3 4 ...
83 // 1 2 3 4 ... .* kFactor
84 // 1 2 3 4 ...
85 //
86 // where the multiplication by kFactor makes it easier to distinguish
87 // between Jacobians of different residuals for the same parameter.
88 if (jacobians[k] != nullptr) {
89 MatrixRef jacobian(jacobians[k],
90 Base::num_residuals(),
91 Base::parameter_block_sizes()[k]);
92 for (int j = 0; j < Base::parameter_block_sizes()[k]; ++j) {
93 jacobian.col(j).setConstant(kFactor * (j + 1));
94 }
95 }
96 }
97 }
98 return succeeds_;
99 }
100
101 private:
102 bool succeeds_;
103};
104
105struct EvaluatorTestOptions {
106 EvaluatorTestOptions(LinearSolverType linear_solver_type,
107 int num_eliminate_blocks,
108 bool dynamic_sparsity = false)
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800109 : linear_solver_type(linear_solver_type),
110 num_eliminate_blocks(num_eliminate_blocks),
111 dynamic_sparsity(dynamic_sparsity) {}
Austin Schuh70cc9552019-01-21 19:46:48 -0800112
113 LinearSolverType linear_solver_type;
114 int num_eliminate_blocks;
115 bool dynamic_sparsity;
116};
117
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800118struct EvaluatorTest : public ::testing::TestWithParam<EvaluatorTestOptions> {
Austin Schuh70cc9552019-01-21 19:46:48 -0800119 Evaluator* CreateEvaluator(Program* program) {
120 // This program is straight from the ProblemImpl, and so has no index/offset
121 // yet; compute it here as required by the evaluator implementations.
122 program->SetParameterOffsetsAndIndex();
123
124 if (VLOG_IS_ON(1)) {
125 string report;
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800126 StringAppendF(&report,
127 "Creating evaluator with type: %d",
Austin Schuh70cc9552019-01-21 19:46:48 -0800128 GetParam().linear_solver_type);
129 if (GetParam().linear_solver_type == SPARSE_NORMAL_CHOLESKY) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800130 StringAppendF(
131 &report, ", dynamic_sparsity: %d", GetParam().dynamic_sparsity);
Austin Schuh70cc9552019-01-21 19:46:48 -0800132 }
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800133 StringAppendF(&report,
134 " and num_eliminate_blocks: %d",
Austin Schuh70cc9552019-01-21 19:46:48 -0800135 GetParam().num_eliminate_blocks);
136 VLOG(1) << report;
137 }
138 Evaluator::Options options;
139 options.linear_solver_type = GetParam().linear_solver_type;
140 options.num_eliminate_blocks = GetParam().num_eliminate_blocks;
141 options.dynamic_sparsity = GetParam().dynamic_sparsity;
142 options.context = problem.context();
143 string error;
144 return Evaluator::Create(options, program, &error);
145 }
146
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800147 void EvaluateAndCompare(ProblemImpl* problem,
Austin Schuh70cc9552019-01-21 19:46:48 -0800148 int expected_num_rows,
149 int expected_num_cols,
150 double expected_cost,
151 const double* expected_residuals,
152 const double* expected_gradient,
153 const double* expected_jacobian) {
154 std::unique_ptr<Evaluator> evaluator(
155 CreateEvaluator(problem->mutable_program()));
156 int num_residuals = expected_num_rows;
157 int num_parameters = expected_num_cols;
158
159 double cost = -1;
160
161 Vector residuals(num_residuals);
162 residuals.setConstant(-2000);
163
164 Vector gradient(num_parameters);
165 gradient.setConstant(-3000);
166
167 std::unique_ptr<SparseMatrix> jacobian(evaluator->CreateJacobian());
168
169 ASSERT_EQ(expected_num_rows, evaluator->NumResiduals());
170 ASSERT_EQ(expected_num_cols, evaluator->NumEffectiveParameters());
171 ASSERT_EQ(expected_num_rows, jacobian->num_rows());
172 ASSERT_EQ(expected_num_cols, jacobian->num_cols());
173
174 vector<double> state(evaluator->NumParameters());
175
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800176 // clang-format off
Austin Schuh70cc9552019-01-21 19:46:48 -0800177 ASSERT_TRUE(evaluator->Evaluate(
178 &state[0],
179 &cost,
180 expected_residuals != nullptr ? &residuals[0] : nullptr,
181 expected_gradient != nullptr ? &gradient[0] : nullptr,
182 expected_jacobian != nullptr ? jacobian.get() : nullptr));
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800183 // clang-format on
Austin Schuh70cc9552019-01-21 19:46:48 -0800184
185 Matrix actual_jacobian;
186 if (expected_jacobian != nullptr) {
187 jacobian->ToDenseMatrix(&actual_jacobian);
188 }
189
190 CompareEvaluations(expected_num_rows,
191 expected_num_cols,
192 expected_cost,
193 expected_residuals,
194 expected_gradient,
195 expected_jacobian,
196 cost,
197 &residuals[0],
198 &gradient[0],
199 actual_jacobian.data());
200 }
201
202 // Try all combinations of parameters for the evaluator.
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800203 void CheckAllEvaluationCombinations(const ExpectedEvaluation& expected) {
Austin Schuh70cc9552019-01-21 19:46:48 -0800204 for (int i = 0; i < 8; ++i) {
205 EvaluateAndCompare(&problem,
206 expected.num_rows,
207 expected.num_cols,
208 expected.cost,
209 (i & 1) ? expected.residuals : nullptr,
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800210 (i & 2) ? expected.gradient : nullptr,
211 (i & 4) ? expected.jacobian : nullptr);
Austin Schuh70cc9552019-01-21 19:46:48 -0800212 }
213 }
214
215 // The values are ignored completely by the cost function.
216 double x[2];
217 double y[3];
218 double z[4];
219
220 ProblemImpl problem;
221};
222
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800223static void SetSparseMatrixConstant(SparseMatrix* sparse_matrix, double value) {
224 VectorRef(sparse_matrix->mutable_values(), sparse_matrix->num_nonzeros())
225 .setConstant(value);
Austin Schuh70cc9552019-01-21 19:46:48 -0800226}
227
228TEST_P(EvaluatorTest, SingleResidualProblem) {
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800229 problem.AddResidualBlock(
230 new ParameterIgnoringCostFunction<1, 3, 2, 3, 4>, nullptr, x, y, z);
Austin Schuh70cc9552019-01-21 19:46:48 -0800231
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800232 // clang-format off
Austin Schuh70cc9552019-01-21 19:46:48 -0800233 ExpectedEvaluation expected = {
234 // Rows/columns
235 3, 9,
236 // Cost
237 7.0,
238 // Residuals
239 { 1.0, 2.0, 3.0 },
240 // Gradient
241 { 6.0, 12.0, // x
242 6.0, 12.0, 18.0, // y
243 6.0, 12.0, 18.0, 24.0, // z
244 },
245 // Jacobian
246 // x y z
247 { 1, 2, 1, 2, 3, 1, 2, 3, 4,
248 1, 2, 1, 2, 3, 1, 2, 3, 4,
249 1, 2, 1, 2, 3, 1, 2, 3, 4
250 }
251 };
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800252 // clang-format on
Austin Schuh70cc9552019-01-21 19:46:48 -0800253 CheckAllEvaluationCombinations(expected);
254}
255
256TEST_P(EvaluatorTest, SingleResidualProblemWithPermutedParameters) {
257 // Add the parameters in explicit order to force the ordering in the program.
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800258 problem.AddParameterBlock(x, 2);
259 problem.AddParameterBlock(y, 3);
260 problem.AddParameterBlock(z, 4);
Austin Schuh70cc9552019-01-21 19:46:48 -0800261
262 // Then use a cost function which is similar to the others, but swap around
263 // the ordering of the parameters to the cost function. This shouldn't affect
264 // the jacobian evaluation, but requires explicit handling in the evaluators.
265 // At one point the compressed row evaluator had a bug that went undetected
266 // for a long time, since by chance most users added parameters to the problem
267 // in the same order that they occurred as parameters to a cost function.
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800268 problem.AddResidualBlock(
269 new ParameterIgnoringCostFunction<1, 3, 4, 3, 2>, nullptr, z, y, x);
Austin Schuh70cc9552019-01-21 19:46:48 -0800270
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800271 // clang-format off
Austin Schuh70cc9552019-01-21 19:46:48 -0800272 ExpectedEvaluation expected = {
273 // Rows/columns
274 3, 9,
275 // Cost
276 7.0,
277 // Residuals
278 { 1.0, 2.0, 3.0 },
279 // Gradient
280 { 6.0, 12.0, // x
281 6.0, 12.0, 18.0, // y
282 6.0, 12.0, 18.0, 24.0, // z
283 },
284 // Jacobian
285 // x y z
286 { 1, 2, 1, 2, 3, 1, 2, 3, 4,
287 1, 2, 1, 2, 3, 1, 2, 3, 4,
288 1, 2, 1, 2, 3, 1, 2, 3, 4
289 }
290 };
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800291 // clang-format on
Austin Schuh70cc9552019-01-21 19:46:48 -0800292 CheckAllEvaluationCombinations(expected);
293}
294
295TEST_P(EvaluatorTest, SingleResidualProblemWithNuisanceParameters) {
296 // These parameters are not used.
297 double a[2];
298 double b[1];
299 double c[1];
300 double d[3];
301
302 // Add the parameters in a mixed order so the Jacobian is "checkered" with the
303 // values from the other parameters.
304 problem.AddParameterBlock(a, 2);
305 problem.AddParameterBlock(x, 2);
306 problem.AddParameterBlock(b, 1);
307 problem.AddParameterBlock(y, 3);
308 problem.AddParameterBlock(c, 1);
309 problem.AddParameterBlock(z, 4);
310 problem.AddParameterBlock(d, 3);
311
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800312 problem.AddResidualBlock(
313 new ParameterIgnoringCostFunction<1, 3, 2, 3, 4>, nullptr, x, y, z);
Austin Schuh70cc9552019-01-21 19:46:48 -0800314
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800315 // clang-format off
Austin Schuh70cc9552019-01-21 19:46:48 -0800316 ExpectedEvaluation expected = {
317 // Rows/columns
318 3, 16,
319 // Cost
320 7.0,
321 // Residuals
322 { 1.0, 2.0, 3.0 },
323 // Gradient
324 { 0.0, 0.0, // a
325 6.0, 12.0, // x
326 0.0, // b
327 6.0, 12.0, 18.0, // y
328 0.0, // c
329 6.0, 12.0, 18.0, 24.0, // z
330 0.0, 0.0, 0.0, // d
331 },
332 // Jacobian
333 // a x b y c z d
334 { 0, 0, 1, 2, 0, 1, 2, 3, 0, 1, 2, 3, 4, 0, 0, 0,
335 0, 0, 1, 2, 0, 1, 2, 3, 0, 1, 2, 3, 4, 0, 0, 0,
336 0, 0, 1, 2, 0, 1, 2, 3, 0, 1, 2, 3, 4, 0, 0, 0
337 }
338 };
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800339 // clang-format on
Austin Schuh70cc9552019-01-21 19:46:48 -0800340 CheckAllEvaluationCombinations(expected);
341}
342
343TEST_P(EvaluatorTest, MultipleResidualProblem) {
344 // Add the parameters in explicit order to force the ordering in the program.
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800345 problem.AddParameterBlock(x, 2);
346 problem.AddParameterBlock(y, 3);
347 problem.AddParameterBlock(z, 4);
Austin Schuh70cc9552019-01-21 19:46:48 -0800348
349 // f(x, y) in R^2
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800350 problem.AddResidualBlock(
351 new ParameterIgnoringCostFunction<1, 2, 2, 3>, nullptr, x, y);
Austin Schuh70cc9552019-01-21 19:46:48 -0800352
353 // g(x, z) in R^3
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800354 problem.AddResidualBlock(
355 new ParameterIgnoringCostFunction<2, 3, 2, 4>, nullptr, x, z);
Austin Schuh70cc9552019-01-21 19:46:48 -0800356
357 // h(y, z) in R^4
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800358 problem.AddResidualBlock(
359 new ParameterIgnoringCostFunction<3, 4, 3, 4>, nullptr, y, z);
Austin Schuh70cc9552019-01-21 19:46:48 -0800360
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800361 // clang-format off
Austin Schuh70cc9552019-01-21 19:46:48 -0800362 ExpectedEvaluation expected = {
363 // Rows/columns
364 9, 9,
365 // Cost
366 // f g h
367 ( 1 + 4 + 1 + 4 + 9 + 1 + 4 + 9 + 16) / 2.0,
368 // Residuals
369 { 1.0, 2.0, // f
370 1.0, 2.0, 3.0, // g
371 1.0, 2.0, 3.0, 4.0 // h
372 },
373 // Gradient
374 { 15.0, 30.0, // x
375 33.0, 66.0, 99.0, // y
376 42.0, 84.0, 126.0, 168.0 // z
377 },
378 // Jacobian
379 // x y z
380 { /* f(x, y) */ 1, 2, 1, 2, 3, 0, 0, 0, 0,
381 1, 2, 1, 2, 3, 0, 0, 0, 0,
382
383 /* g(x, z) */ 2, 4, 0, 0, 0, 2, 4, 6, 8,
384 2, 4, 0, 0, 0, 2, 4, 6, 8,
385 2, 4, 0, 0, 0, 2, 4, 6, 8,
386
387 /* h(y, z) */ 0, 0, 3, 6, 9, 3, 6, 9, 12,
388 0, 0, 3, 6, 9, 3, 6, 9, 12,
389 0, 0, 3, 6, 9, 3, 6, 9, 12,
390 0, 0, 3, 6, 9, 3, 6, 9, 12
391 }
392 };
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800393 // clang-format on
Austin Schuh70cc9552019-01-21 19:46:48 -0800394 CheckAllEvaluationCombinations(expected);
395}
396
397TEST_P(EvaluatorTest, MultipleResidualsWithLocalParameterizations) {
398 // Add the parameters in explicit order to force the ordering in the program.
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800399 problem.AddParameterBlock(x, 2);
Austin Schuh70cc9552019-01-21 19:46:48 -0800400
401 // Fix y's first dimension.
402 vector<int> y_fixed;
403 y_fixed.push_back(0);
404 problem.AddParameterBlock(y, 3, new SubsetParameterization(3, y_fixed));
405
406 // Fix z's second dimension.
407 vector<int> z_fixed;
408 z_fixed.push_back(1);
409 problem.AddParameterBlock(z, 4, new SubsetParameterization(4, z_fixed));
410
411 // f(x, y) in R^2
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800412 problem.AddResidualBlock(
413 new ParameterIgnoringCostFunction<1, 2, 2, 3>, nullptr, x, y);
Austin Schuh70cc9552019-01-21 19:46:48 -0800414
415 // g(x, z) in R^3
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800416 problem.AddResidualBlock(
417 new ParameterIgnoringCostFunction<2, 3, 2, 4>, nullptr, x, z);
Austin Schuh70cc9552019-01-21 19:46:48 -0800418
419 // h(y, z) in R^4
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800420 problem.AddResidualBlock(
421 new ParameterIgnoringCostFunction<3, 4, 3, 4>, nullptr, y, z);
Austin Schuh70cc9552019-01-21 19:46:48 -0800422
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800423 // clang-format off
Austin Schuh70cc9552019-01-21 19:46:48 -0800424 ExpectedEvaluation expected = {
425 // Rows/columns
426 9, 7,
427 // Cost
428 // f g h
429 ( 1 + 4 + 1 + 4 + 9 + 1 + 4 + 9 + 16) / 2.0,
430 // Residuals
431 { 1.0, 2.0, // f
432 1.0, 2.0, 3.0, // g
433 1.0, 2.0, 3.0, 4.0 // h
434 },
435 // Gradient
436 { 15.0, 30.0, // x
437 66.0, 99.0, // y
438 42.0, 126.0, 168.0 // z
439 },
440 // Jacobian
441 // x y z
442 { /* f(x, y) */ 1, 2, 2, 3, 0, 0, 0,
443 1, 2, 2, 3, 0, 0, 0,
444
445 /* g(x, z) */ 2, 4, 0, 0, 2, 6, 8,
446 2, 4, 0, 0, 2, 6, 8,
447 2, 4, 0, 0, 2, 6, 8,
448
449 /* h(y, z) */ 0, 0, 6, 9, 3, 9, 12,
450 0, 0, 6, 9, 3, 9, 12,
451 0, 0, 6, 9, 3, 9, 12,
452 0, 0, 6, 9, 3, 9, 12
453 }
454 };
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800455 // clang-format on
Austin Schuh70cc9552019-01-21 19:46:48 -0800456 CheckAllEvaluationCombinations(expected);
457}
458
459TEST_P(EvaluatorTest, MultipleResidualProblemWithSomeConstantParameters) {
460 // The values are ignored completely by the cost function.
461 double x[2];
462 double y[3];
463 double z[4];
464
465 // Add the parameters in explicit order to force the ordering in the program.
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800466 problem.AddParameterBlock(x, 2);
467 problem.AddParameterBlock(y, 3);
468 problem.AddParameterBlock(z, 4);
Austin Schuh70cc9552019-01-21 19:46:48 -0800469
470 // f(x, y) in R^2
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800471 problem.AddResidualBlock(
472 new ParameterIgnoringCostFunction<1, 2, 2, 3>, nullptr, x, y);
Austin Schuh70cc9552019-01-21 19:46:48 -0800473
474 // g(x, z) in R^3
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800475 problem.AddResidualBlock(
476 new ParameterIgnoringCostFunction<2, 3, 2, 4>, nullptr, x, z);
Austin Schuh70cc9552019-01-21 19:46:48 -0800477
478 // h(y, z) in R^4
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800479 problem.AddResidualBlock(
480 new ParameterIgnoringCostFunction<3, 4, 3, 4>, nullptr, y, z);
Austin Schuh70cc9552019-01-21 19:46:48 -0800481
482 // For this test, "z" is constant.
483 problem.SetParameterBlockConstant(z);
484
485 // Create the reduced program which is missing the fixed "z" variable.
486 // Normally, the preprocessing of the program that happens in solver_impl
487 // takes care of this, but we don't want to invoke the solver here.
488 Program reduced_program;
489 vector<ParameterBlock*>* parameter_blocks =
490 problem.mutable_program()->mutable_parameter_blocks();
491
492 // "z" is the last parameter; save it for later and pop it off temporarily.
493 // Note that "z" will still get read during evaluation, so it cannot be
494 // deleted at this point.
495 ParameterBlock* parameter_block_z = parameter_blocks->back();
496 parameter_blocks->pop_back();
497
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800498 // clang-format off
Austin Schuh70cc9552019-01-21 19:46:48 -0800499 ExpectedEvaluation expected = {
500 // Rows/columns
501 9, 5,
502 // Cost
503 // f g h
504 ( 1 + 4 + 1 + 4 + 9 + 1 + 4 + 9 + 16) / 2.0,
505 // Residuals
506 { 1.0, 2.0, // f
507 1.0, 2.0, 3.0, // g
508 1.0, 2.0, 3.0, 4.0 // h
509 },
510 // Gradient
511 { 15.0, 30.0, // x
512 33.0, 66.0, 99.0, // y
513 },
514 // Jacobian
515 // x y
516 { /* f(x, y) */ 1, 2, 1, 2, 3,
517 1, 2, 1, 2, 3,
518
519 /* g(x, z) */ 2, 4, 0, 0, 0,
520 2, 4, 0, 0, 0,
521 2, 4, 0, 0, 0,
522
523 /* h(y, z) */ 0, 0, 3, 6, 9,
524 0, 0, 3, 6, 9,
525 0, 0, 3, 6, 9,
526 0, 0, 3, 6, 9
527 }
528 };
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800529 // clang-format on
Austin Schuh70cc9552019-01-21 19:46:48 -0800530 CheckAllEvaluationCombinations(expected);
531
532 // Restore parameter block z, so it will get freed in a consistent way.
533 parameter_blocks->push_back(parameter_block_z);
534}
535
536TEST_P(EvaluatorTest, EvaluatorAbortsForResidualsThatFailToEvaluate) {
537 // Switch the return value to failure.
538 problem.AddResidualBlock(
539 new ParameterIgnoringCostFunction<20, 3, 2, 3, 4>(false),
540 nullptr,
541 x,
542 y,
543 z);
544
545 // The values are ignored.
546 double state[9];
547
548 std::unique_ptr<Evaluator> evaluator(
549 CreateEvaluator(problem.mutable_program()));
550 std::unique_ptr<SparseMatrix> jacobian(evaluator->CreateJacobian());
551 double cost;
552 EXPECT_FALSE(evaluator->Evaluate(state, &cost, nullptr, nullptr, nullptr));
553}
554
555// In the pairs, the first argument is the linear solver type, and the second
556// argument is num_eliminate_blocks. Changing the num_eliminate_blocks only
557// makes sense for the schur-based solvers.
558//
559// Try all values of num_eliminate_blocks that make sense given that in the
560// tests a maximum of 4 parameter blocks are present.
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800561INSTANTIATE_TEST_SUITE_P(
Austin Schuh70cc9552019-01-21 19:46:48 -0800562 LinearSolvers,
563 EvaluatorTest,
564 ::testing::Values(EvaluatorTestOptions(DENSE_QR, 0),
565 EvaluatorTestOptions(DENSE_SCHUR, 0),
566 EvaluatorTestOptions(DENSE_SCHUR, 1),
567 EvaluatorTestOptions(DENSE_SCHUR, 2),
568 EvaluatorTestOptions(DENSE_SCHUR, 3),
569 EvaluatorTestOptions(DENSE_SCHUR, 4),
570 EvaluatorTestOptions(SPARSE_SCHUR, 0),
571 EvaluatorTestOptions(SPARSE_SCHUR, 1),
572 EvaluatorTestOptions(SPARSE_SCHUR, 2),
573 EvaluatorTestOptions(SPARSE_SCHUR, 3),
574 EvaluatorTestOptions(SPARSE_SCHUR, 4),
575 EvaluatorTestOptions(ITERATIVE_SCHUR, 0),
576 EvaluatorTestOptions(ITERATIVE_SCHUR, 1),
577 EvaluatorTestOptions(ITERATIVE_SCHUR, 2),
578 EvaluatorTestOptions(ITERATIVE_SCHUR, 3),
579 EvaluatorTestOptions(ITERATIVE_SCHUR, 4),
580 EvaluatorTestOptions(SPARSE_NORMAL_CHOLESKY, 0, false),
581 EvaluatorTestOptions(SPARSE_NORMAL_CHOLESKY, 0, true)));
582
583// Simple cost function used to check if the evaluator is sensitive to
584// state changes.
585class ParameterSensitiveCostFunction : public SizedCostFunction<2, 2> {
586 public:
Austin Schuh1d1e6ea2020-12-23 21:56:30 -0800587 bool Evaluate(double const* const* parameters,
588 double* residuals,
589 double** jacobians) const final {
Austin Schuh70cc9552019-01-21 19:46:48 -0800590 double x1 = parameters[0][0];
591 double x2 = parameters[0][1];
592 residuals[0] = x1 * x1;
593 residuals[1] = x2 * x2;
594
595 if (jacobians != nullptr) {
596 double* jacobian = jacobians[0];
597 if (jacobian != nullptr) {
598 jacobian[0] = 2.0 * x1;
599 jacobian[1] = 0.0;
600 jacobian[2] = 0.0;
601 jacobian[3] = 2.0 * x2;
602 }
603 }
604 return true;
605 }
606};
607
608TEST(Evaluator, EvaluatorRespectsParameterChanges) {
609 ProblemImpl problem;
610
611 double x[2];
612 x[0] = 1.0;
613 x[1] = 1.0;
614
615 problem.AddResidualBlock(new ParameterSensitiveCostFunction(), nullptr, x);
616 Program* program = problem.mutable_program();
617 program->SetParameterOffsetsAndIndex();
618
619 Evaluator::Options options;
620 options.linear_solver_type = DENSE_QR;
621 options.num_eliminate_blocks = 0;
622 options.context = problem.context();
623 string error;
624 std::unique_ptr<Evaluator> evaluator(
625 Evaluator::Create(options, program, &error));
626 std::unique_ptr<SparseMatrix> jacobian(evaluator->CreateJacobian());
627
628 ASSERT_EQ(2, jacobian->num_rows());
629 ASSERT_EQ(2, jacobian->num_cols());
630
631 double state[2];
632 state[0] = 2.0;
633 state[1] = 3.0;
634
635 // The original state of a residual block comes from the user's
636 // state. So the original state is 1.0, 1.0, and the only way we get
637 // the 2.0, 3.0 results in the following tests is if it respects the
638 // values in the state vector.
639
640 // Cost only; no residuals and no jacobian.
641 {
642 double cost = -1;
643 ASSERT_TRUE(evaluator->Evaluate(state, &cost, nullptr, nullptr, nullptr));
644 EXPECT_EQ(48.5, cost);
645 }
646
647 // Cost and residuals, no jacobian.
648 {
649 double cost = -1;
650 double residuals[2] = {-2, -2};
651 ASSERT_TRUE(evaluator->Evaluate(state, &cost, residuals, nullptr, nullptr));
652 EXPECT_EQ(48.5, cost);
653 EXPECT_EQ(4, residuals[0]);
654 EXPECT_EQ(9, residuals[1]);
655 }
656
657 // Cost, residuals, and jacobian.
658 {
659 double cost = -1;
660 double residuals[2] = {-2, -2};
661 SetSparseMatrixConstant(jacobian.get(), -1);
662 ASSERT_TRUE(
663 evaluator->Evaluate(state, &cost, residuals, nullptr, jacobian.get()));
664 EXPECT_EQ(48.5, cost);
665 EXPECT_EQ(4, residuals[0]);
666 EXPECT_EQ(9, residuals[1]);
667 Matrix actual_jacobian;
668 jacobian->ToDenseMatrix(&actual_jacobian);
669
670 Matrix expected_jacobian(2, 2);
671 expected_jacobian << 2 * state[0], 0, 0, 2 * state[1];
672
673 EXPECT_TRUE((actual_jacobian.array() == expected_jacobian.array()).all())
674 << "Actual:\n"
675 << actual_jacobian << "\nExpected:\n"
676 << expected_jacobian;
677 }
678}
679
680} // namespace internal
681} // namespace ceres