Skip to content

Commit

Permalink
Review update.
Browse files Browse the repository at this point in the history
+ Update some forgotten names.
+ Generic improvements.
  • Loading branch information
pratikvn committed Oct 7, 2019
1 parent b8bd477 commit 91cf1e3
Show file tree
Hide file tree
Showing 13 changed files with 39 additions and 40 deletions.
16 changes: 7 additions & 9 deletions core/matrix/dense.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -76,7 +76,7 @@ GKO_REGISTER_OPERATION(convert_to_csr, dense::convert_to_csr);
GKO_REGISTER_OPERATION(convert_to_ell, dense::convert_to_ell);
GKO_REGISTER_OPERATION(convert_to_hybrid, dense::convert_to_hybrid);
GKO_REGISTER_OPERATION(convert_to_sellp, dense::convert_to_sellp);
GKO_REGISTER_OPERATION(convert_to_sparsity, dense::convert_to_sparsity);
GKO_REGISTER_OPERATION(convert_to_sparsity_csr, dense::convert_to_sparsity_csr);


} // namespace dense
Expand Down Expand Up @@ -444,10 +444,9 @@ void Dense<ValueType>::move_to(Sellp<ValueType, int64> *result)
template <typename ValueType>
void Dense<ValueType>::convert_to(SparsityCsr<ValueType, int32> *result) const
{
conversion_helper(
result, this,
dense::template make_convert_to_sparsity<decltype(result),
const Dense<ValueType> *&>);
conversion_helper(result, this,
dense::template make_convert_to_sparsity_csr<
decltype(result), const Dense<ValueType> *&>);
}


Expand All @@ -461,10 +460,9 @@ void Dense<ValueType>::move_to(SparsityCsr<ValueType, int32> *result)
template <typename ValueType>
void Dense<ValueType>::convert_to(SparsityCsr<ValueType, int64> *result) const
{
conversion_helper(
result, this,
dense::template make_convert_to_sparsity<decltype(result),
const Dense<ValueType> *&>);
conversion_helper(result, this,
dense::template make_convert_to_sparsity_csr<
decltype(result), const Dense<ValueType> *&>);
}


Expand Down
8 changes: 4 additions & 4 deletions core/matrix/dense_kernels.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -98,10 +98,10 @@ namespace kernels {
matrix::Sellp<_type, _prec> *other, \
const matrix::Dense<_type> *source)

#define GKO_DECLARE_DENSE_CONVERT_TO_SPARSITY_CSR_KERNEL(_type, _prec) \
void convert_to_sparsity(std::shared_ptr<const DefaultExecutor> exec, \
matrix::SparsityCsr<_type, _prec> *other, \
const matrix::Dense<_type> *source)
#define GKO_DECLARE_DENSE_CONVERT_TO_SPARSITY_CSR_KERNEL(_type, _prec) \
void convert_to_sparsity_csr(std::shared_ptr<const DefaultExecutor> exec, \
matrix::SparsityCsr<_type, _prec> *other, \
const matrix::Dense<_type> *source)

#define GKO_DECLARE_DENSE_COUNT_NONZEROS_KERNEL(_type) \
void count_nonzeros(std::shared_ptr<const DefaultExecutor> exec, \
Expand Down
2 changes: 1 addition & 1 deletion core/matrix/sparsity_csr.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -171,7 +171,7 @@ SparsityCsr<ValueType, IndexType>::to_adjacency_matrix() const
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));
this, &num_diagonal_elements));
auto adj_mat =
SparsityCsr::create(exec, this->get_size(),
this->get_num_nonzeros() - num_diagonal_elements);
Expand Down
2 changes: 1 addition & 1 deletion core/matrix/sparsity_csr_kernels.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -68,7 +68,7 @@ namespace kernels {
void count_num_diagonal_elements( \
std::shared_ptr<const DefaultExecutor> exec, \
const matrix::SparsityCsr<ValueType, IndexType> *matrix, \
size_type &num_diagonal_elements)
size_type *num_diagonal_elements)

#define GKO_DECLARE_SPARSITY_CSR_TRANSPOSE_KERNEL(ValueType, IndexType) \
void transpose(std::shared_ptr<const DefaultExecutor> exec, \
Expand Down
6 changes: 3 additions & 3 deletions cuda/matrix/dense_kernels.cu
Original file line number Diff line number Diff line change
Expand Up @@ -732,9 +732,9 @@ GKO_INSTANTIATE_FOR_EACH_VALUE_AND_INDEX_TYPE(


template <typename ValueType, typename IndexType>
void convert_to_sparsity(std::shared_ptr<const CudaExecutor> exec,
matrix::SparsityCsr<ValueType, IndexType> *result,
const matrix::Dense<ValueType> *source)
void convert_to_sparsity_csr(std::shared_ptr<const CudaExecutor> exec,
matrix::SparsityCsr<ValueType, IndexType> *result,
const matrix::Dense<ValueType> *source)
GKO_NOT_IMPLEMENTED;

GKO_INSTANTIATE_FOR_EACH_VALUE_AND_INDEX_TYPE(
Expand Down
2 changes: 1 addition & 1 deletion cuda/matrix/sparsity_csr_kernels.cu
Original file line number Diff line number Diff line change
Expand Up @@ -92,7 +92,7 @@ template <typename ValueType, typename IndexType>
void count_num_diagonal_elements(
std::shared_ptr<const CudaExecutor> exec,
const matrix::SparsityCsr<ValueType, IndexType> *matrix,
size_type &num_diagonal_elements) GKO_NOT_IMPLEMENTED;
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);
Expand Down
6 changes: 3 additions & 3 deletions omp/matrix/dense_kernels.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -457,9 +457,9 @@ GKO_INSTANTIATE_FOR_EACH_VALUE_AND_INDEX_TYPE(


template <typename ValueType, typename IndexType>
void convert_to_sparsity(std::shared_ptr<const OmpExecutor> exec,
matrix::SparsityCsr<ValueType, IndexType> *result,
const matrix::Dense<ValueType> *source)
void convert_to_sparsity_csr(std::shared_ptr<const OmpExecutor> exec,
matrix::SparsityCsr<ValueType, IndexType> *result,
const matrix::Dense<ValueType> *source)
{
auto num_rows = result->get_size()[0];
auto num_cols = result->get_size()[1];
Expand Down
8 changes: 4 additions & 4 deletions omp/matrix/sparsity_csr_kernels.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -68,6 +68,7 @@ void spmv(std::shared_ptr<const OmpExecutor> exec,
{
auto row_ptrs = a->get_const_row_ptrs();
auto col_idxs = a->get_const_col_idxs();
auto val = a->get_const_value()[0];

#pragma omp parallel for
for (size_type row = 0; row < a->get_size()[0]; ++row) {
Expand All @@ -76,7 +77,6 @@ void spmv(std::shared_ptr<const OmpExecutor> exec,
}
for (size_type k = row_ptrs[row];
k < static_cast<size_type>(row_ptrs[row + 1]); ++k) {
auto val = one<ValueType>();
auto col = col_idxs[k];
for (size_type j = 0; j < c->get_size()[1]; ++j) {
c->at(row, j) += val * b->at(col, j);
Expand All @@ -101,6 +101,7 @@ void advanced_spmv(std::shared_ptr<const OmpExecutor> exec,
auto col_idxs = a->get_const_col_idxs();
auto valpha = alpha->at(0, 0);
auto vbeta = beta->at(0, 0);
auto val = a->get_const_value()[0];

#pragma omp parallel for
for (size_type row = 0; row < a->get_size()[0]; ++row) {
Expand All @@ -109,7 +110,6 @@ void advanced_spmv(std::shared_ptr<const OmpExecutor> exec,
}
for (size_type k = row_ptrs[row];
k < static_cast<size_type>(row_ptrs[row + 1]); ++k) {
auto val = one<ValueType>();
auto col = col_idxs[k];
for (size_type j = 0; j < c->get_size()[1]; ++j) {
c->at(row, j) += valpha * val * b->at(col, j);
Expand All @@ -126,7 +126,7 @@ 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)
size_type *num_diagonal_elements)
{
auto num_rows = matrix->get_size()[0];
auto row_ptrs = matrix->get_const_row_ptrs();
Expand All @@ -139,7 +139,7 @@ void count_num_diagonal_elements(
}
}
}
num_diagonal_elements = num_diag;
*num_diagonal_elements = num_diag;
}

GKO_INSTANTIATE_FOR_EACH_VALUE_AND_INDEX_TYPE(
Expand Down
7 changes: 4 additions & 3 deletions omp/test/matrix/sparsity_csr_kernels.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -225,9 +225,9 @@ TEST_F(SparsityCsr, CountsNumberOfDiagElementsIsEqualToRef)
gko::size_type d_num_diags = 0;

gko::kernels::reference::sparsity_csr::count_num_diagonal_elements(
ref, mtx.get(), num_diags);
ref, mtx.get(), &num_diags);
gko::kernels::omp::sparsity_csr::count_num_diagonal_elements(
omp, dmtx.get(), d_num_diags);
omp, dmtx.get(), &d_num_diags);

ASSERT_EQ(d_num_diags, num_diags);
}
Expand All @@ -238,11 +238,12 @@ TEST_F(SparsityCsr, RemovesDiagElementsKernelIsEquivalentToRef)
set_up_apply_data();
gko::size_type num_diags = 0;
gko::kernels::reference::sparsity_csr::count_num_diagonal_elements(
ref, mtx.get(), num_diags);
ref, mtx.get(), &num_diags);
auto tmp =
Mtx::create(ref, mtx->get_size(), mtx->get_num_nonzeros() - num_diags);
auto d_tmp = Mtx::create(omp, dmtx->get_size(),
dmtx->get_num_nonzeros() - num_diags);

gko::kernels::reference::sparsity_csr::remove_diagonal_elements(
ref, tmp.get(), mtx->get_const_row_ptrs(), mtx->get_const_col_idxs());
gko::kernels::omp::sparsity_csr::remove_diagonal_elements(
Expand Down
6 changes: 3 additions & 3 deletions reference/matrix/dense_kernels.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -416,9 +416,9 @@ GKO_INSTANTIATE_FOR_EACH_VALUE_AND_INDEX_TYPE(


template <typename ValueType, typename IndexType>
void convert_to_sparsity(std::shared_ptr<const ReferenceExecutor> exec,
matrix::SparsityCsr<ValueType, IndexType> *result,
const matrix::Dense<ValueType> *source)
void convert_to_sparsity_csr(std::shared_ptr<const ReferenceExecutor> exec,
matrix::SparsityCsr<ValueType, IndexType> *result,
const matrix::Dense<ValueType> *source)
{
auto num_rows = result->get_size()[0];
auto num_cols = result->get_size()[1];
Expand Down
8 changes: 4 additions & 4 deletions reference/matrix/sparsity_csr_kernels.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -65,14 +65,14 @@ void spmv(std::shared_ptr<const ReferenceExecutor> exec,
{
auto row_ptrs = a->get_const_row_ptrs();
auto col_idxs = a->get_const_col_idxs();
auto val = a->get_const_value()[0];

for (size_type row = 0; row < a->get_size()[0]; ++row) {
for (size_type j = 0; j < c->get_size()[1]; ++j) {
c->at(row, j) = zero<ValueType>();
}
for (size_type k = row_ptrs[row];
k < static_cast<size_type>(row_ptrs[row + 1]); ++k) {
auto val = one<ValueType>();
auto col = col_idxs[k];
for (size_type j = 0; j < c->get_size()[1]; ++j) {
c->at(row, j) += val * b->at(col, j);
Expand All @@ -97,14 +97,14 @@ void advanced_spmv(std::shared_ptr<const ReferenceExecutor> exec,
auto col_idxs = a->get_const_col_idxs();
auto valpha = alpha->at(0, 0);
auto vbeta = beta->at(0, 0);
auto val = a->get_const_value()[0];

for (size_type row = 0; row < a->get_size()[0]; ++row) {
for (size_type j = 0; j < c->get_size()[1]; ++j) {
c->at(row, j) *= vbeta;
}
for (size_type k = row_ptrs[row];
k < static_cast<size_type>(row_ptrs[row + 1]); ++k) {
auto val = one<ValueType>();
auto col = col_idxs[k];
for (size_type j = 0; j < c->get_size()[1]; ++j) {
c->at(row, j) += valpha * val * b->at(col, j);
Expand All @@ -121,7 +121,7 @@ template <typename ValueType, typename IndexType>
void count_num_diagonal_elements(
std::shared_ptr<const ReferenceExecutor> exec,
const matrix::SparsityCsr<ValueType, IndexType> *matrix,
size_type &num_diagonal_elements)
size_type *num_diagonal_elements)
{
auto num_rows = matrix->get_size()[0];
auto row_ptrs = matrix->get_const_row_ptrs();
Expand All @@ -134,7 +134,7 @@ void count_num_diagonal_elements(
}
}
}
num_diagonal_elements = num_diag;
*num_diagonal_elements = num_diag;
}

GKO_INSTANTIATE_FOR_EACH_VALUE_AND_INDEX_TYPE(
Expand Down
4 changes: 2 additions & 2 deletions reference/test/matrix/sparsity_csr.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -86,7 +86,7 @@ TEST_F(SparsityCsr, CanBeCreatedFromExistingCsrMatrix)

auto mtx = Mtx::create(exec, std::move(csr_mtx));

GKO_ASSERT_MTX_NEAR(comp_mtx.get(), mtx.get(), 0.0)
GKO_ASSERT_MTX_NEAR(comp_mtx.get(), mtx.get(), 0.0);
}


Expand All @@ -99,7 +99,7 @@ TEST_F(SparsityCsr, CanBeCreatedFromExistingDenseMatrix)

auto mtx = Mtx::create(exec, std::move(dense_mtx));

GKO_ASSERT_MTX_NEAR(comp_mtx.get(), mtx.get(), 0.0)
GKO_ASSERT_MTX_NEAR(comp_mtx.get(), mtx.get(), 0.0);
}


Expand Down
4 changes: 2 additions & 2 deletions reference/test/matrix/sparsity_csr_kernels.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -288,9 +288,9 @@ TEST_F(SparsityCsr, CountsCorrectNumberOfDiagonalElements)
gko::size_type ms_num_diags = 0;

gko::kernels::reference::sparsity_csr::count_num_diagonal_elements(
exec, mtx2.get(), m2_num_diags);
exec, mtx2.get(), &m2_num_diags);
gko::kernels::reference::sparsity_csr::count_num_diagonal_elements(
exec, mtx_s.get(), ms_num_diags);
exec, mtx_s.get(), &ms_num_diags);

ASSERT_EQ(m2_num_diags, 3);
ASSERT_EQ(ms_num_diags, 2);
Expand Down

0 comments on commit 91cf1e3

Please sign in to comment.