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