schur_eliminator_impl.h 29 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718
  1. // Ceres Solver - A fast non-linear least squares minimizer
  2. // Copyright 2022 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. // TODO(sameeragarwal): row_block_counter can perhaps be replaced by
  32. // Chunk::start ?
  33. #ifndef CERES_INTERNAL_SCHUR_ELIMINATOR_IMPL_H_
  34. #define CERES_INTERNAL_SCHUR_ELIMINATOR_IMPL_H_
  35. // Eigen has an internal threshold switching between different matrix
  36. // multiplication algorithms. In particular for matrices larger than
  37. // EIGEN_CACHEFRIENDLY_PRODUCT_THRESHOLD it uses a cache friendly
  38. // matrix matrix product algorithm that has a higher setup cost. For
  39. // matrix sizes close to this threshold, especially when the matrices
  40. // are thin and long, the default choice may not be optimal. This is
  41. // the case for us, as the default choice causes a 30% performance
  42. // regression when we moved from Eigen2 to Eigen3.
  43. #define EIGEN_CACHEFRIENDLY_PRODUCT_THRESHOLD 10
  44. // This include must come before any #ifndef check on Ceres compile options.
  45. // clang-format off
  46. #include "ceres/internal/config.h"
  47. // clang-format on
  48. #include <algorithm>
  49. #include <map>
  50. #include "Eigen/Dense"
  51. #include "ceres/block_random_access_matrix.h"
  52. #include "ceres/block_sparse_matrix.h"
  53. #include "ceres/block_structure.h"
  54. #include "ceres/internal/eigen.h"
  55. #include "ceres/internal/fixed_array.h"
  56. #include "ceres/invert_psd_matrix.h"
  57. #include "ceres/map_util.h"
  58. #include "ceres/parallel_for.h"
  59. #include "ceres/schur_eliminator.h"
  60. #include "ceres/scoped_thread_token.h"
  61. #include "ceres/small_blas.h"
  62. #include "ceres/stl_util.h"
  63. #include "ceres/thread_token_provider.h"
  64. #include "glog/logging.h"
  65. namespace ceres::internal {
  66. template <int kRowBlockSize, int kEBlockSize, int kFBlockSize>
  67. SchurEliminator<kRowBlockSize, kEBlockSize, kFBlockSize>::~SchurEliminator() {
  68. STLDeleteElements(&rhs_locks_);
  69. }
  70. template <int kRowBlockSize, int kEBlockSize, int kFBlockSize>
  71. void SchurEliminator<kRowBlockSize, kEBlockSize, kFBlockSize>::Init(
  72. int num_eliminate_blocks,
  73. bool assume_full_rank_ete,
  74. const CompressedRowBlockStructure* bs) {
  75. CHECK_GT(num_eliminate_blocks, 0)
  76. << "SchurComplementSolver cannot be initialized with "
  77. << "num_eliminate_blocks = 0.";
  78. num_eliminate_blocks_ = num_eliminate_blocks;
  79. assume_full_rank_ete_ = assume_full_rank_ete;
  80. const int num_col_blocks = bs->cols.size();
  81. const int num_row_blocks = bs->rows.size();
  82. buffer_size_ = 1;
  83. chunks_.clear();
  84. lhs_row_layout_.clear();
  85. int lhs_num_rows = 0;
  86. // Add a map object for each block in the reduced linear system
  87. // and build the row/column block structure of the reduced linear
  88. // system.
  89. lhs_row_layout_.resize(num_col_blocks - num_eliminate_blocks_);
  90. for (int i = num_eliminate_blocks_; i < num_col_blocks; ++i) {
  91. lhs_row_layout_[i - num_eliminate_blocks_] = lhs_num_rows;
  92. lhs_num_rows += bs->cols[i].size;
  93. }
  94. // TODO(sameeragarwal): Now that we may have subset block structure,
  95. // we need to make sure that we account for the fact that some
  96. // point blocks only have a "diagonal" row and nothing more.
  97. //
  98. // This likely requires a slightly different algorithm, which works
  99. // off of the number of elimination blocks.
  100. int r = 0;
  101. // Iterate over the row blocks of A, and detect the chunks. The
  102. // matrix should already have been ordered so that all rows
  103. // containing the same y block are vertically contiguous. Along
  104. // the way also compute the amount of space each chunk will need
  105. // to perform the elimination.
  106. while (r < num_row_blocks) {
  107. const int chunk_block_id = bs->rows[r].cells.front().block_id;
  108. if (chunk_block_id >= num_eliminate_blocks_) {
  109. break;
  110. }
  111. chunks_.push_back(Chunk(r));
  112. Chunk& chunk = chunks_.back();
  113. int buffer_size = 0;
  114. const int e_block_size = bs->cols[chunk_block_id].size;
  115. // Add to the chunk until the first block in the row is
  116. // different than the one in the first row for the chunk.
  117. while (r + chunk.size < num_row_blocks) {
  118. const CompressedRow& row = bs->rows[r + chunk.size];
  119. if (row.cells.front().block_id != chunk_block_id) {
  120. break;
  121. }
  122. // Iterate over the blocks in the row, ignoring the first
  123. // block since it is the one to be eliminated.
  124. for (int c = 1; c < row.cells.size(); ++c) {
  125. const Cell& cell = row.cells[c];
  126. if (InsertIfNotPresent(
  127. &(chunk.buffer_layout), cell.block_id, buffer_size)) {
  128. buffer_size += e_block_size * bs->cols[cell.block_id].size;
  129. }
  130. }
  131. buffer_size_ = std::max(buffer_size, buffer_size_);
  132. ++chunk.size;
  133. }
  134. CHECK_GT(chunk.size, 0); // This check will need to be resolved.
  135. r += chunk.size;
  136. }
  137. const Chunk& chunk = chunks_.back();
  138. uneliminated_row_begins_ = chunk.start + chunk.size;
  139. buffer_ = std::make_unique<double[]>(buffer_size_ * num_threads_);
  140. // chunk_outer_product_buffer_ only needs to store e_block_size *
  141. // f_block_size, which is always less than buffer_size_, so we just
  142. // allocate buffer_size_ per thread.
  143. chunk_outer_product_buffer_ =
  144. std::make_unique<double[]>(buffer_size_ * num_threads_);
  145. STLDeleteElements(&rhs_locks_);
  146. rhs_locks_.resize(num_col_blocks - num_eliminate_blocks_);
  147. for (int i = 0; i < num_col_blocks - num_eliminate_blocks_; ++i) {
  148. rhs_locks_[i] = new std::mutex;
  149. }
  150. }
  151. template <int kRowBlockSize, int kEBlockSize, int kFBlockSize>
  152. void SchurEliminator<kRowBlockSize, kEBlockSize, kFBlockSize>::Eliminate(
  153. const BlockSparseMatrixData& A,
  154. const double* b,
  155. const double* D,
  156. BlockRandomAccessMatrix* lhs,
  157. double* rhs) {
  158. if (lhs->num_rows() > 0) {
  159. lhs->SetZero();
  160. if (rhs) {
  161. VectorRef(rhs, lhs->num_rows()).setZero();
  162. }
  163. }
  164. const CompressedRowBlockStructure* bs = A.block_structure();
  165. const int num_col_blocks = bs->cols.size();
  166. // Add the diagonal to the schur complement.
  167. if (D != nullptr) {
  168. ParallelFor(context_,
  169. num_eliminate_blocks_,
  170. num_col_blocks,
  171. num_threads_,
  172. [&](int i) {
  173. const int block_id = i - num_eliminate_blocks_;
  174. int r, c, row_stride, col_stride;
  175. CellInfo* cell_info = lhs->GetCell(
  176. block_id, block_id, &r, &c, &row_stride, &col_stride);
  177. if (cell_info != nullptr) {
  178. const int block_size = bs->cols[i].size;
  179. typename EigenTypes<Eigen::Dynamic>::ConstVectorRef diag(
  180. D + bs->cols[i].position, block_size);
  181. MatrixRef m(cell_info->values, row_stride, col_stride);
  182. m.block(r, c, block_size, block_size).diagonal() +=
  183. diag.array().square().matrix();
  184. }
  185. });
  186. }
  187. // Eliminate y blocks one chunk at a time. For each chunk, compute
  188. // the entries of the normal equations and the gradient vector block
  189. // corresponding to the y block and then apply Gaussian elimination
  190. // to them. The matrix ete stores the normal matrix corresponding to
  191. // the block being eliminated and array buffer_ contains the
  192. // non-zero blocks in the row corresponding to this y block in the
  193. // normal equations. This computation is done in
  194. // ChunkDiagonalBlockAndGradient. UpdateRhs then applies gaussian
  195. // elimination to the rhs of the normal equations, updating the rhs
  196. // of the reduced linear system by modifying rhs blocks for all the
  197. // z blocks that share a row block/residual term with the y
  198. // block. EliminateRowOuterProduct does the corresponding operation
  199. // for the lhs of the reduced linear system.
  200. ParallelFor(
  201. context_,
  202. 0,
  203. int(chunks_.size()),
  204. num_threads_,
  205. [&](int thread_id, int i) {
  206. double* buffer = buffer_.get() + thread_id * buffer_size_;
  207. const Chunk& chunk = chunks_[i];
  208. const int e_block_id = bs->rows[chunk.start].cells.front().block_id;
  209. const int e_block_size = bs->cols[e_block_id].size;
  210. VectorRef(buffer, buffer_size_).setZero();
  211. typename EigenTypes<kEBlockSize, kEBlockSize>::Matrix ete(e_block_size,
  212. e_block_size);
  213. if (D != nullptr) {
  214. const typename EigenTypes<kEBlockSize>::ConstVectorRef diag(
  215. D + bs->cols[e_block_id].position, e_block_size);
  216. ete = diag.array().square().matrix().asDiagonal();
  217. } else {
  218. ete.setZero();
  219. }
  220. FixedArray<double, 8> g(e_block_size);
  221. typename EigenTypes<kEBlockSize>::VectorRef gref(g.data(),
  222. e_block_size);
  223. gref.setZero();
  224. // We are going to be computing
  225. //
  226. // S += F'F - F'E(E'E)^{-1}E'F
  227. //
  228. // for each Chunk. The computation is broken down into a number of
  229. // function calls as below.
  230. // Compute the outer product of the e_blocks with themselves (ete
  231. // = E'E). Compute the product of the e_blocks with the
  232. // corresponding f_blocks (buffer = E'F), the gradient of the terms
  233. // in this chunk (g) and add the outer product of the f_blocks to
  234. // Schur complement (S += F'F).
  235. ChunkDiagonalBlockAndGradient(
  236. chunk, A, b, chunk.start, &ete, g.data(), buffer, lhs);
  237. // Normally one wouldn't compute the inverse explicitly, but
  238. // e_block_size will typically be a small number like 3, in
  239. // which case its much faster to compute the inverse once and
  240. // use it to multiply other matrices/vectors instead of doing a
  241. // Solve call over and over again.
  242. typename EigenTypes<kEBlockSize, kEBlockSize>::Matrix inverse_ete =
  243. InvertPSDMatrix<kEBlockSize>(assume_full_rank_ete_, ete);
  244. // For the current chunk compute and update the rhs of the reduced
  245. // linear system.
  246. //
  247. // rhs = F'b - F'E(E'E)^(-1) E'b
  248. if (rhs) {
  249. FixedArray<double, 8> inverse_ete_g(e_block_size);
  250. MatrixVectorMultiply<kEBlockSize, kEBlockSize, 0>(
  251. inverse_ete.data(),
  252. e_block_size,
  253. e_block_size,
  254. g.data(),
  255. inverse_ete_g.data());
  256. UpdateRhs(chunk, A, b, chunk.start, inverse_ete_g.data(), rhs);
  257. }
  258. // S -= F'E(E'E)^{-1}E'F
  259. ChunkOuterProduct(
  260. thread_id, bs, inverse_ete, buffer, chunk.buffer_layout, lhs);
  261. });
  262. // For rows with no e_blocks, the Schur complement update reduces to
  263. // S += F'F.
  264. NoEBlockRowsUpdate(A, b, uneliminated_row_begins_, lhs, rhs);
  265. }
  266. template <int kRowBlockSize, int kEBlockSize, int kFBlockSize>
  267. void SchurEliminator<kRowBlockSize, kEBlockSize, kFBlockSize>::BackSubstitute(
  268. const BlockSparseMatrixData& A,
  269. const double* b,
  270. const double* D,
  271. const double* z,
  272. double* y) {
  273. const CompressedRowBlockStructure* bs = A.block_structure();
  274. const double* values = A.values();
  275. ParallelFor(context_, 0, int(chunks_.size()), num_threads_, [&](int i) {
  276. const Chunk& chunk = chunks_[i];
  277. const int e_block_id = bs->rows[chunk.start].cells.front().block_id;
  278. const int e_block_size = bs->cols[e_block_id].size;
  279. double* y_ptr = y + bs->cols[e_block_id].position;
  280. typename EigenTypes<kEBlockSize>::VectorRef y_block(y_ptr, e_block_size);
  281. typename EigenTypes<kEBlockSize, kEBlockSize>::Matrix ete(e_block_size,
  282. e_block_size);
  283. if (D != nullptr) {
  284. const typename EigenTypes<kEBlockSize>::ConstVectorRef diag(
  285. D + bs->cols[e_block_id].position, e_block_size);
  286. ete = diag.array().square().matrix().asDiagonal();
  287. } else {
  288. ete.setZero();
  289. }
  290. for (int j = 0; j < chunk.size; ++j) {
  291. const CompressedRow& row = bs->rows[chunk.start + j];
  292. const Cell& e_cell = row.cells.front();
  293. DCHECK_EQ(e_block_id, e_cell.block_id);
  294. FixedArray<double, 8> sj(row.block.size);
  295. typename EigenTypes<kRowBlockSize>::VectorRef(sj.data(), row.block.size) =
  296. typename EigenTypes<kRowBlockSize>::ConstVectorRef(
  297. b + bs->rows[chunk.start + j].block.position, row.block.size);
  298. for (int c = 1; c < row.cells.size(); ++c) {
  299. const int f_block_id = row.cells[c].block_id;
  300. const int f_block_size = bs->cols[f_block_id].size;
  301. const int r_block = f_block_id - num_eliminate_blocks_;
  302. // clang-format off
  303. MatrixVectorMultiply<kRowBlockSize, kFBlockSize, -1>(
  304. values + row.cells[c].position, row.block.size, f_block_size,
  305. z + lhs_row_layout_[r_block],
  306. sj.data());
  307. }
  308. MatrixTransposeVectorMultiply<kRowBlockSize, kEBlockSize, 1>(
  309. values + e_cell.position, row.block.size, e_block_size,
  310. sj.data(),
  311. y_ptr);
  312. MatrixTransposeMatrixMultiply
  313. <kRowBlockSize, kEBlockSize, kRowBlockSize, kEBlockSize, 1>(
  314. values + e_cell.position, row.block.size, e_block_size,
  315. values + e_cell.position, row.block.size, e_block_size,
  316. ete.data(), 0, 0, e_block_size, e_block_size);
  317. // clang-format on
  318. }
  319. y_block =
  320. InvertPSDMatrix<kEBlockSize>(assume_full_rank_ete_, ete) * y_block;
  321. });
  322. }
  323. // Update the rhs of the reduced linear system. Compute
  324. //
  325. // F'b - F'E(E'E)^(-1) E'b
  326. template <int kRowBlockSize, int kEBlockSize, int kFBlockSize>
  327. void SchurEliminator<kRowBlockSize, kEBlockSize, kFBlockSize>::UpdateRhs(
  328. const Chunk& chunk,
  329. const BlockSparseMatrixData& A,
  330. const double* b,
  331. int row_block_counter,
  332. const double* inverse_ete_g,
  333. double* rhs) {
  334. const CompressedRowBlockStructure* bs = A.block_structure();
  335. const double* values = A.values();
  336. const int e_block_id = bs->rows[chunk.start].cells.front().block_id;
  337. const int e_block_size = bs->cols[e_block_id].size;
  338. int b_pos = bs->rows[row_block_counter].block.position;
  339. for (int j = 0; j < chunk.size; ++j) {
  340. const CompressedRow& row = bs->rows[row_block_counter + j];
  341. const Cell& e_cell = row.cells.front();
  342. typename EigenTypes<kRowBlockSize>::Vector sj =
  343. typename EigenTypes<kRowBlockSize>::ConstVectorRef(b + b_pos,
  344. row.block.size);
  345. // clang-format off
  346. MatrixVectorMultiply<kRowBlockSize, kEBlockSize, -1>(
  347. values + e_cell.position, row.block.size, e_block_size,
  348. inverse_ete_g, sj.data());
  349. // clang-format on
  350. for (int c = 1; c < row.cells.size(); ++c) {
  351. const int block_id = row.cells[c].block_id;
  352. const int block_size = bs->cols[block_id].size;
  353. const int block = block_id - num_eliminate_blocks_;
  354. auto lock = MakeConditionalLock(num_threads_, *rhs_locks_[block]);
  355. // clang-format off
  356. MatrixTransposeVectorMultiply<kRowBlockSize, kFBlockSize, 1>(
  357. values + row.cells[c].position,
  358. row.block.size, block_size,
  359. sj.data(), rhs + lhs_row_layout_[block]);
  360. // clang-format on
  361. }
  362. b_pos += row.block.size;
  363. }
  364. }
  365. // Given a Chunk - set of rows with the same e_block, e.g. in the
  366. // following Chunk with two rows.
  367. //
  368. // E F
  369. // [ y11 0 0 0 | z11 0 0 0 z51]
  370. // [ y12 0 0 0 | z12 z22 0 0 0]
  371. //
  372. // this function computes twp matrices. The diagonal block matrix
  373. //
  374. // ete = y11 * y11' + y12 * y12'
  375. //
  376. // and the off diagonal blocks in the Gauss Newton Hessian.
  377. //
  378. // buffer = [y11'(z11 + z12), y12' * z22, y11' * z51]
  379. //
  380. // which are zero compressed versions of the block sparse matrices E'E
  381. // and E'F.
  382. //
  383. // and the gradient of the e_block, E'b.
  384. template <int kRowBlockSize, int kEBlockSize, int kFBlockSize>
  385. void SchurEliminator<kRowBlockSize, kEBlockSize, kFBlockSize>::
  386. ChunkDiagonalBlockAndGradient(
  387. const Chunk& chunk,
  388. const BlockSparseMatrixData& A,
  389. const double* b,
  390. int row_block_counter,
  391. typename EigenTypes<kEBlockSize, kEBlockSize>::Matrix* ete,
  392. double* g,
  393. double* buffer,
  394. BlockRandomAccessMatrix* lhs) {
  395. const CompressedRowBlockStructure* bs = A.block_structure();
  396. const double* values = A.values();
  397. int b_pos = bs->rows[row_block_counter].block.position;
  398. const int e_block_size = ete->rows();
  399. // Iterate over the rows in this chunk, for each row, compute the
  400. // contribution of its F blocks to the Schur complement, the
  401. // contribution of its E block to the matrix EE' (ete), and the
  402. // corresponding block in the gradient vector.
  403. for (int j = 0; j < chunk.size; ++j) {
  404. const CompressedRow& row = bs->rows[row_block_counter + j];
  405. if (row.cells.size() > 1) {
  406. EBlockRowOuterProduct(A, row_block_counter + j, lhs);
  407. }
  408. // Extract the e_block, ETE += E_i' E_i
  409. const Cell& e_cell = row.cells.front();
  410. // clang-format off
  411. MatrixTransposeMatrixMultiply
  412. <kRowBlockSize, kEBlockSize, kRowBlockSize, kEBlockSize, 1>(
  413. values + e_cell.position, row.block.size, e_block_size,
  414. values + e_cell.position, row.block.size, e_block_size,
  415. ete->data(), 0, 0, e_block_size, e_block_size);
  416. // clang-format on
  417. if (b) {
  418. // g += E_i' b_i
  419. // clang-format off
  420. MatrixTransposeVectorMultiply<kRowBlockSize, kEBlockSize, 1>(
  421. values + e_cell.position, row.block.size, e_block_size,
  422. b + b_pos,
  423. g);
  424. // clang-format on
  425. }
  426. // buffer = E'F. This computation is done by iterating over the
  427. // f_blocks for each row in the chunk.
  428. for (int c = 1; c < row.cells.size(); ++c) {
  429. const int f_block_id = row.cells[c].block_id;
  430. const int f_block_size = bs->cols[f_block_id].size;
  431. double* buffer_ptr = buffer + FindOrDie(chunk.buffer_layout, f_block_id);
  432. // clang-format off
  433. MatrixTransposeMatrixMultiply
  434. <kRowBlockSize, kEBlockSize, kRowBlockSize, kFBlockSize, 1>(
  435. values + e_cell.position, row.block.size, e_block_size,
  436. values + row.cells[c].position, row.block.size, f_block_size,
  437. buffer_ptr, 0, 0, e_block_size, f_block_size);
  438. // clang-format on
  439. }
  440. b_pos += row.block.size;
  441. }
  442. }
  443. // Compute the outer product F'E(E'E)^{-1}E'F and subtract it from the
  444. // Schur complement matrix, i.e
  445. //
  446. // S -= F'E(E'E)^{-1}E'F.
  447. template <int kRowBlockSize, int kEBlockSize, int kFBlockSize>
  448. void SchurEliminator<kRowBlockSize, kEBlockSize, kFBlockSize>::
  449. ChunkOuterProduct(int thread_id,
  450. const CompressedRowBlockStructure* bs,
  451. const Matrix& inverse_ete,
  452. const double* buffer,
  453. const BufferLayoutType& buffer_layout,
  454. BlockRandomAccessMatrix* lhs) {
  455. // This is the most computationally expensive part of this
  456. // code. Profiling experiments reveal that the bottleneck is not the
  457. // computation of the right-hand matrix product, but memory
  458. // references to the left hand side.
  459. const int e_block_size = inverse_ete.rows();
  460. auto it1 = buffer_layout.begin();
  461. double* b1_transpose_inverse_ete =
  462. chunk_outer_product_buffer_.get() + thread_id * buffer_size_;
  463. // S(i,j) -= bi' * ete^{-1} b_j
  464. for (; it1 != buffer_layout.end(); ++it1) {
  465. const int block1 = it1->first - num_eliminate_blocks_;
  466. const int block1_size = bs->cols[it1->first].size;
  467. // clang-format off
  468. MatrixTransposeMatrixMultiply
  469. <kEBlockSize, kFBlockSize, kEBlockSize, kEBlockSize, 0>(
  470. buffer + it1->second, e_block_size, block1_size,
  471. inverse_ete.data(), e_block_size, e_block_size,
  472. b1_transpose_inverse_ete, 0, 0, block1_size, e_block_size);
  473. // clang-format on
  474. auto it2 = it1;
  475. for (; it2 != buffer_layout.end(); ++it2) {
  476. const int block2 = it2->first - num_eliminate_blocks_;
  477. int r, c, row_stride, col_stride;
  478. CellInfo* cell_info =
  479. lhs->GetCell(block1, block2, &r, &c, &row_stride, &col_stride);
  480. if (cell_info != nullptr) {
  481. const int block2_size = bs->cols[it2->first].size;
  482. auto lock = MakeConditionalLock(num_threads_, cell_info->m);
  483. // clang-format off
  484. MatrixMatrixMultiply
  485. <kFBlockSize, kEBlockSize, kEBlockSize, kFBlockSize, -1>(
  486. b1_transpose_inverse_ete, block1_size, e_block_size,
  487. buffer + it2->second, e_block_size, block2_size,
  488. cell_info->values, r, c, row_stride, col_stride);
  489. // clang-format on
  490. }
  491. }
  492. }
  493. }
  494. // For rows with no e_blocks, the Schur complement update reduces to S
  495. // += F'F. This function iterates over the rows of A with no e_block,
  496. // and calls NoEBlockRowOuterProduct on each row.
  497. template <int kRowBlockSize, int kEBlockSize, int kFBlockSize>
  498. void SchurEliminator<kRowBlockSize, kEBlockSize, kFBlockSize>::
  499. NoEBlockRowsUpdate(const BlockSparseMatrixData& A,
  500. const double* b,
  501. int row_block_counter,
  502. BlockRandomAccessMatrix* lhs,
  503. double* rhs) {
  504. const CompressedRowBlockStructure* bs = A.block_structure();
  505. const double* values = A.values();
  506. for (; row_block_counter < bs->rows.size(); ++row_block_counter) {
  507. NoEBlockRowOuterProduct(A, row_block_counter, lhs);
  508. if (!rhs) {
  509. continue;
  510. }
  511. const CompressedRow& row = bs->rows[row_block_counter];
  512. for (int c = 0; c < row.cells.size(); ++c) {
  513. const int block_id = row.cells[c].block_id;
  514. const int block_size = bs->cols[block_id].size;
  515. const int block = block_id - num_eliminate_blocks_;
  516. // clang-format off
  517. MatrixTransposeVectorMultiply<Eigen::Dynamic, Eigen::Dynamic, 1>(
  518. values + row.cells[c].position, row.block.size, block_size,
  519. b + row.block.position,
  520. rhs + lhs_row_layout_[block]);
  521. // clang-format on
  522. }
  523. }
  524. }
  525. // A row r of A, which has no e_blocks gets added to the Schur
  526. // complement as S += r r'. This function is responsible for computing
  527. // the contribution of a single row r to the Schur complement. It is
  528. // very similar in structure to EBlockRowOuterProduct except for
  529. // one difference. It does not use any of the template
  530. // parameters. This is because the algorithm used for detecting the
  531. // static structure of the matrix A only pays attention to rows with
  532. // e_blocks. This is because rows without e_blocks are rare and
  533. // typically arise from regularization terms in the original
  534. // optimization problem, and have a very different structure than the
  535. // rows with e_blocks. Including them in the static structure
  536. // detection will lead to most template parameters being set to
  537. // dynamic. Since the number of rows without e_blocks is small, the
  538. // lack of templating is not an issue.
  539. template <int kRowBlockSize, int kEBlockSize, int kFBlockSize>
  540. void SchurEliminator<kRowBlockSize, kEBlockSize, kFBlockSize>::
  541. NoEBlockRowOuterProduct(const BlockSparseMatrixData& A,
  542. int row_block_index,
  543. BlockRandomAccessMatrix* lhs) {
  544. const CompressedRowBlockStructure* bs = A.block_structure();
  545. const double* values = A.values();
  546. const CompressedRow& row = bs->rows[row_block_index];
  547. for (int i = 0; i < row.cells.size(); ++i) {
  548. const int block1 = row.cells[i].block_id - num_eliminate_blocks_;
  549. DCHECK_GE(block1, 0);
  550. const int block1_size = bs->cols[row.cells[i].block_id].size;
  551. int r, c, row_stride, col_stride;
  552. CellInfo* cell_info =
  553. lhs->GetCell(block1, block1, &r, &c, &row_stride, &col_stride);
  554. if (cell_info != nullptr) {
  555. auto lock = MakeConditionalLock(num_threads_, cell_info->m);
  556. // This multiply currently ignores the fact that this is a
  557. // symmetric outer product.
  558. // clang-format off
  559. MatrixTransposeMatrixMultiply
  560. <Eigen::Dynamic, Eigen::Dynamic, Eigen::Dynamic, Eigen::Dynamic, 1>(
  561. values + row.cells[i].position, row.block.size, block1_size,
  562. values + row.cells[i].position, row.block.size, block1_size,
  563. cell_info->values, r, c, row_stride, col_stride);
  564. // clang-format on
  565. }
  566. for (int j = i + 1; j < row.cells.size(); ++j) {
  567. const int block2 = row.cells[j].block_id - num_eliminate_blocks_;
  568. DCHECK_GE(block2, 0);
  569. DCHECK_LT(block1, block2);
  570. int r, c, row_stride, col_stride;
  571. CellInfo* cell_info =
  572. lhs->GetCell(block1, block2, &r, &c, &row_stride, &col_stride);
  573. if (cell_info != nullptr) {
  574. const int block2_size = bs->cols[row.cells[j].block_id].size;
  575. auto lock = MakeConditionalLock(num_threads_, cell_info->m);
  576. // clang-format off
  577. MatrixTransposeMatrixMultiply
  578. <Eigen::Dynamic, Eigen::Dynamic, Eigen::Dynamic, Eigen::Dynamic, 1>(
  579. values + row.cells[i].position, row.block.size, block1_size,
  580. values + row.cells[j].position, row.block.size, block2_size,
  581. cell_info->values, r, c, row_stride, col_stride);
  582. // clang-format on
  583. }
  584. }
  585. }
  586. }
  587. // For a row with an e_block, compute the contribution S += F'F. This
  588. // function has the same structure as NoEBlockRowOuterProduct, except
  589. // that this function uses the template parameters.
  590. template <int kRowBlockSize, int kEBlockSize, int kFBlockSize>
  591. void SchurEliminator<kRowBlockSize, kEBlockSize, kFBlockSize>::
  592. EBlockRowOuterProduct(const BlockSparseMatrixData& A,
  593. int row_block_index,
  594. BlockRandomAccessMatrix* lhs) {
  595. const CompressedRowBlockStructure* bs = A.block_structure();
  596. const double* values = A.values();
  597. const CompressedRow& row = bs->rows[row_block_index];
  598. for (int i = 1; i < row.cells.size(); ++i) {
  599. const int block1 = row.cells[i].block_id - num_eliminate_blocks_;
  600. DCHECK_GE(block1, 0);
  601. const int block1_size = bs->cols[row.cells[i].block_id].size;
  602. int r, c, row_stride, col_stride;
  603. CellInfo* cell_info =
  604. lhs->GetCell(block1, block1, &r, &c, &row_stride, &col_stride);
  605. if (cell_info != nullptr) {
  606. auto lock = MakeConditionalLock(num_threads_, cell_info->m);
  607. // block += b1.transpose() * b1;
  608. // clang-format off
  609. MatrixTransposeMatrixMultiply
  610. <kRowBlockSize, kFBlockSize, kRowBlockSize, kFBlockSize, 1>(
  611. values + row.cells[i].position, row.block.size, block1_size,
  612. values + row.cells[i].position, row.block.size, block1_size,
  613. cell_info->values, r, c, row_stride, col_stride);
  614. // clang-format on
  615. }
  616. for (int j = i + 1; j < row.cells.size(); ++j) {
  617. const int block2 = row.cells[j].block_id - num_eliminate_blocks_;
  618. DCHECK_GE(block2, 0);
  619. DCHECK_LT(block1, block2);
  620. const int block2_size = bs->cols[row.cells[j].block_id].size;
  621. int r, c, row_stride, col_stride;
  622. CellInfo* cell_info =
  623. lhs->GetCell(block1, block2, &r, &c, &row_stride, &col_stride);
  624. if (cell_info != nullptr) {
  625. // block += b1.transpose() * b2;
  626. auto lock = MakeConditionalLock(num_threads_, cell_info->m);
  627. // clang-format off
  628. MatrixTransposeMatrixMultiply
  629. <kRowBlockSize, kFBlockSize, kRowBlockSize, kFBlockSize, 1>(
  630. values + row.cells[i].position, row.block.size, block1_size,
  631. values + row.cells[j].position, row.block.size, block2_size,
  632. cell_info->values, r, c, row_stride, col_stride);
  633. // clang-format on
  634. }
  635. }
  636. }
  637. }
  638. } // namespace ceres::internal
  639. #endif // CERES_INTERNAL_SCHUR_ELIMINATOR_IMPL_H_