blob: dea07236e7a003d8cb700e01a58314cb8990b360 [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/covariance.h"
32
33#include <algorithm>
34#include <cstdint>
35#include <cmath>
36#include <map>
37#include <memory>
38#include <utility>
39
40#include "ceres/compressed_row_sparse_matrix.h"
41#include "ceres/cost_function.h"
42#include "ceres/covariance_impl.h"
43#include "ceres/local_parameterization.h"
44#include "ceres/map_util.h"
45#include "ceres/problem_impl.h"
46#include "gtest/gtest.h"
47
48namespace ceres {
49namespace internal {
50
51using std::make_pair;
52using std::map;
53using std::pair;
54using std::vector;
55
56class UnaryCostFunction: public CostFunction {
57 public:
58 UnaryCostFunction(const int num_residuals,
59 const int32_t parameter_block_size,
60 const double* jacobian)
61 : jacobian_(jacobian, jacobian + num_residuals * parameter_block_size) {
62 set_num_residuals(num_residuals);
63 mutable_parameter_block_sizes()->push_back(parameter_block_size);
64 }
65
66 virtual bool Evaluate(double const* const* parameters,
67 double* residuals,
68 double** jacobians) const {
69 for (int i = 0; i < num_residuals(); ++i) {
70 residuals[i] = 1;
71 }
72
73 if (jacobians == NULL) {
74 return true;
75 }
76
77 if (jacobians[0] != NULL) {
78 copy(jacobian_.begin(), jacobian_.end(), jacobians[0]);
79 }
80
81 return true;
82 }
83
84 private:
85 vector<double> jacobian_;
86};
87
88
89class BinaryCostFunction: public CostFunction {
90 public:
91 BinaryCostFunction(const int num_residuals,
92 const int32_t parameter_block1_size,
93 const int32_t parameter_block2_size,
94 const double* jacobian1,
95 const double* jacobian2)
96 : jacobian1_(jacobian1,
97 jacobian1 + num_residuals * parameter_block1_size),
98 jacobian2_(jacobian2,
99 jacobian2 + num_residuals * parameter_block2_size) {
100 set_num_residuals(num_residuals);
101 mutable_parameter_block_sizes()->push_back(parameter_block1_size);
102 mutable_parameter_block_sizes()->push_back(parameter_block2_size);
103 }
104
105 virtual bool Evaluate(double const* const* parameters,
106 double* residuals,
107 double** jacobians) const {
108 for (int i = 0; i < num_residuals(); ++i) {
109 residuals[i] = 2;
110 }
111
112 if (jacobians == NULL) {
113 return true;
114 }
115
116 if (jacobians[0] != NULL) {
117 copy(jacobian1_.begin(), jacobian1_.end(), jacobians[0]);
118 }
119
120 if (jacobians[1] != NULL) {
121 copy(jacobian2_.begin(), jacobian2_.end(), jacobians[1]);
122 }
123
124 return true;
125 }
126
127 private:
128 vector<double> jacobian1_;
129 vector<double> jacobian2_;
130};
131
132// x_plus_delta = delta * x;
133class PolynomialParameterization : public LocalParameterization {
134 public:
135 virtual ~PolynomialParameterization() {}
136
137 virtual bool Plus(const double* x,
138 const double* delta,
139 double* x_plus_delta) const {
140 x_plus_delta[0] = delta[0] * x[0];
141 x_plus_delta[1] = delta[0] * x[1];
142 return true;
143 }
144
145 virtual bool ComputeJacobian(const double* x, double* jacobian) const {
146 jacobian[0] = x[0];
147 jacobian[1] = x[1];
148 return true;
149 }
150
151 virtual int GlobalSize() const { return 2; }
152 virtual int LocalSize() const { return 1; }
153};
154
155TEST(CovarianceImpl, ComputeCovarianceSparsity) {
156 double parameters[10];
157
158 double* block1 = parameters;
159 double* block2 = block1 + 1;
160 double* block3 = block2 + 2;
161 double* block4 = block3 + 3;
162
163 ProblemImpl problem;
164
165 // Add in random order
166 Vector junk_jacobian = Vector::Zero(10);
167 problem.AddResidualBlock(
168 new UnaryCostFunction(1, 1, junk_jacobian.data()), NULL, block1);
169 problem.AddResidualBlock(
170 new UnaryCostFunction(1, 4, junk_jacobian.data()), NULL, block4);
171 problem.AddResidualBlock(
172 new UnaryCostFunction(1, 3, junk_jacobian.data()), NULL, block3);
173 problem.AddResidualBlock(
174 new UnaryCostFunction(1, 2, junk_jacobian.data()), NULL, block2);
175
176 // Sparsity pattern
177 //
178 // Note that the problem structure does not imply this sparsity
179 // pattern since all the residual blocks are unary. But the
180 // ComputeCovarianceSparsity function in its current incarnation
181 // does not pay attention to this fact and only looks at the
182 // parameter block pairs that the user provides.
183 //
184 // X . . . . . X X X X
185 // . X X X X X . . . .
186 // . X X X X X . . . .
187 // . . . X X X . . . .
188 // . . . X X X . . . .
189 // . . . X X X . . . .
190 // . . . . . . X X X X
191 // . . . . . . X X X X
192 // . . . . . . X X X X
193 // . . . . . . X X X X
194
195 int expected_rows[] = {0, 5, 10, 15, 18, 21, 24, 28, 32, 36, 40};
196 int expected_cols[] = {0, 6, 7, 8, 9,
197 1, 2, 3, 4, 5,
198 1, 2, 3, 4, 5,
199 3, 4, 5,
200 3, 4, 5,
201 3, 4, 5,
202 6, 7, 8, 9,
203 6, 7, 8, 9,
204 6, 7, 8, 9,
205 6, 7, 8, 9};
206
207
208 vector<pair<const double*, const double*>> covariance_blocks;
209 covariance_blocks.push_back(make_pair(block1, block1));
210 covariance_blocks.push_back(make_pair(block4, block4));
211 covariance_blocks.push_back(make_pair(block2, block2));
212 covariance_blocks.push_back(make_pair(block3, block3));
213 covariance_blocks.push_back(make_pair(block2, block3));
214 covariance_blocks.push_back(make_pair(block4, block1)); // reversed
215
216 Covariance::Options options;
217 CovarianceImpl covariance_impl(options);
218 EXPECT_TRUE(covariance_impl
219 .ComputeCovarianceSparsity(covariance_blocks, &problem));
220
221 const CompressedRowSparseMatrix* crsm = covariance_impl.covariance_matrix();
222
223 EXPECT_EQ(crsm->num_rows(), 10);
224 EXPECT_EQ(crsm->num_cols(), 10);
225 EXPECT_EQ(crsm->num_nonzeros(), 40);
226
227 const int* rows = crsm->rows();
228 for (int r = 0; r < crsm->num_rows() + 1; ++r) {
229 EXPECT_EQ(rows[r], expected_rows[r])
230 << r << " "
231 << rows[r] << " "
232 << expected_rows[r];
233 }
234
235 const int* cols = crsm->cols();
236 for (int c = 0; c < crsm->num_nonzeros(); ++c) {
237 EXPECT_EQ(cols[c], expected_cols[c])
238 << c << " "
239 << cols[c] << " "
240 << expected_cols[c];
241 }
242}
243
244TEST(CovarianceImpl, ComputeCovarianceSparsityWithConstantParameterBlock) {
245 double parameters[10];
246
247 double* block1 = parameters;
248 double* block2 = block1 + 1;
249 double* block3 = block2 + 2;
250 double* block4 = block3 + 3;
251
252 ProblemImpl problem;
253
254 // Add in random order
255 Vector junk_jacobian = Vector::Zero(10);
256 problem.AddResidualBlock(
257 new UnaryCostFunction(1, 1, junk_jacobian.data()), NULL, block1);
258 problem.AddResidualBlock(
259 new UnaryCostFunction(1, 4, junk_jacobian.data()), NULL, block4);
260 problem.AddResidualBlock(
261 new UnaryCostFunction(1, 3, junk_jacobian.data()), NULL, block3);
262 problem.AddResidualBlock(
263 new UnaryCostFunction(1, 2, junk_jacobian.data()), NULL, block2);
264 problem.SetParameterBlockConstant(block3);
265
266 // Sparsity pattern
267 //
268 // Note that the problem structure does not imply this sparsity
269 // pattern since all the residual blocks are unary. But the
270 // ComputeCovarianceSparsity function in its current incarnation
271 // does not pay attention to this fact and only looks at the
272 // parameter block pairs that the user provides.
273 //
274 // X . . X X X X
275 // . X X . . . .
276 // . X X . . . .
277 // . . . X X X X
278 // . . . X X X X
279 // . . . X X X X
280 // . . . X X X X
281
282 int expected_rows[] = {0, 5, 7, 9, 13, 17, 21, 25};
283 int expected_cols[] = {0, 3, 4, 5, 6,
284 1, 2,
285 1, 2,
286 3, 4, 5, 6,
287 3, 4, 5, 6,
288 3, 4, 5, 6,
289 3, 4, 5, 6};
290
291 vector<pair<const double*, const double*>> covariance_blocks;
292 covariance_blocks.push_back(make_pair(block1, block1));
293 covariance_blocks.push_back(make_pair(block4, block4));
294 covariance_blocks.push_back(make_pair(block2, block2));
295 covariance_blocks.push_back(make_pair(block3, block3));
296 covariance_blocks.push_back(make_pair(block2, block3));
297 covariance_blocks.push_back(make_pair(block4, block1)); // reversed
298
299 Covariance::Options options;
300 CovarianceImpl covariance_impl(options);
301 EXPECT_TRUE(covariance_impl
302 .ComputeCovarianceSparsity(covariance_blocks, &problem));
303
304 const CompressedRowSparseMatrix* crsm = covariance_impl.covariance_matrix();
305
306 EXPECT_EQ(crsm->num_rows(), 7);
307 EXPECT_EQ(crsm->num_cols(), 7);
308 EXPECT_EQ(crsm->num_nonzeros(), 25);
309
310 const int* rows = crsm->rows();
311 for (int r = 0; r < crsm->num_rows() + 1; ++r) {
312 EXPECT_EQ(rows[r], expected_rows[r])
313 << r << " "
314 << rows[r] << " "
315 << expected_rows[r];
316 }
317
318 const int* cols = crsm->cols();
319 for (int c = 0; c < crsm->num_nonzeros(); ++c) {
320 EXPECT_EQ(cols[c], expected_cols[c])
321 << c << " "
322 << cols[c] << " "
323 << expected_cols[c];
324 }
325}
326
327TEST(CovarianceImpl, ComputeCovarianceSparsityWithFreeParameterBlock) {
328 double parameters[10];
329
330 double* block1 = parameters;
331 double* block2 = block1 + 1;
332 double* block3 = block2 + 2;
333 double* block4 = block3 + 3;
334
335 ProblemImpl problem;
336
337 // Add in random order
338 Vector junk_jacobian = Vector::Zero(10);
339 problem.AddResidualBlock(
340 new UnaryCostFunction(1, 1, junk_jacobian.data()), NULL, block1);
341 problem.AddResidualBlock(
342 new UnaryCostFunction(1, 4, junk_jacobian.data()), NULL, block4);
343 problem.AddParameterBlock(block3, 3);
344 problem.AddResidualBlock(
345 new UnaryCostFunction(1, 2, junk_jacobian.data()), NULL, block2);
346
347 // Sparsity pattern
348 //
349 // Note that the problem structure does not imply this sparsity
350 // pattern since all the residual blocks are unary. But the
351 // ComputeCovarianceSparsity function in its current incarnation
352 // does not pay attention to this fact and only looks at the
353 // parameter block pairs that the user provides.
354 //
355 // X . . X X X X
356 // . X X . . . .
357 // . X X . . . .
358 // . . . X X X X
359 // . . . X X X X
360 // . . . X X X X
361 // . . . X X X X
362
363 int expected_rows[] = {0, 5, 7, 9, 13, 17, 21, 25};
364 int expected_cols[] = {0, 3, 4, 5, 6,
365 1, 2,
366 1, 2,
367 3, 4, 5, 6,
368 3, 4, 5, 6,
369 3, 4, 5, 6,
370 3, 4, 5, 6};
371
372 vector<pair<const double*, const double*>> covariance_blocks;
373 covariance_blocks.push_back(make_pair(block1, block1));
374 covariance_blocks.push_back(make_pair(block4, block4));
375 covariance_blocks.push_back(make_pair(block2, block2));
376 covariance_blocks.push_back(make_pair(block3, block3));
377 covariance_blocks.push_back(make_pair(block2, block3));
378 covariance_blocks.push_back(make_pair(block4, block1)); // reversed
379
380 Covariance::Options options;
381 CovarianceImpl covariance_impl(options);
382 EXPECT_TRUE(covariance_impl
383 .ComputeCovarianceSparsity(covariance_blocks, &problem));
384
385 const CompressedRowSparseMatrix* crsm = covariance_impl.covariance_matrix();
386
387 EXPECT_EQ(crsm->num_rows(), 7);
388 EXPECT_EQ(crsm->num_cols(), 7);
389 EXPECT_EQ(crsm->num_nonzeros(), 25);
390
391 const int* rows = crsm->rows();
392 for (int r = 0; r < crsm->num_rows() + 1; ++r) {
393 EXPECT_EQ(rows[r], expected_rows[r])
394 << r << " "
395 << rows[r] << " "
396 << expected_rows[r];
397 }
398
399 const int* cols = crsm->cols();
400 for (int c = 0; c < crsm->num_nonzeros(); ++c) {
401 EXPECT_EQ(cols[c], expected_cols[c])
402 << c << " "
403 << cols[c] << " "
404 << expected_cols[c];
405 }
406}
407
408class CovarianceTest : public ::testing::Test {
409 protected:
410 typedef map<const double*, pair<int, int>> BoundsMap;
411
412 virtual void SetUp() {
413 double* x = parameters_;
414 double* y = x + 2;
415 double* z = y + 3;
416
417 x[0] = 1;
418 x[1] = 1;
419 y[0] = 2;
420 y[1] = 2;
421 y[2] = 2;
422 z[0] = 3;
423
424 {
425 double jacobian[] = { 1.0, 0.0, 0.0, 1.0};
426 problem_.AddResidualBlock(new UnaryCostFunction(2, 2, jacobian), NULL, x);
427 }
428
429 {
430 double jacobian[] = { 2.0, 0.0, 0.0, 0.0, 2.0, 0.0, 0.0, 0.0, 2.0 };
431 problem_.AddResidualBlock(new UnaryCostFunction(3, 3, jacobian), NULL, y);
432 }
433
434 {
435 double jacobian = 5.0;
436 problem_.AddResidualBlock(new UnaryCostFunction(1, 1, &jacobian),
437 NULL,
438 z);
439 }
440
441 {
442 double jacobian1[] = { 1.0, 2.0, 3.0 };
443 double jacobian2[] = { -5.0, -6.0 };
444 problem_.AddResidualBlock(
445 new BinaryCostFunction(1, 3, 2, jacobian1, jacobian2),
446 NULL,
447 y,
448 x);
449 }
450
451 {
452 double jacobian1[] = {2.0 };
453 double jacobian2[] = { 3.0, -2.0 };
454 problem_.AddResidualBlock(
455 new BinaryCostFunction(1, 1, 2, jacobian1, jacobian2),
456 NULL,
457 z,
458 x);
459 }
460
461 all_covariance_blocks_.push_back(make_pair(x, x));
462 all_covariance_blocks_.push_back(make_pair(y, y));
463 all_covariance_blocks_.push_back(make_pair(z, z));
464 all_covariance_blocks_.push_back(make_pair(x, y));
465 all_covariance_blocks_.push_back(make_pair(x, z));
466 all_covariance_blocks_.push_back(make_pair(y, z));
467
468 column_bounds_[x] = make_pair(0, 2);
469 column_bounds_[y] = make_pair(2, 5);
470 column_bounds_[z] = make_pair(5, 6);
471 }
472
473 // Computes covariance in ambient space.
474 void ComputeAndCompareCovarianceBlocks(const Covariance::Options& options,
475 const double* expected_covariance) {
476 ComputeAndCompareCovarianceBlocksInTangentOrAmbientSpace(
477 options,
478 true, // ambient
479 expected_covariance);
480 }
481
482 // Computes covariance in tangent space.
483 void ComputeAndCompareCovarianceBlocksInTangentSpace(
484 const Covariance::Options& options,
485 const double* expected_covariance) {
486 ComputeAndCompareCovarianceBlocksInTangentOrAmbientSpace(
487 options,
488 false, // tangent
489 expected_covariance);
490 }
491
492 void ComputeAndCompareCovarianceBlocksInTangentOrAmbientSpace(
493 const Covariance::Options& options,
494 bool lift_covariance_to_ambient_space,
495 const double* expected_covariance) {
496 // Generate all possible combination of block pairs and check if the
497 // covariance computation is correct.
498 for (int i = 0; i <= 64; ++i) {
499 vector<pair<const double*, const double*>> covariance_blocks;
500 if (i & 1) {
501 covariance_blocks.push_back(all_covariance_blocks_[0]);
502 }
503
504 if (i & 2) {
505 covariance_blocks.push_back(all_covariance_blocks_[1]);
506 }
507
508 if (i & 4) {
509 covariance_blocks.push_back(all_covariance_blocks_[2]);
510 }
511
512 if (i & 8) {
513 covariance_blocks.push_back(all_covariance_blocks_[3]);
514 }
515
516 if (i & 16) {
517 covariance_blocks.push_back(all_covariance_blocks_[4]);
518 }
519
520 if (i & 32) {
521 covariance_blocks.push_back(all_covariance_blocks_[5]);
522 }
523
524 Covariance covariance(options);
525 EXPECT_TRUE(covariance.Compute(covariance_blocks, &problem_));
526
527 for (int i = 0; i < covariance_blocks.size(); ++i) {
528 const double* block1 = covariance_blocks[i].first;
529 const double* block2 = covariance_blocks[i].second;
530 // block1, block2
531 GetCovarianceBlockAndCompare(block1,
532 block2,
533 lift_covariance_to_ambient_space,
534 covariance,
535 expected_covariance);
536 // block2, block1
537 GetCovarianceBlockAndCompare(block2,
538 block1,
539 lift_covariance_to_ambient_space,
540 covariance,
541 expected_covariance);
542 }
543 }
544 }
545
546 void GetCovarianceBlockAndCompare(const double* block1,
547 const double* block2,
548 bool lift_covariance_to_ambient_space,
549 const Covariance& covariance,
550 const double* expected_covariance) {
551 const BoundsMap& column_bounds = lift_covariance_to_ambient_space ?
552 column_bounds_ : local_column_bounds_;
553 const int row_begin = FindOrDie(column_bounds, block1).first;
554 const int row_end = FindOrDie(column_bounds, block1).second;
555 const int col_begin = FindOrDie(column_bounds, block2).first;
556 const int col_end = FindOrDie(column_bounds, block2).second;
557
558 Matrix actual(row_end - row_begin, col_end - col_begin);
559 if (lift_covariance_to_ambient_space) {
560 EXPECT_TRUE(covariance.GetCovarianceBlock(block1,
561 block2,
562 actual.data()));
563 } else {
564 EXPECT_TRUE(covariance.GetCovarianceBlockInTangentSpace(block1,
565 block2,
566 actual.data()));
567 }
568
569 int dof = 0; // degrees of freedom = sum of LocalSize()s
570 for (const auto& bound : column_bounds) {
571 dof = std::max(dof, bound.second.second);
572 }
573 ConstMatrixRef expected(expected_covariance, dof, dof);
574 double diff_norm = (expected.block(row_begin,
575 col_begin,
576 row_end - row_begin,
577 col_end - col_begin) - actual).norm();
578 diff_norm /= (row_end - row_begin) * (col_end - col_begin);
579
580 const double kTolerance = 1e-5;
581 EXPECT_NEAR(diff_norm, 0.0, kTolerance)
582 << "rows: " << row_begin << " " << row_end << " "
583 << "cols: " << col_begin << " " << col_end << " "
584 << "\n\n expected: \n " << expected.block(row_begin,
585 col_begin,
586 row_end - row_begin,
587 col_end - col_begin)
588 << "\n\n actual: \n " << actual
589 << "\n\n full expected: \n" << expected;
590 }
591
592 double parameters_[6];
593 Problem problem_;
594 vector<pair<const double*, const double*>> all_covariance_blocks_;
595 BoundsMap column_bounds_;
596 BoundsMap local_column_bounds_;
597};
598
599
600TEST_F(CovarianceTest, NormalBehavior) {
601 // J
602 //
603 // 1 0 0 0 0 0
604 // 0 1 0 0 0 0
605 // 0 0 2 0 0 0
606 // 0 0 0 2 0 0
607 // 0 0 0 0 2 0
608 // 0 0 0 0 0 5
609 // -5 -6 1 2 3 0
610 // 3 -2 0 0 0 2
611
612 // J'J
613 //
614 // 35 24 -5 -10 -15 6
615 // 24 41 -6 -12 -18 -4
616 // -5 -6 5 2 3 0
617 // -10 -12 2 8 6 0
618 // -15 -18 3 6 13 0
619 // 6 -4 0 0 0 29
620
621 // inv(J'J) computed using octave.
622 double expected_covariance[] = {
623 7.0747e-02, -8.4923e-03, 1.6821e-02, 3.3643e-02, 5.0464e-02, -1.5809e-02, // NOLINT
624 -8.4923e-03, 8.1352e-02, 2.4758e-02, 4.9517e-02, 7.4275e-02, 1.2978e-02, // NOLINT
625 1.6821e-02, 2.4758e-02, 2.4904e-01, -1.9271e-03, -2.8906e-03, -6.5325e-05, // NOLINT
626 3.3643e-02, 4.9517e-02, -1.9271e-03, 2.4615e-01, -5.7813e-03, -1.3065e-04, // NOLINT
627 5.0464e-02, 7.4275e-02, -2.8906e-03, -5.7813e-03, 2.4133e-01, -1.9598e-04, // NOLINT
628 -1.5809e-02, 1.2978e-02, -6.5325e-05, -1.3065e-04, -1.9598e-04, 3.9544e-02, // NOLINT
629 };
630
631 Covariance::Options options;
632
633#ifndef CERES_NO_SUITESPARSE
634 options.algorithm_type = SPARSE_QR;
635 options.sparse_linear_algebra_library_type = SUITE_SPARSE;
636 ComputeAndCompareCovarianceBlocks(options, expected_covariance);
637#endif
638
639 options.algorithm_type = DENSE_SVD;
640 ComputeAndCompareCovarianceBlocks(options, expected_covariance);
641
642 options.algorithm_type = SPARSE_QR;
643 options.sparse_linear_algebra_library_type = EIGEN_SPARSE;
644 ComputeAndCompareCovarianceBlocks(options, expected_covariance);
645}
646
647#ifdef CERES_USE_OPENMP
648
649TEST_F(CovarianceTest, ThreadedNormalBehavior) {
650 // J
651 //
652 // 1 0 0 0 0 0
653 // 0 1 0 0 0 0
654 // 0 0 2 0 0 0
655 // 0 0 0 2 0 0
656 // 0 0 0 0 2 0
657 // 0 0 0 0 0 5
658 // -5 -6 1 2 3 0
659 // 3 -2 0 0 0 2
660
661 // J'J
662 //
663 // 35 24 -5 -10 -15 6
664 // 24 41 -6 -12 -18 -4
665 // -5 -6 5 2 3 0
666 // -10 -12 2 8 6 0
667 // -15 -18 3 6 13 0
668 // 6 -4 0 0 0 29
669
670 // inv(J'J) computed using octave.
671 double expected_covariance[] = {
672 7.0747e-02, -8.4923e-03, 1.6821e-02, 3.3643e-02, 5.0464e-02, -1.5809e-02, // NOLINT
673 -8.4923e-03, 8.1352e-02, 2.4758e-02, 4.9517e-02, 7.4275e-02, 1.2978e-02, // NOLINT
674 1.6821e-02, 2.4758e-02, 2.4904e-01, -1.9271e-03, -2.8906e-03, -6.5325e-05, // NOLINT
675 3.3643e-02, 4.9517e-02, -1.9271e-03, 2.4615e-01, -5.7813e-03, -1.3065e-04, // NOLINT
676 5.0464e-02, 7.4275e-02, -2.8906e-03, -5.7813e-03, 2.4133e-01, -1.9598e-04, // NOLINT
677 -1.5809e-02, 1.2978e-02, -6.5325e-05, -1.3065e-04, -1.9598e-04, 3.9544e-02, // NOLINT
678 };
679
680 Covariance::Options options;
681 options.num_threads = 4;
682
683#ifndef CERES_NO_SUITESPARSE
684 options.algorithm_type = SPARSE_QR;
685 options.sparse_linear_algebra_library_type = SUITE_SPARSE;
686 ComputeAndCompareCovarianceBlocks(options, expected_covariance);
687#endif
688
689 options.algorithm_type = DENSE_SVD;
690 ComputeAndCompareCovarianceBlocks(options, expected_covariance);
691
692 options.algorithm_type = SPARSE_QR;
693 options.sparse_linear_algebra_library_type = EIGEN_SPARSE;
694 ComputeAndCompareCovarianceBlocks(options, expected_covariance);
695}
696
697#endif // CERES_USE_OPENMP
698
699TEST_F(CovarianceTest, ConstantParameterBlock) {
700 problem_.SetParameterBlockConstant(parameters_);
701
702 // J
703 //
704 // 0 0 0 0 0 0
705 // 0 0 0 0 0 0
706 // 0 0 2 0 0 0
707 // 0 0 0 2 0 0
708 // 0 0 0 0 2 0
709 // 0 0 0 0 0 5
710 // 0 0 1 2 3 0
711 // 0 0 0 0 0 2
712
713 // J'J
714 //
715 // 0 0 0 0 0 0
716 // 0 0 0 0 0 0
717 // 0 0 5 2 3 0
718 // 0 0 2 8 6 0
719 // 0 0 3 6 13 0
720 // 0 0 0 0 0 29
721
722 // pinv(J'J) computed using octave.
723 double expected_covariance[] = {
724 0, 0, 0, 0, 0, 0, // NOLINT
725 0, 0, 0, 0, 0, 0, // NOLINT
726 0, 0, 0.23611, -0.02778, -0.04167, -0.00000, // NOLINT
727 0, 0, -0.02778, 0.19444, -0.08333, -0.00000, // NOLINT
728 0, 0, -0.04167, -0.08333, 0.12500, -0.00000, // NOLINT
729 0, 0, -0.00000, -0.00000, -0.00000, 0.03448 // NOLINT
730 };
731
732 Covariance::Options options;
733
734#ifndef CERES_NO_SUITESPARSE
735 options.algorithm_type = SPARSE_QR;
736 options.sparse_linear_algebra_library_type = SUITE_SPARSE;
737 ComputeAndCompareCovarianceBlocks(options, expected_covariance);
738#endif
739
740 options.algorithm_type = DENSE_SVD;
741 ComputeAndCompareCovarianceBlocks(options, expected_covariance);
742
743 options.algorithm_type = SPARSE_QR;
744 options.sparse_linear_algebra_library_type = EIGEN_SPARSE;
745 ComputeAndCompareCovarianceBlocks(options, expected_covariance);
746}
747
748TEST_F(CovarianceTest, LocalParameterization) {
749 double* x = parameters_;
750 double* y = x + 2;
751
752 problem_.SetParameterization(x, new PolynomialParameterization);
753
754 vector<int> subset;
755 subset.push_back(2);
756 problem_.SetParameterization(y, new SubsetParameterization(3, subset));
757
758 // Raw Jacobian: J
759 //
760 // 1 0 0 0 0 0
761 // 0 1 0 0 0 0
762 // 0 0 2 0 0 0
763 // 0 0 0 2 0 0
764 // 0 0 0 0 2 0
765 // 0 0 0 0 0 5
766 // -5 -6 1 2 3 0
767 // 3 -2 0 0 0 2
768
769 // Local to global jacobian: A
770 //
771 // 1 0 0 0
772 // 1 0 0 0
773 // 0 1 0 0
774 // 0 0 1 0
775 // 0 0 0 0
776 // 0 0 0 1
777
778 // A * inv((J*A)'*(J*A)) * A'
779 // Computed using octave.
780 double expected_covariance[] = {
781 0.01766, 0.01766, 0.02158, 0.04316, 0.00000, -0.00122,
782 0.01766, 0.01766, 0.02158, 0.04316, 0.00000, -0.00122,
783 0.02158, 0.02158, 0.24860, -0.00281, 0.00000, -0.00149,
784 0.04316, 0.04316, -0.00281, 0.24439, 0.00000, -0.00298,
785 0.00000, 0.00000, 0.00000, 0.00000, 0.00000, 0.00000,
786 -0.00122, -0.00122, -0.00149, -0.00298, 0.00000, 0.03457
787 };
788
789 Covariance::Options options;
790
791#ifndef CERES_NO_SUITESPARSE
792 options.algorithm_type = SPARSE_QR;
793 options.sparse_linear_algebra_library_type = SUITE_SPARSE;
794 ComputeAndCompareCovarianceBlocks(options, expected_covariance);
795#endif
796
797 options.algorithm_type = DENSE_SVD;
798 ComputeAndCompareCovarianceBlocks(options, expected_covariance);
799
800 options.algorithm_type = SPARSE_QR;
801 options.sparse_linear_algebra_library_type = EIGEN_SPARSE;
802 ComputeAndCompareCovarianceBlocks(options, expected_covariance);
803}
804
805TEST_F(CovarianceTest, LocalParameterizationInTangentSpace) {
806 double* x = parameters_;
807 double* y = x + 2;
808 double* z = y + 3;
809
810 problem_.SetParameterization(x, new PolynomialParameterization);
811
812 vector<int> subset;
813 subset.push_back(2);
814 problem_.SetParameterization(y, new SubsetParameterization(3, subset));
815
816 local_column_bounds_[x] = make_pair(0, 1);
817 local_column_bounds_[y] = make_pair(1, 3);
818 local_column_bounds_[z] = make_pair(3, 4);
819
820 // Raw Jacobian: J
821 //
822 // 1 0 0 0 0 0
823 // 0 1 0 0 0 0
824 // 0 0 2 0 0 0
825 // 0 0 0 2 0 0
826 // 0 0 0 0 2 0
827 // 0 0 0 0 0 5
828 // -5 -6 1 2 3 0
829 // 3 -2 0 0 0 2
830
831 // Local to global jacobian: A
832 //
833 // 1 0 0 0
834 // 1 0 0 0
835 // 0 1 0 0
836 // 0 0 1 0
837 // 0 0 0 0
838 // 0 0 0 1
839
840 // inv((J*A)'*(J*A))
841 // Computed using octave.
842 double expected_covariance[] = {
843 0.01766, 0.02158, 0.04316, -0.00122,
844 0.02158, 0.24860, -0.00281, -0.00149,
845 0.04316, -0.00281, 0.24439, -0.00298,
846 -0.00122, -0.00149, -0.00298, 0.03457 // NOLINT
847 };
848
849 Covariance::Options options;
850
851#ifndef CERES_NO_SUITESPARSE
852 options.algorithm_type = SPARSE_QR;
853 options.sparse_linear_algebra_library_type = SUITE_SPARSE;
854
855 ComputeAndCompareCovarianceBlocksInTangentSpace(options, expected_covariance);
856#endif
857
858 options.algorithm_type = DENSE_SVD;
859 ComputeAndCompareCovarianceBlocksInTangentSpace(options, expected_covariance);
860
861 options.algorithm_type = SPARSE_QR;
862 options.sparse_linear_algebra_library_type = EIGEN_SPARSE;
863 ComputeAndCompareCovarianceBlocksInTangentSpace(options, expected_covariance);
864}
865
866TEST_F(CovarianceTest, LocalParameterizationInTangentSpaceWithConstantBlocks) {
867 double* x = parameters_;
868 double* y = x + 2;
869 double* z = y + 3;
870
871 problem_.SetParameterization(x, new PolynomialParameterization);
872 problem_.SetParameterBlockConstant(x);
873
874 vector<int> subset;
875 subset.push_back(2);
876 problem_.SetParameterization(y, new SubsetParameterization(3, subset));
877 problem_.SetParameterBlockConstant(y);
878
879 local_column_bounds_[x] = make_pair(0, 1);
880 local_column_bounds_[y] = make_pair(1, 3);
881 local_column_bounds_[z] = make_pair(3, 4);
882
883 // Raw Jacobian: J
884 //
885 // 1 0 0 0 0 0
886 // 0 1 0 0 0 0
887 // 0 0 2 0 0 0
888 // 0 0 0 2 0 0
889 // 0 0 0 0 2 0
890 // 0 0 0 0 0 5
891 // -5 -6 1 2 3 0
892 // 3 -2 0 0 0 2
893
894 // Local to global jacobian: A
895 //
896 // 0 0 0 0
897 // 0 0 0 0
898 // 0 0 0 0
899 // 0 0 0 0
900 // 0 0 0 0
901 // 0 0 0 1
902
903 // pinv((J*A)'*(J*A))
904 // Computed using octave.
905 double expected_covariance[] = {
906 0.0, 0.0, 0.0, 0.0,
907 0.0, 0.0, 0.0, 0.0,
908 0.0, 0.0, 0.0, 0.0,
909 0.0, 0.0, 0.0, 0.034482 // NOLINT
910 };
911
912 Covariance::Options options;
913
914#ifndef CERES_NO_SUITESPARSE
915 options.algorithm_type = SPARSE_QR;
916 options.sparse_linear_algebra_library_type = SUITE_SPARSE;
917
918 ComputeAndCompareCovarianceBlocksInTangentSpace(options, expected_covariance);
919#endif
920
921 options.algorithm_type = DENSE_SVD;
922 ComputeAndCompareCovarianceBlocksInTangentSpace(options, expected_covariance);
923
924 options.algorithm_type = SPARSE_QR;
925 options.sparse_linear_algebra_library_type = EIGEN_SPARSE;
926 ComputeAndCompareCovarianceBlocksInTangentSpace(options, expected_covariance);
927}
928
929TEST_F(CovarianceTest, TruncatedRank) {
930 // J
931 //
932 // 1 0 0 0 0 0
933 // 0 1 0 0 0 0
934 // 0 0 2 0 0 0
935 // 0 0 0 2 0 0
936 // 0 0 0 0 2 0
937 // 0 0 0 0 0 5
938 // -5 -6 1 2 3 0
939 // 3 -2 0 0 0 2
940
941 // J'J
942 //
943 // 35 24 -5 -10 -15 6
944 // 24 41 -6 -12 -18 -4
945 // -5 -6 5 2 3 0
946 // -10 -12 2 8 6 0
947 // -15 -18 3 6 13 0
948 // 6 -4 0 0 0 29
949
950 // 3.4142 is the smallest eigen value of J'J. The following matrix
951 // was obtained by dropping the eigenvector corresponding to this
952 // eigenvalue.
953 double expected_covariance[] = {
954 5.4135e-02, -3.5121e-02, 1.7257e-04, 3.4514e-04, 5.1771e-04, -1.6076e-02, // NOLINT
955 -3.5121e-02, 3.8667e-02, -1.9288e-03, -3.8576e-03, -5.7864e-03, 1.2549e-02, // NOLINT
956 1.7257e-04, -1.9288e-03, 2.3235e-01, -3.5297e-02, -5.2946e-02, -3.3329e-04, // NOLINT
957 3.4514e-04, -3.8576e-03, -3.5297e-02, 1.7941e-01, -1.0589e-01, -6.6659e-04, // NOLINT
958 5.1771e-04, -5.7864e-03, -5.2946e-02, -1.0589e-01, 9.1162e-02, -9.9988e-04, // NOLINT
959 -1.6076e-02, 1.2549e-02, -3.3329e-04, -6.6659e-04, -9.9988e-04, 3.9539e-02 // NOLINT
960 };
961
962
963 {
964 Covariance::Options options;
965 options.algorithm_type = DENSE_SVD;
966 // Force dropping of the smallest eigenvector.
967 options.null_space_rank = 1;
968 ComputeAndCompareCovarianceBlocks(options, expected_covariance);
969 }
970
971 {
972 Covariance::Options options;
973 options.algorithm_type = DENSE_SVD;
974 // Force dropping of the smallest eigenvector via the ratio but
975 // automatic truncation.
976 options.min_reciprocal_condition_number = 0.044494;
977 options.null_space_rank = -1;
978 ComputeAndCompareCovarianceBlocks(options, expected_covariance);
979 }
980}
981
982TEST_F(CovarianceTest, DenseCovarianceMatrixFromSetOfParameters) {
983 Covariance::Options options;
984 Covariance covariance(options);
985 double* x = parameters_;
986 double* y = x + 2;
987 double* z = y + 3;
988 vector<const double*> parameter_blocks;
989 parameter_blocks.push_back(x);
990 parameter_blocks.push_back(y);
991 parameter_blocks.push_back(z);
992 covariance.Compute(parameter_blocks, &problem_);
993 double expected_covariance[36];
994 covariance.GetCovarianceMatrix(parameter_blocks, expected_covariance);
995
996#ifndef CERES_NO_SUITESPARSE
997 options.algorithm_type = SPARSE_QR;
998 options.sparse_linear_algebra_library_type = SUITE_SPARSE;
999 ComputeAndCompareCovarianceBlocks(options, expected_covariance);
1000#endif
1001
1002 options.algorithm_type = DENSE_SVD;
1003 ComputeAndCompareCovarianceBlocks(options, expected_covariance);
1004
1005 options.algorithm_type = SPARSE_QR;
1006 options.sparse_linear_algebra_library_type = EIGEN_SPARSE;
1007 ComputeAndCompareCovarianceBlocks(options, expected_covariance);
1008}
1009
1010TEST_F(CovarianceTest, DenseCovarianceMatrixFromSetOfParametersThreaded) {
1011 Covariance::Options options;
1012 options.num_threads = 4;
1013 Covariance covariance(options);
1014 double* x = parameters_;
1015 double* y = x + 2;
1016 double* z = y + 3;
1017 vector<const double*> parameter_blocks;
1018 parameter_blocks.push_back(x);
1019 parameter_blocks.push_back(y);
1020 parameter_blocks.push_back(z);
1021 covariance.Compute(parameter_blocks, &problem_);
1022 double expected_covariance[36];
1023 covariance.GetCovarianceMatrix(parameter_blocks, expected_covariance);
1024
1025#ifndef CERES_NO_SUITESPARSE
1026 options.algorithm_type = SPARSE_QR;
1027 options.sparse_linear_algebra_library_type = SUITE_SPARSE;
1028 ComputeAndCompareCovarianceBlocks(options, expected_covariance);
1029#endif
1030
1031 options.algorithm_type = DENSE_SVD;
1032 ComputeAndCompareCovarianceBlocks(options, expected_covariance);
1033
1034 options.algorithm_type = SPARSE_QR;
1035 options.sparse_linear_algebra_library_type = EIGEN_SPARSE;
1036 ComputeAndCompareCovarianceBlocks(options, expected_covariance);
1037}
1038
1039TEST_F(CovarianceTest, DenseCovarianceMatrixFromSetOfParametersInTangentSpace) {
1040 Covariance::Options options;
1041 Covariance covariance(options);
1042 double* x = parameters_;
1043 double* y = x + 2;
1044 double* z = y + 3;
1045
1046 problem_.SetParameterization(x, new PolynomialParameterization);
1047
1048 vector<int> subset;
1049 subset.push_back(2);
1050 problem_.SetParameterization(y, new SubsetParameterization(3, subset));
1051
1052 local_column_bounds_[x] = make_pair(0, 1);
1053 local_column_bounds_[y] = make_pair(1, 3);
1054 local_column_bounds_[z] = make_pair(3, 4);
1055
1056 vector<const double*> parameter_blocks;
1057 parameter_blocks.push_back(x);
1058 parameter_blocks.push_back(y);
1059 parameter_blocks.push_back(z);
1060 covariance.Compute(parameter_blocks, &problem_);
1061 double expected_covariance[16];
1062 covariance.GetCovarianceMatrixInTangentSpace(parameter_blocks,
1063 expected_covariance);
1064
1065#ifndef CERES_NO_SUITESPARSE
1066 options.algorithm_type = SPARSE_QR;
1067 options.sparse_linear_algebra_library_type = SUITE_SPARSE;
1068
1069 ComputeAndCompareCovarianceBlocksInTangentSpace(options, expected_covariance);
1070#endif
1071
1072 options.algorithm_type = DENSE_SVD;
1073 ComputeAndCompareCovarianceBlocksInTangentSpace(options, expected_covariance);
1074
1075 options.algorithm_type = SPARSE_QR;
1076 options.sparse_linear_algebra_library_type = EIGEN_SPARSE;
1077 ComputeAndCompareCovarianceBlocksInTangentSpace(options, expected_covariance);
1078}
1079
1080TEST_F(CovarianceTest, ComputeCovarianceFailure) {
1081 Covariance::Options options;
1082 Covariance covariance(options);
1083 double* x = parameters_;
1084 double* y = x + 2;
1085 vector<const double*> parameter_blocks;
1086 parameter_blocks.push_back(x);
1087 parameter_blocks.push_back(x);
1088 parameter_blocks.push_back(y);
1089 parameter_blocks.push_back(y);
1090 EXPECT_DEATH_IF_SUPPORTED(covariance.Compute(parameter_blocks, &problem_),
1091 "Covariance::Compute called with duplicate blocks "
1092 "at indices \\(0, 1\\) and \\(2, 3\\)");
1093 vector<pair<const double*, const double*>> covariance_blocks;
1094 covariance_blocks.push_back(make_pair(x, x));
1095 covariance_blocks.push_back(make_pair(x, x));
1096 covariance_blocks.push_back(make_pair(y, y));
1097 covariance_blocks.push_back(make_pair(y, y));
1098 EXPECT_DEATH_IF_SUPPORTED(covariance.Compute(covariance_blocks, &problem_),
1099 "Covariance::Compute called with duplicate blocks "
1100 "at indices \\(0, 1\\) and \\(2, 3\\)");
1101}
1102
1103class RankDeficientCovarianceTest : public CovarianceTest {
1104 protected:
1105 virtual void SetUp() {
1106 double* x = parameters_;
1107 double* y = x + 2;
1108 double* z = y + 3;
1109
1110 {
1111 double jacobian[] = { 1.0, 0.0, 0.0, 1.0};
1112 problem_.AddResidualBlock(new UnaryCostFunction(2, 2, jacobian), NULL, x);
1113 }
1114
1115 {
1116 double jacobian[] = { 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0 };
1117 problem_.AddResidualBlock(new UnaryCostFunction(3, 3, jacobian), NULL, y);
1118 }
1119
1120 {
1121 double jacobian = 5.0;
1122 problem_.AddResidualBlock(new UnaryCostFunction(1, 1, &jacobian),
1123 NULL,
1124 z);
1125 }
1126
1127 {
1128 double jacobian1[] = { 0.0, 0.0, 0.0 };
1129 double jacobian2[] = { -5.0, -6.0 };
1130 problem_.AddResidualBlock(
1131 new BinaryCostFunction(1, 3, 2, jacobian1, jacobian2),
1132 NULL,
1133 y,
1134 x);
1135 }
1136
1137 {
1138 double jacobian1[] = {2.0 };
1139 double jacobian2[] = { 3.0, -2.0 };
1140 problem_.AddResidualBlock(
1141 new BinaryCostFunction(1, 1, 2, jacobian1, jacobian2),
1142 NULL,
1143 z,
1144 x);
1145 }
1146
1147 all_covariance_blocks_.push_back(make_pair(x, x));
1148 all_covariance_blocks_.push_back(make_pair(y, y));
1149 all_covariance_blocks_.push_back(make_pair(z, z));
1150 all_covariance_blocks_.push_back(make_pair(x, y));
1151 all_covariance_blocks_.push_back(make_pair(x, z));
1152 all_covariance_blocks_.push_back(make_pair(y, z));
1153
1154 column_bounds_[x] = make_pair(0, 2);
1155 column_bounds_[y] = make_pair(2, 5);
1156 column_bounds_[z] = make_pair(5, 6);
1157 }
1158};
1159
1160TEST_F(RankDeficientCovarianceTest, AutomaticTruncation) {
1161 // J
1162 //
1163 // 1 0 0 0 0 0
1164 // 0 1 0 0 0 0
1165 // 0 0 0 0 0 0
1166 // 0 0 0 0 0 0
1167 // 0 0 0 0 0 0
1168 // 0 0 0 0 0 5
1169 // -5 -6 0 0 0 0
1170 // 3 -2 0 0 0 2
1171
1172 // J'J
1173 //
1174 // 35 24 0 0 0 6
1175 // 24 41 0 0 0 -4
1176 // 0 0 0 0 0 0
1177 // 0 0 0 0 0 0
1178 // 0 0 0 0 0 0
1179 // 6 -4 0 0 0 29
1180
1181 // pinv(J'J) computed using octave.
1182 double expected_covariance[] = {
1183 0.053998, -0.033145, 0.000000, 0.000000, 0.000000, -0.015744,
1184 -0.033145, 0.045067, 0.000000, 0.000000, 0.000000, 0.013074,
1185 0.000000, 0.000000, 0.000000, 0.000000, 0.000000, 0.000000,
1186 0.000000, 0.000000, 0.000000, 0.000000, 0.000000, 0.000000,
1187 0.000000, 0.000000, 0.000000, 0.000000, 0.000000, 0.000000,
1188 -0.015744, 0.013074, 0.000000, 0.000000, 0.000000, 0.039543
1189 };
1190
1191 Covariance::Options options;
1192 options.algorithm_type = DENSE_SVD;
1193 options.null_space_rank = -1;
1194 ComputeAndCompareCovarianceBlocks(options, expected_covariance);
1195}
1196
1197class LargeScaleCovarianceTest : public ::testing::Test {
1198 protected:
1199 virtual void SetUp() {
1200 num_parameter_blocks_ = 2000;
1201 parameter_block_size_ = 5;
1202 parameters_.reset(
1203 new double[parameter_block_size_ * num_parameter_blocks_]);
1204
1205 Matrix jacobian(parameter_block_size_, parameter_block_size_);
1206 for (int i = 0; i < num_parameter_blocks_; ++i) {
1207 jacobian.setIdentity();
1208 jacobian *= (i + 1);
1209
1210 double* block_i = parameters_.get() + i * parameter_block_size_;
1211 problem_.AddResidualBlock(new UnaryCostFunction(parameter_block_size_,
1212 parameter_block_size_,
1213 jacobian.data()),
1214 NULL,
1215 block_i);
1216 for (int j = i; j < num_parameter_blocks_; ++j) {
1217 double* block_j = parameters_.get() + j * parameter_block_size_;
1218 all_covariance_blocks_.push_back(make_pair(block_i, block_j));
1219 }
1220 }
1221 }
1222
1223 void ComputeAndCompare(
1224 CovarianceAlgorithmType algorithm_type,
1225 SparseLinearAlgebraLibraryType sparse_linear_algebra_library_type,
1226 int num_threads) {
1227 Covariance::Options options;
1228 options.algorithm_type = algorithm_type;
1229 options.sparse_linear_algebra_library_type =
1230 sparse_linear_algebra_library_type;
1231 options.num_threads = num_threads;
1232 Covariance covariance(options);
1233 EXPECT_TRUE(covariance.Compute(all_covariance_blocks_, &problem_));
1234
1235 Matrix expected(parameter_block_size_, parameter_block_size_);
1236 Matrix actual(parameter_block_size_, parameter_block_size_);
1237 const double kTolerance = 1e-16;
1238
1239 for (int i = 0; i < num_parameter_blocks_; ++i) {
1240 expected.setIdentity();
1241 expected /= (i + 1.0) * (i + 1.0);
1242
1243 double* block_i = parameters_.get() + i * parameter_block_size_;
1244 covariance.GetCovarianceBlock(block_i, block_i, actual.data());
1245 EXPECT_NEAR((expected - actual).norm(), 0.0, kTolerance)
1246 << "block: " << i << ", " << i << "\n"
1247 << "expected: \n" << expected << "\n"
1248 << "actual: \n" << actual;
1249
1250 expected.setZero();
1251 for (int j = i + 1; j < num_parameter_blocks_; ++j) {
1252 double* block_j = parameters_.get() + j * parameter_block_size_;
1253 covariance.GetCovarianceBlock(block_i, block_j, actual.data());
1254 EXPECT_NEAR((expected - actual).norm(), 0.0, kTolerance)
1255 << "block: " << i << ", " << j << "\n"
1256 << "expected: \n" << expected << "\n"
1257 << "actual: \n" << actual;
1258 }
1259 }
1260 }
1261
1262 std::unique_ptr<double[]> parameters_;
1263 int parameter_block_size_;
1264 int num_parameter_blocks_;
1265
1266 Problem problem_;
1267 vector<pair<const double*, const double*>> all_covariance_blocks_;
1268};
1269
1270#if !defined(CERES_NO_SUITESPARSE) && defined(CERES_USE_OPENMP)
1271
1272TEST_F(LargeScaleCovarianceTest, Parallel) {
1273 ComputeAndCompare(SPARSE_QR, SUITE_SPARSE, 4);
1274}
1275
1276#endif // !defined(CERES_NO_SUITESPARSE) && defined(CERES_USE_OPENMP)
1277
1278} // namespace internal
1279} // namespace ceres