evaluator_test.cc 24 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726
  1. // Ceres Solver - A fast non-linear least squares minimizer
  2. // Copyright 2023 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. #include "ceres/evaluator.h"
  34. #include <memory>
  35. #include <string>
  36. #include <vector>
  37. #include "ceres/casts.h"
  38. #include "ceres/cost_function.h"
  39. #include "ceres/crs_matrix.h"
  40. #include "ceres/evaluator_test_utils.h"
  41. #include "ceres/internal/eigen.h"
  42. #include "ceres/manifold.h"
  43. #include "ceres/problem_impl.h"
  44. #include "ceres/program.h"
  45. #include "ceres/sized_cost_function.h"
  46. #include "ceres/sparse_matrix.h"
  47. #include "ceres/stringprintf.h"
  48. #include "ceres/types.h"
  49. #include "gtest/gtest.h"
  50. namespace ceres {
  51. namespace internal {
  52. // TODO(keir): Consider pushing this into a common test utils file.
  53. template <int kFactor, int kNumResiduals, int... Ns>
  54. class ParameterIgnoringCostFunction
  55. : public SizedCostFunction<kNumResiduals, Ns...> {
  56. using Base = SizedCostFunction<kNumResiduals, Ns...>;
  57. public:
  58. explicit ParameterIgnoringCostFunction(bool succeeds = true)
  59. : succeeds_(succeeds) {}
  60. bool Evaluate(double const* const* parameters,
  61. double* residuals,
  62. double** jacobians) const final {
  63. for (int i = 0; i < Base::num_residuals(); ++i) {
  64. residuals[i] = i + 1;
  65. }
  66. if (jacobians) {
  67. for (int k = 0; k < Base::parameter_block_sizes().size(); ++k) {
  68. // The jacobians here are full sized, but they are transformed in the
  69. // evaluator into the "local" jacobian. In the tests, the "subset
  70. // constant" manifold is used, which should pick out columns from these
  71. // jacobians. Put values in the jacobian that make this obvious; in
  72. // particular, make the jacobians like this:
  73. //
  74. // 1 2 3 4 ...
  75. // 1 2 3 4 ... .* kFactor
  76. // 1 2 3 4 ...
  77. //
  78. // where the multiplication by kFactor makes it easier to distinguish
  79. // between Jacobians of different residuals for the same parameter.
  80. if (jacobians[k] != nullptr) {
  81. MatrixRef jacobian(jacobians[k],
  82. Base::num_residuals(),
  83. Base::parameter_block_sizes()[k]);
  84. for (int j = 0; j < Base::parameter_block_sizes()[k]; ++j) {
  85. jacobian.col(j).setConstant(kFactor * (j + 1));
  86. }
  87. }
  88. }
  89. }
  90. return succeeds_;
  91. }
  92. private:
  93. bool succeeds_;
  94. };
  95. struct EvaluatorTestOptions {
  96. EvaluatorTestOptions(LinearSolverType linear_solver_type,
  97. int num_eliminate_blocks,
  98. bool dynamic_sparsity = false)
  99. : linear_solver_type(linear_solver_type),
  100. num_eliminate_blocks(num_eliminate_blocks),
  101. dynamic_sparsity(dynamic_sparsity) {}
  102. LinearSolverType linear_solver_type;
  103. int num_eliminate_blocks;
  104. bool dynamic_sparsity;
  105. };
  106. struct EvaluatorTest : public ::testing::TestWithParam<EvaluatorTestOptions> {
  107. std::unique_ptr<Evaluator> CreateEvaluator(Program* program) {
  108. // This program is straight from the ProblemImpl, and so has no index/offset
  109. // yet; compute it here as required by the evaluator implementations.
  110. program->SetParameterOffsetsAndIndex();
  111. if (VLOG_IS_ON(1)) {
  112. std::string report;
  113. StringAppendF(&report,
  114. "Creating evaluator with type: %d",
  115. GetParam().linear_solver_type);
  116. if (GetParam().linear_solver_type == SPARSE_NORMAL_CHOLESKY) {
  117. StringAppendF(
  118. &report, ", dynamic_sparsity: %d", GetParam().dynamic_sparsity);
  119. }
  120. StringAppendF(&report,
  121. " and num_eliminate_blocks: %d",
  122. GetParam().num_eliminate_blocks);
  123. VLOG(1) << report;
  124. }
  125. Evaluator::Options options;
  126. options.linear_solver_type = GetParam().linear_solver_type;
  127. options.num_eliminate_blocks = GetParam().num_eliminate_blocks;
  128. options.dynamic_sparsity = GetParam().dynamic_sparsity;
  129. options.context = problem.context();
  130. std::string error;
  131. return Evaluator::Create(options, program, &error);
  132. }
  133. void EvaluateAndCompare(ProblemImpl* problem,
  134. int expected_num_rows,
  135. int expected_num_cols,
  136. double expected_cost,
  137. const double* expected_residuals,
  138. const double* expected_gradient,
  139. const double* expected_jacobian) {
  140. std::unique_ptr<Evaluator> evaluator =
  141. CreateEvaluator(problem->mutable_program());
  142. int num_residuals = expected_num_rows;
  143. int num_parameters = expected_num_cols;
  144. double cost = -1;
  145. Vector residuals(num_residuals);
  146. residuals.setConstant(-2000);
  147. Vector gradient(num_parameters);
  148. gradient.setConstant(-3000);
  149. std::unique_ptr<SparseMatrix> jacobian(evaluator->CreateJacobian());
  150. ASSERT_EQ(expected_num_rows, evaluator->NumResiduals());
  151. ASSERT_EQ(expected_num_cols, evaluator->NumEffectiveParameters());
  152. ASSERT_EQ(expected_num_rows, jacobian->num_rows());
  153. ASSERT_EQ(expected_num_cols, jacobian->num_cols());
  154. std::vector<double> state(evaluator->NumParameters());
  155. // clang-format off
  156. ASSERT_TRUE(evaluator->Evaluate(
  157. &state[0],
  158. &cost,
  159. expected_residuals != nullptr ? &residuals[0] : nullptr,
  160. expected_gradient != nullptr ? &gradient[0] : nullptr,
  161. expected_jacobian != nullptr ? jacobian.get() : nullptr));
  162. // clang-format on
  163. Matrix actual_jacobian;
  164. if (expected_jacobian != nullptr) {
  165. jacobian->ToDenseMatrix(&actual_jacobian);
  166. }
  167. CompareEvaluations(expected_num_rows,
  168. expected_num_cols,
  169. expected_cost,
  170. expected_residuals,
  171. expected_gradient,
  172. expected_jacobian,
  173. cost,
  174. &residuals[0],
  175. &gradient[0],
  176. actual_jacobian.data());
  177. }
  178. // Try all combinations of parameters for the evaluator.
  179. void CheckAllEvaluationCombinations(const ExpectedEvaluation& expected) {
  180. for (int i = 0; i < 8; ++i) {
  181. EvaluateAndCompare(&problem,
  182. expected.num_rows,
  183. expected.num_cols,
  184. expected.cost,
  185. (i & 1) ? expected.residuals : nullptr,
  186. (i & 2) ? expected.gradient : nullptr,
  187. (i & 4) ? expected.jacobian : nullptr);
  188. }
  189. }
  190. // The values are ignored completely by the cost function.
  191. double x[2];
  192. double y[3];
  193. double z[4];
  194. ProblemImpl problem;
  195. };
  196. static void SetSparseMatrixConstant(SparseMatrix* sparse_matrix, double value) {
  197. VectorRef(sparse_matrix->mutable_values(), sparse_matrix->num_nonzeros())
  198. .setConstant(value);
  199. }
  200. TEST_P(EvaluatorTest, SingleResidualProblem) {
  201. problem.AddResidualBlock(
  202. new ParameterIgnoringCostFunction<1, 3, 2, 3, 4>, nullptr, x, y, z);
  203. // clang-format off
  204. ExpectedEvaluation expected = {
  205. // Rows/columns
  206. 3, 9,
  207. // Cost
  208. 7.0,
  209. // Residuals
  210. { 1.0, 2.0, 3.0 },
  211. // Gradient
  212. { 6.0, 12.0, // x
  213. 6.0, 12.0, 18.0, // y
  214. 6.0, 12.0, 18.0, 24.0, // z
  215. },
  216. // Jacobian
  217. // x y z
  218. { 1, 2, 1, 2, 3, 1, 2, 3, 4,
  219. 1, 2, 1, 2, 3, 1, 2, 3, 4,
  220. 1, 2, 1, 2, 3, 1, 2, 3, 4
  221. }
  222. };
  223. // clang-format on
  224. CheckAllEvaluationCombinations(expected);
  225. }
  226. TEST_P(EvaluatorTest, SingleResidualProblemWithPermutedParameters) {
  227. // Add the parameters in explicit order to force the ordering in the program.
  228. problem.AddParameterBlock(x, 2);
  229. problem.AddParameterBlock(y, 3);
  230. problem.AddParameterBlock(z, 4);
  231. // Then use a cost function which is similar to the others, but swap around
  232. // the ordering of the parameters to the cost function. This shouldn't affect
  233. // the jacobian evaluation, but requires explicit handling in the evaluators.
  234. // At one point the compressed row evaluator had a bug that went undetected
  235. // for a long time, since by chance most users added parameters to the problem
  236. // in the same order that they occurred as parameters to a cost function.
  237. problem.AddResidualBlock(
  238. new ParameterIgnoringCostFunction<1, 3, 4, 3, 2>, nullptr, z, y, x);
  239. // clang-format off
  240. ExpectedEvaluation expected = {
  241. // Rows/columns
  242. 3, 9,
  243. // Cost
  244. 7.0,
  245. // Residuals
  246. { 1.0, 2.0, 3.0 },
  247. // Gradient
  248. { 6.0, 12.0, // x
  249. 6.0, 12.0, 18.0, // y
  250. 6.0, 12.0, 18.0, 24.0, // z
  251. },
  252. // Jacobian
  253. // x y z
  254. { 1, 2, 1, 2, 3, 1, 2, 3, 4,
  255. 1, 2, 1, 2, 3, 1, 2, 3, 4,
  256. 1, 2, 1, 2, 3, 1, 2, 3, 4
  257. }
  258. };
  259. // clang-format on
  260. CheckAllEvaluationCombinations(expected);
  261. }
  262. TEST_P(EvaluatorTest, SingleResidualProblemWithNuisanceParameters) {
  263. // These parameters are not used.
  264. double a[2];
  265. double b[1];
  266. double c[1];
  267. double d[3];
  268. // Add the parameters in a mixed order so the Jacobian is "checkered" with the
  269. // values from the other parameters.
  270. problem.AddParameterBlock(a, 2);
  271. problem.AddParameterBlock(x, 2);
  272. problem.AddParameterBlock(b, 1);
  273. problem.AddParameterBlock(y, 3);
  274. problem.AddParameterBlock(c, 1);
  275. problem.AddParameterBlock(z, 4);
  276. problem.AddParameterBlock(d, 3);
  277. problem.AddResidualBlock(
  278. new ParameterIgnoringCostFunction<1, 3, 2, 3, 4>, nullptr, x, y, z);
  279. // clang-format off
  280. ExpectedEvaluation expected = {
  281. // Rows/columns
  282. 3, 16,
  283. // Cost
  284. 7.0,
  285. // Residuals
  286. { 1.0, 2.0, 3.0 },
  287. // Gradient
  288. { 0.0, 0.0, // a
  289. 6.0, 12.0, // x
  290. 0.0, // b
  291. 6.0, 12.0, 18.0, // y
  292. 0.0, // c
  293. 6.0, 12.0, 18.0, 24.0, // z
  294. 0.0, 0.0, 0.0, // d
  295. },
  296. // Jacobian
  297. // a x b y c z d
  298. { 0, 0, 1, 2, 0, 1, 2, 3, 0, 1, 2, 3, 4, 0, 0, 0,
  299. 0, 0, 1, 2, 0, 1, 2, 3, 0, 1, 2, 3, 4, 0, 0, 0,
  300. 0, 0, 1, 2, 0, 1, 2, 3, 0, 1, 2, 3, 4, 0, 0, 0
  301. }
  302. };
  303. // clang-format on
  304. CheckAllEvaluationCombinations(expected);
  305. }
  306. TEST_P(EvaluatorTest, MultipleResidualProblem) {
  307. // Add the parameters in explicit order to force the ordering in the program.
  308. problem.AddParameterBlock(x, 2);
  309. problem.AddParameterBlock(y, 3);
  310. problem.AddParameterBlock(z, 4);
  311. // f(x, y) in R^2
  312. problem.AddResidualBlock(
  313. new ParameterIgnoringCostFunction<1, 2, 2, 3>, nullptr, x, y);
  314. // g(x, z) in R^3
  315. problem.AddResidualBlock(
  316. new ParameterIgnoringCostFunction<2, 3, 2, 4>, nullptr, x, z);
  317. // h(y, z) in R^4
  318. problem.AddResidualBlock(
  319. new ParameterIgnoringCostFunction<3, 4, 3, 4>, nullptr, y, z);
  320. // clang-format off
  321. ExpectedEvaluation expected = {
  322. // Rows/columns
  323. 9, 9,
  324. // Cost
  325. // f g h
  326. ( 1 + 4 + 1 + 4 + 9 + 1 + 4 + 9 + 16) / 2.0,
  327. // Residuals
  328. { 1.0, 2.0, // f
  329. 1.0, 2.0, 3.0, // g
  330. 1.0, 2.0, 3.0, 4.0 // h
  331. },
  332. // Gradient
  333. { 15.0, 30.0, // x
  334. 33.0, 66.0, 99.0, // y
  335. 42.0, 84.0, 126.0, 168.0 // z
  336. },
  337. // Jacobian
  338. // x y z
  339. { /* f(x, y) */ 1, 2, 1, 2, 3, 0, 0, 0, 0,
  340. 1, 2, 1, 2, 3, 0, 0, 0, 0,
  341. /* g(x, z) */ 2, 4, 0, 0, 0, 2, 4, 6, 8,
  342. 2, 4, 0, 0, 0, 2, 4, 6, 8,
  343. 2, 4, 0, 0, 0, 2, 4, 6, 8,
  344. /* h(y, z) */ 0, 0, 3, 6, 9, 3, 6, 9, 12,
  345. 0, 0, 3, 6, 9, 3, 6, 9, 12,
  346. 0, 0, 3, 6, 9, 3, 6, 9, 12,
  347. 0, 0, 3, 6, 9, 3, 6, 9, 12
  348. }
  349. };
  350. // clang-format on
  351. CheckAllEvaluationCombinations(expected);
  352. }
  353. TEST_P(EvaluatorTest, MultipleResidualsWithManifolds) {
  354. // Add the parameters in explicit order to force the ordering in the program.
  355. problem.AddParameterBlock(x, 2);
  356. // Fix y's first dimension.
  357. std::vector<int> y_fixed;
  358. y_fixed.push_back(0);
  359. problem.AddParameterBlock(y, 3, new SubsetManifold(3, y_fixed));
  360. // Fix z's second dimension.
  361. std::vector<int> z_fixed;
  362. z_fixed.push_back(1);
  363. problem.AddParameterBlock(z, 4, new SubsetManifold(4, z_fixed));
  364. // f(x, y) in R^2
  365. problem.AddResidualBlock(
  366. new ParameterIgnoringCostFunction<1, 2, 2, 3>, nullptr, x, y);
  367. // g(x, z) in R^3
  368. problem.AddResidualBlock(
  369. new ParameterIgnoringCostFunction<2, 3, 2, 4>, nullptr, x, z);
  370. // h(y, z) in R^4
  371. problem.AddResidualBlock(
  372. new ParameterIgnoringCostFunction<3, 4, 3, 4>, nullptr, y, z);
  373. // clang-format off
  374. ExpectedEvaluation expected = {
  375. // Rows/columns
  376. 9, 7,
  377. // Cost
  378. // f g h
  379. ( 1 + 4 + 1 + 4 + 9 + 1 + 4 + 9 + 16) / 2.0,
  380. // Residuals
  381. { 1.0, 2.0, // f
  382. 1.0, 2.0, 3.0, // g
  383. 1.0, 2.0, 3.0, 4.0 // h
  384. },
  385. // Gradient
  386. { 15.0, 30.0, // x
  387. 66.0, 99.0, // y
  388. 42.0, 126.0, 168.0 // z
  389. },
  390. // Jacobian
  391. // x y z
  392. { /* f(x, y) */ 1, 2, 2, 3, 0, 0, 0,
  393. 1, 2, 2, 3, 0, 0, 0,
  394. /* g(x, z) */ 2, 4, 0, 0, 2, 6, 8,
  395. 2, 4, 0, 0, 2, 6, 8,
  396. 2, 4, 0, 0, 2, 6, 8,
  397. /* h(y, z) */ 0, 0, 6, 9, 3, 9, 12,
  398. 0, 0, 6, 9, 3, 9, 12,
  399. 0, 0, 6, 9, 3, 9, 12,
  400. 0, 0, 6, 9, 3, 9, 12
  401. }
  402. };
  403. // clang-format on
  404. CheckAllEvaluationCombinations(expected);
  405. }
  406. TEST_P(EvaluatorTest, MultipleResidualProblemWithSomeConstantParameters) {
  407. // The values are ignored completely by the cost function.
  408. double x[2];
  409. double y[3];
  410. double z[4];
  411. // Add the parameters in explicit order to force the ordering in the program.
  412. problem.AddParameterBlock(x, 2);
  413. problem.AddParameterBlock(y, 3);
  414. problem.AddParameterBlock(z, 4);
  415. // f(x, y) in R^2
  416. problem.AddResidualBlock(
  417. new ParameterIgnoringCostFunction<1, 2, 2, 3>, nullptr, x, y);
  418. // g(x, z) in R^3
  419. problem.AddResidualBlock(
  420. new ParameterIgnoringCostFunction<2, 3, 2, 4>, nullptr, x, z);
  421. // h(y, z) in R^4
  422. problem.AddResidualBlock(
  423. new ParameterIgnoringCostFunction<3, 4, 3, 4>, nullptr, y, z);
  424. // For this test, "z" is constant.
  425. problem.SetParameterBlockConstant(z);
  426. // Create the reduced program which is missing the fixed "z" variable.
  427. // Normally, the preprocessing of the program that happens in solver_impl
  428. // takes care of this, but we don't want to invoke the solver here.
  429. Program reduced_program;
  430. std::vector<ParameterBlock*>* parameter_blocks =
  431. problem.mutable_program()->mutable_parameter_blocks();
  432. // "z" is the last parameter; save it for later and pop it off temporarily.
  433. // Note that "z" will still get read during evaluation, so it cannot be
  434. // deleted at this point.
  435. ParameterBlock* parameter_block_z = parameter_blocks->back();
  436. parameter_blocks->pop_back();
  437. // clang-format off
  438. ExpectedEvaluation expected = {
  439. // Rows/columns
  440. 9, 5,
  441. // Cost
  442. // f g h
  443. ( 1 + 4 + 1 + 4 + 9 + 1 + 4 + 9 + 16) / 2.0,
  444. // Residuals
  445. { 1.0, 2.0, // f
  446. 1.0, 2.0, 3.0, // g
  447. 1.0, 2.0, 3.0, 4.0 // h
  448. },
  449. // Gradient
  450. { 15.0, 30.0, // x
  451. 33.0, 66.0, 99.0, // y
  452. },
  453. // Jacobian
  454. // x y
  455. { /* f(x, y) */ 1, 2, 1, 2, 3,
  456. 1, 2, 1, 2, 3,
  457. /* g(x, z) */ 2, 4, 0, 0, 0,
  458. 2, 4, 0, 0, 0,
  459. 2, 4, 0, 0, 0,
  460. /* h(y, z) */ 0, 0, 3, 6, 9,
  461. 0, 0, 3, 6, 9,
  462. 0, 0, 3, 6, 9,
  463. 0, 0, 3, 6, 9
  464. }
  465. };
  466. // clang-format on
  467. CheckAllEvaluationCombinations(expected);
  468. // Restore parameter block z, so it will get freed in a consistent way.
  469. parameter_blocks->push_back(parameter_block_z);
  470. }
  471. TEST_P(EvaluatorTest, EvaluatorAbortsForResidualsThatFailToEvaluate) {
  472. // Switch the return value to failure.
  473. problem.AddResidualBlock(
  474. new ParameterIgnoringCostFunction<20, 3, 2, 3, 4>(false),
  475. nullptr,
  476. x,
  477. y,
  478. z);
  479. // The values are ignored.
  480. double state[9];
  481. std::unique_ptr<Evaluator> evaluator =
  482. CreateEvaluator(problem.mutable_program());
  483. std::unique_ptr<SparseMatrix> jacobian(evaluator->CreateJacobian());
  484. double cost;
  485. EXPECT_FALSE(evaluator->Evaluate(state, &cost, nullptr, nullptr, nullptr));
  486. }
  487. // In the pairs, the first argument is the linear solver type, and the second
  488. // argument is num_eliminate_blocks. Changing the num_eliminate_blocks only
  489. // makes sense for the schur-based solvers.
  490. //
  491. // Try all values of num_eliminate_blocks that make sense given that in the
  492. // tests a maximum of 4 parameter blocks are present.
  493. INSTANTIATE_TEST_SUITE_P(
  494. LinearSolvers,
  495. EvaluatorTest,
  496. ::testing::Values(EvaluatorTestOptions(DENSE_QR, 0),
  497. EvaluatorTestOptions(DENSE_SCHUR, 0),
  498. EvaluatorTestOptions(DENSE_SCHUR, 1),
  499. EvaluatorTestOptions(DENSE_SCHUR, 2),
  500. EvaluatorTestOptions(DENSE_SCHUR, 3),
  501. EvaluatorTestOptions(DENSE_SCHUR, 4),
  502. EvaluatorTestOptions(SPARSE_SCHUR, 0),
  503. EvaluatorTestOptions(SPARSE_SCHUR, 1),
  504. EvaluatorTestOptions(SPARSE_SCHUR, 2),
  505. EvaluatorTestOptions(SPARSE_SCHUR, 3),
  506. EvaluatorTestOptions(SPARSE_SCHUR, 4),
  507. EvaluatorTestOptions(ITERATIVE_SCHUR, 0),
  508. EvaluatorTestOptions(ITERATIVE_SCHUR, 1),
  509. EvaluatorTestOptions(ITERATIVE_SCHUR, 2),
  510. EvaluatorTestOptions(ITERATIVE_SCHUR, 3),
  511. EvaluatorTestOptions(ITERATIVE_SCHUR, 4),
  512. EvaluatorTestOptions(SPARSE_NORMAL_CHOLESKY, 0, false),
  513. EvaluatorTestOptions(SPARSE_NORMAL_CHOLESKY, 0, true)));
  514. // Simple cost function used to check if the evaluator is sensitive to
  515. // state changes.
  516. class ParameterSensitiveCostFunction : public SizedCostFunction<2, 2> {
  517. public:
  518. bool Evaluate(double const* const* parameters,
  519. double* residuals,
  520. double** jacobians) const final {
  521. double x1 = parameters[0][0];
  522. double x2 = parameters[0][1];
  523. residuals[0] = x1 * x1;
  524. residuals[1] = x2 * x2;
  525. if (jacobians != nullptr) {
  526. double* jacobian = jacobians[0];
  527. if (jacobian != nullptr) {
  528. jacobian[0] = 2.0 * x1;
  529. jacobian[1] = 0.0;
  530. jacobian[2] = 0.0;
  531. jacobian[3] = 2.0 * x2;
  532. }
  533. }
  534. return true;
  535. }
  536. };
  537. TEST(Evaluator, EvaluatorRespectsParameterChanges) {
  538. ProblemImpl problem;
  539. double x[2];
  540. x[0] = 1.0;
  541. x[1] = 1.0;
  542. problem.AddResidualBlock(new ParameterSensitiveCostFunction(), nullptr, x);
  543. Program* program = problem.mutable_program();
  544. program->SetParameterOffsetsAndIndex();
  545. Evaluator::Options options;
  546. options.linear_solver_type = DENSE_QR;
  547. options.num_eliminate_blocks = 0;
  548. options.context = problem.context();
  549. std::string error;
  550. std::unique_ptr<Evaluator> evaluator(
  551. Evaluator::Create(options, program, &error));
  552. std::unique_ptr<SparseMatrix> jacobian(evaluator->CreateJacobian());
  553. ASSERT_EQ(2, jacobian->num_rows());
  554. ASSERT_EQ(2, jacobian->num_cols());
  555. double state[2];
  556. state[0] = 2.0;
  557. state[1] = 3.0;
  558. // The original state of a residual block comes from the user's
  559. // state. So the original state is 1.0, 1.0, and the only way we get
  560. // the 2.0, 3.0 results in the following tests is if it respects the
  561. // values in the state vector.
  562. // Cost only; no residuals and no jacobian.
  563. {
  564. double cost = -1;
  565. ASSERT_TRUE(evaluator->Evaluate(state, &cost, nullptr, nullptr, nullptr));
  566. EXPECT_EQ(48.5, cost);
  567. }
  568. // Cost and residuals, no jacobian.
  569. {
  570. double cost = -1;
  571. double residuals[2] = {-2, -2};
  572. ASSERT_TRUE(evaluator->Evaluate(state, &cost, residuals, nullptr, nullptr));
  573. EXPECT_EQ(48.5, cost);
  574. EXPECT_EQ(4, residuals[0]);
  575. EXPECT_EQ(9, residuals[1]);
  576. }
  577. // Cost, residuals, and jacobian.
  578. {
  579. double cost = -1;
  580. double residuals[2] = {-2, -2};
  581. SetSparseMatrixConstant(jacobian.get(), -1);
  582. ASSERT_TRUE(
  583. evaluator->Evaluate(state, &cost, residuals, nullptr, jacobian.get()));
  584. EXPECT_EQ(48.5, cost);
  585. EXPECT_EQ(4, residuals[0]);
  586. EXPECT_EQ(9, residuals[1]);
  587. Matrix actual_jacobian;
  588. jacobian->ToDenseMatrix(&actual_jacobian);
  589. Matrix expected_jacobian(2, 2);
  590. expected_jacobian << 2 * state[0], 0, 0, 2 * state[1];
  591. EXPECT_TRUE((actual_jacobian.array() == expected_jacobian.array()).all())
  592. << "Actual:\n"
  593. << actual_jacobian << "\nExpected:\n"
  594. << expected_jacobian;
  595. }
  596. }
  597. class HugeCostFunction : public SizedCostFunction<46341, 46345> {
  598. bool Evaluate(double const* const* parameters,
  599. double* residuals,
  600. double** jacobians) const override {
  601. return true;
  602. }
  603. };
  604. TEST(Evaluator, LargeProblemDoesNotCauseCrashBlockJacobianWriter) {
  605. ProblemImpl problem;
  606. std::vector<double> x(46345);
  607. problem.AddResidualBlock(new HugeCostFunction, nullptr, x.data());
  608. Evaluator::Options options;
  609. options.linear_solver_type = SPARSE_NORMAL_CHOLESKY;
  610. options.context = problem.context();
  611. options.num_eliminate_blocks = 0;
  612. options.dynamic_sparsity = false;
  613. std::string error;
  614. auto program = problem.mutable_program();
  615. program->SetParameterOffsetsAndIndex();
  616. auto evaluator = Evaluator::Create(options, program, &error);
  617. auto jacobian = evaluator->CreateJacobian();
  618. EXPECT_EQ(jacobian, nullptr);
  619. }
  620. TEST(Evaluator, LargeProblemDoesNotCauseCrashCompressedRowJacobianWriter) {
  621. ProblemImpl problem;
  622. std::vector<double> x(46345);
  623. problem.AddResidualBlock(new HugeCostFunction, nullptr, x.data());
  624. Evaluator::Options options;
  625. // CGNR on CUDA_SPARSE is the only combination that triggers a
  626. // CompressedRowJacobianWriter.
  627. options.linear_solver_type = CGNR;
  628. options.sparse_linear_algebra_library_type = CUDA_SPARSE;
  629. options.context = problem.context();
  630. options.num_eliminate_blocks = 0;
  631. std::string error;
  632. auto program = problem.mutable_program();
  633. program->SetParameterOffsetsAndIndex();
  634. auto evaluator = Evaluator::Create(options, program, &error);
  635. auto jacobian = evaluator->CreateJacobian();
  636. EXPECT_EQ(jacobian, nullptr);
  637. }
  638. } // namespace internal
  639. } // namespace ceres