1 // 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/compressed_col_sparse_matrix_utils.h"
32
33 #include <vector>
34 #include <algorithm>
35 #include "ceres/internal/port.h"
36 #include "glog/logging.h"
37
38 namespace ceres {
39 namespace internal {
40
41 using std::vector;
42
CompressedColumnScalarMatrixToBlockMatrix(const int * scalar_rows,const int * scalar_cols,const vector<int> & row_blocks,const vector<int> & col_blocks,vector<int> * block_rows,vector<int> * block_cols)43 void CompressedColumnScalarMatrixToBlockMatrix(
44 const int* scalar_rows,
45 const int* scalar_cols,
46 const vector<int>& row_blocks,
47 const vector<int>& col_blocks,
48 vector<int>* block_rows,
49 vector<int>* block_cols) {
50 CHECK(block_rows != nullptr);
51 CHECK(block_cols != nullptr);
52 block_rows->clear();
53 block_cols->clear();
54 const int num_row_blocks = row_blocks.size();
55 const int num_col_blocks = col_blocks.size();
56
57 vector<int> row_block_starts(num_row_blocks);
58 for (int i = 0, cursor = 0; i < num_row_blocks; ++i) {
59 row_block_starts[i] = cursor;
60 cursor += row_blocks[i];
61 }
62
63 // This loop extracts the block sparsity of the scalar sparse matrix
64 // It does so by iterating over the columns, but only considering
65 // the columns corresponding to the first element of each column
66 // block. Within each column, the inner loop iterates over the rows,
67 // and detects the presence of a row block by checking for the
68 // presence of a non-zero entry corresponding to its first element.
69 block_cols->push_back(0);
70 int c = 0;
71 for (int col_block = 0; col_block < num_col_blocks; ++col_block) {
72 int column_size = 0;
73 for (int idx = scalar_cols[c]; idx < scalar_cols[c + 1]; ++idx) {
74 vector<int>::const_iterator it =
75 std::lower_bound(row_block_starts.begin(),
76 row_block_starts.end(),
77 scalar_rows[idx]);
78 // Since we are using lower_bound, it will return the row id
79 // where the row block starts. For everything but the first row
80 // of the block, where these values will be the same, we can
81 // skip, as we only need the first row to detect the presence of
82 // the block.
83 //
84 // For rows all but the first row in the last row block,
85 // lower_bound will return row_block_starts.end(), but those can
86 // be skipped like the rows in other row blocks too.
87 if (it == row_block_starts.end() || *it != scalar_rows[idx]) {
88 continue;
89 }
90
91 block_rows->push_back(it - row_block_starts.begin());
92 ++column_size;
93 }
94 block_cols->push_back(block_cols->back() + column_size);
95 c += col_blocks[col_block];
96 }
97 }
98
BlockOrderingToScalarOrdering(const vector<int> & blocks,const vector<int> & block_ordering,vector<int> * scalar_ordering)99 void BlockOrderingToScalarOrdering(const vector<int>& blocks,
100 const vector<int>& block_ordering,
101 vector<int>* scalar_ordering) {
102 CHECK_EQ(blocks.size(), block_ordering.size());
103 const int num_blocks = blocks.size();
104
105 // block_starts = [0, block1, block1 + block2 ..]
106 vector<int> block_starts(num_blocks);
107 for (int i = 0, cursor = 0; i < num_blocks ; ++i) {
108 block_starts[i] = cursor;
109 cursor += blocks[i];
110 }
111
112 scalar_ordering->resize(block_starts.back() + blocks.back());
113 int cursor = 0;
114 for (int i = 0; i < num_blocks; ++i) {
115 const int block_id = block_ordering[i];
116 const int block_size = blocks[block_id];
117 int block_position = block_starts[block_id];
118 for (int j = 0; j < block_size; ++j) {
119 (*scalar_ordering)[cursor++] = block_position++;
120 }
121 }
122 }
123 } // namespace internal
124 } // namespace ceres
125