blob: 784ffeb03adb695945760a027e87ef700e5d813e [file] [log] [blame]
// Ceres Solver - A fast non-linear least squares minimizer
// Copyright 2022 Google Inc. All rights reserved.
// http://ceres-solver.org/
//
// Redistribution and use in source and binary forms, with or without
// modification, are permitted provided that the following conditions are met:
//
// * Redistributions of source code must retain the above copyright notice,
// this list of conditions and the following disclaimer.
// * Redistributions in binary form must reproduce the above copyright notice,
// this list of conditions and the following disclaimer in the documentation
// and/or other materials provided with the distribution.
// * Neither the name of Google Inc. nor the names of its contributors may be
// used to endorse or promote products derived from this software without
// specific prior written permission.
//
// THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
// AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
// ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE
// LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
// CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
// SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
// INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
// CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
// ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
// POSSIBILITY OF SUCH DAMAGE.
//
// Author: sameeragarwal@google.com (Sameer Agarwal)
#include "ceres/block_random_access_sparse_matrix.h"
#include <algorithm>
#include <memory>
#include <set>
#include <utility>
#include <vector>
#include "ceres/internal/export.h"
#include "ceres/parallel_for.h"
#include "ceres/triplet_sparse_matrix.h"
#include "ceres/types.h"
#include "glog/logging.h"
namespace ceres::internal {
BlockRandomAccessSparseMatrix::BlockRandomAccessSparseMatrix(
const std::vector<Block>& blocks,
const std::set<std::pair<int, int>>& block_pairs,
ContextImpl* context,
int num_threads)
: kMaxRowBlocks(10 * 1000 * 1000),
blocks_(blocks),
context_(context),
num_threads_(num_threads) {
CHECK_LT(blocks.size(), kMaxRowBlocks);
const int num_cols = NumScalarEntries(blocks);
// Count the number of scalar non-zero entries and build the layout
// object for looking into the values array of the
// TripletSparseMatrix.
int num_nonzeros = 0;
for (const auto& block_pair : block_pairs) {
const int row_block_size = blocks_[block_pair.first].size;
const int col_block_size = blocks_[block_pair.second].size;
num_nonzeros += row_block_size * col_block_size;
}
VLOG(1) << "Matrix Size [" << num_cols << "," << num_cols << "] "
<< num_nonzeros;
tsm_ =
std::make_unique<TripletSparseMatrix>(num_cols, num_cols, num_nonzeros);
tsm_->set_num_nonzeros(num_nonzeros);
int* rows = tsm_->mutable_rows();
int* cols = tsm_->mutable_cols();
double* values = tsm_->mutable_values();
int pos = 0;
for (const auto& block_pair : block_pairs) {
const int row_block_size = blocks_[block_pair.first].size;
const int col_block_size = blocks_[block_pair.second].size;
cell_values_.emplace_back(block_pair, values + pos);
layout_[IntPairToLong(block_pair.first, block_pair.second)] =
std::make_unique<CellInfo>(values + pos);
pos += row_block_size * col_block_size;
}
// Fill the sparsity pattern of the underlying matrix.
for (const auto& block_pair : block_pairs) {
const int row_block_id = block_pair.first;
const int col_block_id = block_pair.second;
const int row_block_size = blocks_[row_block_id].size;
const int col_block_size = blocks_[col_block_id].size;
int pos =
layout_[IntPairToLong(row_block_id, col_block_id)]->values - values;
for (int r = 0; r < row_block_size; ++r) {
for (int c = 0; c < col_block_size; ++c, ++pos) {
rows[pos] = blocks_[row_block_id].position + r;
cols[pos] = blocks_[col_block_id].position + c;
values[pos] = 1.0;
DCHECK_LT(rows[pos], tsm_->num_rows());
DCHECK_LT(cols[pos], tsm_->num_rows());
}
}
}
}
CellInfo* BlockRandomAccessSparseMatrix::GetCell(int row_block_id,
int col_block_id,
int* row,
int* col,
int* row_stride,
int* col_stride) {
const auto it = layout_.find(IntPairToLong(row_block_id, col_block_id));
if (it == layout_.end()) {
return nullptr;
}
// Each cell is stored contiguously as its own little dense matrix.
*row = 0;
*col = 0;
*row_stride = blocks_[row_block_id].size;
*col_stride = blocks_[col_block_id].size;
return it->second.get();
}
// Assume that the user does not hold any locks on any cell blocks
// when they are calling SetZero.
void BlockRandomAccessSparseMatrix::SetZero() {
ParallelSetZero(
context_, num_threads_, tsm_->mutable_values(), tsm_->num_nonzeros());
}
void BlockRandomAccessSparseMatrix::SymmetricRightMultiplyAndAccumulate(
const double* x, double* y) const {
for (const auto& cell_position_and_data : cell_values_) {
const int row = cell_position_and_data.first.first;
const int row_block_size = blocks_[row].size;
const int row_block_pos = blocks_[row].position;
const int col = cell_position_and_data.first.second;
const int col_block_size = blocks_[col].size;
const int col_block_pos = blocks_[col].position;
MatrixVectorMultiply<Eigen::Dynamic, Eigen::Dynamic, 1>(
cell_position_and_data.second,
row_block_size,
col_block_size,
x + col_block_pos,
y + row_block_pos);
// Since the matrix is symmetric, but only the upper triangular
// part is stored, if the block being accessed is not a diagonal
// block, then use the same block to do the corresponding lower
// triangular multiply also.
if (row != col) {
MatrixTransposeVectorMultiply<Eigen::Dynamic, Eigen::Dynamic, 1>(
cell_position_and_data.second,
row_block_size,
col_block_size,
x + row_block_pos,
y + col_block_pos);
}
}
}
} // namespace ceres::internal