Skip to content

Commit

Permalink
add sparsity_csr remove_diagonal_elements kernel
Browse files Browse the repository at this point in the history
  • Loading branch information
upsj committed Mar 13, 2023
1 parent a4bf2d4 commit a81c60f
Show file tree
Hide file tree
Showing 14 changed files with 346 additions and 170 deletions.
1 change: 1 addition & 0 deletions common/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -13,6 +13,7 @@ set(UNIFIED_SOURCES
matrix/ell_kernels.cpp
matrix/hybrid_kernels.cpp
matrix/sellp_kernels.cpp
matrix/sparsity_csr_kernels.cpp
matrix/diagonal_kernels.cpp
multigrid/pgm_kernels.cpp
preconditioner/jacobi_kernels.cpp
Expand Down
29 changes: 0 additions & 29 deletions common/cuda_hip/matrix/sparsity_csr_kernels.hpp.inc
Original file line number Diff line number Diff line change
Expand Up @@ -110,35 +110,6 @@ __global__ __launch_bounds__(spmv_block_size) void abstract_classical_spmv(
} // namespace kernel


template <typename ValueType, typename IndexType>
void fill_in_dense(std::shared_ptr<const DefaultExecutor> exec,
const matrix::SparsityCsr<ValueType, IndexType>* input,
matrix::Dense<ValueType>* output) GKO_NOT_IMPLEMENTED;

GKO_INSTANTIATE_FOR_EACH_VALUE_AND_INDEX_TYPE(
GKO_DECLARE_SPARSITY_CSR_FILL_IN_DENSE_KERNEL);


template <typename ValueType, typename IndexType>
void count_num_diagonal_elements(
std::shared_ptr<const DefaultExecutor> exec,
const matrix::SparsityCsr<ValueType, IndexType>* matrix,
size_type* num_diagonal_elements) GKO_NOT_IMPLEMENTED;

GKO_INSTANTIATE_FOR_EACH_VALUE_AND_INDEX_TYPE(
GKO_DECLARE_SPARSITY_CSR_COUNT_NUM_DIAGONAL_ELEMENTS_KERNEL);


template <typename ValueType, typename IndexType>
void remove_diagonal_elements(
std::shared_ptr<const DefaultExecutor> exec, const IndexType* row_ptrs,
const IndexType* col_idxs,
matrix::SparsityCsr<ValueType, IndexType>* matrix) GKO_NOT_IMPLEMENTED;

GKO_INSTANTIATE_FOR_EACH_VALUE_AND_INDEX_TYPE(
GKO_DECLARE_SPARSITY_CSR_REMOVE_DIAGONAL_ELEMENTS_KERNEL);


template <typename ValueType, typename IndexType>
void transpose(std::shared_ptr<const DefaultExecutor> exec,
const matrix::SparsityCsr<ValueType, IndexType>* orig,
Expand Down
146 changes: 146 additions & 0 deletions common/unified/matrix/sparsity_csr_kernels.cpp
Original file line number Diff line number Diff line change
@@ -0,0 +1,146 @@
/*******************************<GINKGO LICENSE>******************************
Copyright (c) 2017-2023, the Ginkgo authors
All rights reserved.
Redistribution and use in source and binary forms, with or without
modification, are permitted provided that the following conditions
are met:
1. Redistributions of source code must retain the above copyright
notice, this list of conditions and the following disclaimer.
2. 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.
3. Neither the name of the copyright holder 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
HOLDER 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.
******************************<GINKGO LICENSE>*******************************/

#include "core/matrix/sparsity_csr_kernels.hpp"


#include <ginkgo/core/base/math.hpp>


#include "common/unified/base/kernel_launch.hpp"
#include "core/components/prefix_sum_kernels.hpp"


namespace gko {
namespace kernels {
namespace GKO_DEVICE_NAMESPACE {
/**
* @brief The SparsityCsr matrix format namespace.
*
* @ingroup sparsity_csr
*/
namespace sparsity_csr {


template <typename ValueType, typename IndexType>
void fill_in_dense(std::shared_ptr<const DefaultExecutor> exec,
const matrix::SparsityCsr<ValueType, IndexType>* input,
matrix::Dense<ValueType>* output)
{
run_kernel(
exec,
[] GKO_KERNEL(auto row, auto row_ptrs, auto col_idxs, auto value,
auto output) {
const auto begin = row_ptrs[row];
const auto end = row_ptrs[row + 1];
const auto val = *value;
for (auto nz = begin; nz < end; nz++) {
output(row, col_idxs[nz]) = val;
}
},
input->get_size()[0], input->get_const_row_ptrs(),
input->get_const_col_idxs(), input->get_const_value(), output);
}

GKO_INSTANTIATE_FOR_EACH_VALUE_AND_INDEX_TYPE(
GKO_DECLARE_SPARSITY_CSR_FILL_IN_DENSE_KERNEL);


template <typename ValueType, typename IndexType>
void diagonal_element_prefix_sum(
std::shared_ptr<const DefaultExecutor> exec,
const matrix::SparsityCsr<ValueType, IndexType>* matrix,
IndexType* prefix_sum)
{
const auto num_rows = matrix->get_size()[0];
run_kernel(
exec,
[] GKO_KERNEL(auto row, auto row_ptrs, auto col_idxs, auto prefix_sum) {
const auto begin = row_ptrs[row];
const auto end = row_ptrs[row + 1];
IndexType count = 0;
for (auto nz = begin; nz < end; nz++) {
if (col_idxs[nz] == row) {
count++;
}
}
prefix_sum[row] = count;
},
num_rows, matrix->get_const_row_ptrs(), matrix->get_const_col_idxs(),
prefix_sum);
components::prefix_sum(exec, prefix_sum, num_rows + 1);
}

GKO_INSTANTIATE_FOR_EACH_VALUE_AND_INDEX_TYPE(
GKO_DECLARE_SPARSITY_CSR_DIAGONAL_ELEMENT_PREFIX_SUM_KERNEL);


template <typename ValueType, typename IndexType>
void remove_diagonal_elements(std::shared_ptr<const DefaultExecutor> exec,
const IndexType* row_ptrs,
const IndexType* col_idxs,
const IndexType* diag_prefix_sum,
matrix::SparsityCsr<ValueType, IndexType>* matrix)
{
const auto num_rows = matrix->get_size()[0];
run_kernel(
exec,
[] GKO_KERNEL(auto row, auto in_row_ptrs, auto in_col_idxs,
auto diag_prefix_sum, auto out_row_ptrs,
auto out_col_idxs) {
const auto in_begin = in_row_ptrs[row];
const auto in_end = in_row_ptrs[row + 1];
auto out_idx = in_begin - diag_prefix_sum[row];
for (auto nz = in_begin; nz < in_end; nz++) {
const auto col = in_col_idxs[nz];
if (col != row) {
out_col_idxs[out_idx] = col;
out_idx++;
}
}
if (row == 0) {
out_row_ptrs[0] = 0;
}
out_row_ptrs[row + 1] = out_idx;
},
num_rows, row_ptrs, col_idxs, diag_prefix_sum, matrix->get_row_ptrs(),
matrix->get_col_idxs());
}

GKO_INSTANTIATE_FOR_EACH_VALUE_AND_INDEX_TYPE(
GKO_DECLARE_SPARSITY_CSR_REMOVE_DIAGONAL_ELEMENTS_KERNEL);


} // namespace sparsity_csr
} // namespace GKO_DEVICE_NAMESPACE
} // namespace kernels
} // namespace gko
2 changes: 1 addition & 1 deletion core/device_hooks/common_kernels.inc.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -499,7 +499,7 @@ GKO_STUB_MIXED_VALUE_AND_INDEX_TYPE(
GKO_DECLARE_SPARSITY_CSR_ADVANCED_SPMV_KERNEL);
GKO_STUB_VALUE_AND_INDEX_TYPE(GKO_DECLARE_SPARSITY_CSR_FILL_IN_DENSE_KERNEL);
GKO_STUB_VALUE_AND_INDEX_TYPE(
GKO_DECLARE_SPARSITY_CSR_COUNT_NUM_DIAGONAL_ELEMENTS_KERNEL);
GKO_DECLARE_SPARSITY_CSR_DIAGONAL_ELEMENT_PREFIX_SUM_KERNEL);
GKO_STUB_VALUE_AND_INDEX_TYPE(
GKO_DECLARE_SPARSITY_CSR_REMOVE_DIAGONAL_ELEMENTS_KERNEL);
GKO_STUB_VALUE_AND_INDEX_TYPE(GKO_DECLARE_SPARSITY_CSR_TRANSPOSE_KERNEL);
Expand Down
16 changes: 10 additions & 6 deletions core/matrix/sparsity_csr.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -56,8 +56,8 @@ namespace {
GKO_REGISTER_OPERATION(spmv, sparsity_csr::spmv);
GKO_REGISTER_OPERATION(advanced_spmv, sparsity_csr::advanced_spmv);
GKO_REGISTER_OPERATION(transpose, sparsity_csr::transpose);
GKO_REGISTER_OPERATION(count_num_diagonal_elements,
sparsity_csr::count_num_diagonal_elements);
GKO_REGISTER_OPERATION(diagonal_element_prefix_sum,
sparsity_csr::diagonal_element_prefix_sum);
GKO_REGISTER_OPERATION(convert_idxs_to_ptrs, components::convert_idxs_to_ptrs);
GKO_REGISTER_OPERATION(fill_in_dense, sparsity_csr::fill_in_dense);
GKO_REGISTER_OPERATION(remove_diagonal_elements,
Expand Down Expand Up @@ -270,15 +270,19 @@ SparsityCsr<ValueType, IndexType>::to_adjacency_matrix() const
auto exec = this->get_executor();
// Adjacency matrix has to be square.
GKO_ASSERT_IS_SQUARE_MATRIX(this);
size_type num_diagonal_elements = 0;
exec->run(sparsity_csr::make_count_num_diagonal_elements(
this, &num_diagonal_elements));
const auto num_rows = this->get_size()[0];
array<IndexType> diag_prefix_sum{exec, num_rows + 1};
exec->run(sparsity_csr::make_diagonal_element_prefix_sum(
this, diag_prefix_sum.get_data()));
const auto num_diagonal_elements = static_cast<size_type>(
exec->copy_val_to_host(diag_prefix_sum.get_const_data() + num_rows));
auto adj_mat =
SparsityCsr::create(exec, this->get_size(),
this->get_num_nonzeros() - num_diagonal_elements);

exec->run(sparsity_csr::make_remove_diagonal_elements(
this->get_const_row_ptrs(), this->get_const_col_idxs(), adj_mat.get()));
this->get_const_row_ptrs(), this->get_const_col_idxs(),
diag_prefix_sum.get_const_data(), adj_mat.get()));
return std::move(adj_mat);
}

Expand Down
9 changes: 5 additions & 4 deletions core/matrix/sparsity_csr_kernels.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -75,14 +75,15 @@ namespace kernels {
void remove_diagonal_elements( \
std::shared_ptr<const DefaultExecutor> exec, \
const IndexType* row_ptrs, const IndexType* col_idxs, \
const IndexType* diag_prefix_sum, \
matrix::SparsityCsr<ValueType, IndexType>* matrix)

#define GKO_DECLARE_SPARSITY_CSR_COUNT_NUM_DIAGONAL_ELEMENTS_KERNEL(ValueType, \
#define GKO_DECLARE_SPARSITY_CSR_DIAGONAL_ELEMENT_PREFIX_SUM_KERNEL(ValueType, \
IndexType) \
void count_num_diagonal_elements( \
void diagonal_element_prefix_sum( \
std::shared_ptr<const DefaultExecutor> exec, \
const matrix::SparsityCsr<ValueType, IndexType>* matrix, \
size_type* num_diagonal_elements)
IndexType* prefix_sum)

#define GKO_DECLARE_SPARSITY_CSR_TRANSPOSE_KERNEL(ValueType, IndexType) \
void transpose(std::shared_ptr<const DefaultExecutor> exec, \
Expand Down Expand Up @@ -116,7 +117,7 @@ namespace kernels {
GKO_DECLARE_SPARSITY_CSR_REMOVE_DIAGONAL_ELEMENTS_KERNEL(ValueType, \
IndexType); \
template <typename ValueType, typename IndexType> \
GKO_DECLARE_SPARSITY_CSR_COUNT_NUM_DIAGONAL_ELEMENTS_KERNEL(ValueType, \
GKO_DECLARE_SPARSITY_CSR_DIAGONAL_ELEMENT_PREFIX_SUM_KERNEL(ValueType, \
IndexType); \
template <typename ValueType, typename IndexType> \
GKO_DECLARE_SPARSITY_CSR_TRANSPOSE_KERNEL(ValueType, IndexType); \
Expand Down
29 changes: 0 additions & 29 deletions dpcpp/matrix/sparsity_csr_kernels.dp.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -290,35 +290,6 @@ GKO_INSTANTIATE_FOR_EACH_MIXED_VALUE_AND_INDEX_TYPE(
GKO_DECLARE_SPARSITY_CSR_ADVANCED_SPMV_KERNEL);


template <typename ValueType, typename IndexType>
void fill_in_dense(std::shared_ptr<const DefaultExecutor> exec,
const matrix::SparsityCsr<ValueType, IndexType>* input,
matrix::Dense<ValueType>* output) GKO_NOT_IMPLEMENTED;

GKO_INSTANTIATE_FOR_EACH_VALUE_AND_INDEX_TYPE(
GKO_DECLARE_SPARSITY_CSR_FILL_IN_DENSE_KERNEL);


template <typename ValueType, typename IndexType>
void count_num_diagonal_elements(
std::shared_ptr<const DpcppExecutor> exec,
const matrix::SparsityCsr<ValueType, IndexType>* matrix,
size_type* num_diagonal_elements) GKO_NOT_IMPLEMENTED;

GKO_INSTANTIATE_FOR_EACH_VALUE_AND_INDEX_TYPE(
GKO_DECLARE_SPARSITY_CSR_COUNT_NUM_DIAGONAL_ELEMENTS_KERNEL);


template <typename ValueType, typename IndexType>
void remove_diagonal_elements(
std::shared_ptr<const DpcppExecutor> exec, const IndexType* row_ptrs,
const IndexType* col_idxs,
matrix::SparsityCsr<ValueType, IndexType>* matrix) GKO_NOT_IMPLEMENTED;

GKO_INSTANTIATE_FOR_EACH_VALUE_AND_INDEX_TYPE(
GKO_DECLARE_SPARSITY_CSR_REMOVE_DIAGONAL_ELEMENTS_KERNEL);


template <typename ValueType, typename IndexType>
void transpose(std::shared_ptr<const DpcppExecutor> exec,
const matrix::SparsityCsr<ValueType, IndexType>* orig,
Expand Down
1 change: 1 addition & 0 deletions include/ginkgo/ginkgo.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -123,6 +123,7 @@ OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
#include <ginkgo/core/preconditioner/isai.hpp>
#include <ginkgo/core/preconditioner/jacobi.hpp>

#include <ginkgo/core/reorder/nested_dissection.hpp>
#include <ginkgo/core/reorder/rcm.hpp>
#include <ginkgo/core/reorder/reordering_base.hpp>
#include <ginkgo/core/reorder/scaled_reordered.hpp>
Expand Down
81 changes: 0 additions & 81 deletions omp/matrix/sparsity_csr_kernels.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -131,87 +131,6 @@ GKO_INSTANTIATE_FOR_EACH_MIXED_VALUE_AND_INDEX_TYPE(
GKO_DECLARE_SPARSITY_CSR_ADVANCED_SPMV_KERNEL);


template <typename ValueType, typename IndexType>
void fill_in_dense(std::shared_ptr<const DefaultExecutor> exec,
const matrix::SparsityCsr<ValueType, IndexType>* input,
matrix::Dense<ValueType>* output)
{
auto row_ptrs = input->get_const_row_ptrs();
auto col_idxs = input->get_const_col_idxs();
auto val = input->get_const_value()[0];
const auto num_rows = input->get_size()[0];

#pragma omp parallel for
for (size_type row = 0; row < num_rows; ++row) {
for (auto k = row_ptrs[row]; k < row_ptrs[row + 1]; ++k) {
auto col = col_idxs[k];
output->at(row, col) = val;
}
}
}

GKO_INSTANTIATE_FOR_EACH_VALUE_AND_INDEX_TYPE(
GKO_DECLARE_SPARSITY_CSR_FILL_IN_DENSE_KERNEL);


template <typename ValueType, typename IndexType>
void count_num_diagonal_elements(
std::shared_ptr<const OmpExecutor> exec,
const matrix::SparsityCsr<ValueType, IndexType>* matrix,
size_type* num_diagonal_elements)
{
auto num_rows = matrix->get_size()[0];
auto row_ptrs = matrix->get_const_row_ptrs();
auto col_idxs = matrix->get_const_col_idxs();
size_type num_diag = 0;
for (auto i = 0; i < num_rows; ++i) {
for (auto j = row_ptrs[i]; j < row_ptrs[i + 1]; ++j) {
if (col_idxs[j] == i) {
num_diag++;
}
}
}
*num_diagonal_elements = num_diag;
}

GKO_INSTANTIATE_FOR_EACH_VALUE_AND_INDEX_TYPE(
GKO_DECLARE_SPARSITY_CSR_COUNT_NUM_DIAGONAL_ELEMENTS_KERNEL);


template <typename ValueType, typename IndexType>
void remove_diagonal_elements(std::shared_ptr<const OmpExecutor> exec,
const IndexType* row_ptrs,
const IndexType* col_idxs,
matrix::SparsityCsr<ValueType, IndexType>* matrix)
{
auto num_rows = matrix->get_size()[0];
auto adj_ptrs = matrix->get_row_ptrs();
auto adj_idxs = matrix->get_col_idxs();
size_type num_diag = 0;
adj_ptrs[0] = row_ptrs[0];
for (auto i = 0; i < num_rows; ++i) {
for (auto j = row_ptrs[i]; j < row_ptrs[i + 1]; ++j) {
if (col_idxs[j] == i) {
num_diag++;
}
}
adj_ptrs[i + 1] = row_ptrs[i + 1] - num_diag;
}
auto nnz = 0;
for (auto i = 0; i < num_rows; ++i) {
for (auto j = row_ptrs[i]; j < row_ptrs[i + 1]; ++j) {
if (col_idxs[j] != i) {
adj_idxs[nnz] = col_idxs[j];
nnz++;
}
}
}
}

GKO_INSTANTIATE_FOR_EACH_VALUE_AND_INDEX_TYPE(
GKO_DECLARE_SPARSITY_CSR_REMOVE_DIAGONAL_ELEMENTS_KERNEL);


template <typename IndexType>
inline void convert_sparsity_to_csc(size_type num_rows,
const IndexType* row_ptrs,
Expand Down
Loading

0 comments on commit a81c60f

Please sign in to comment.