Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

fixed a mistake in the calculate_total_cols kernels. #295

Merged
merged 1 commit into from
Apr 24, 2019
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
10 changes: 6 additions & 4 deletions cuda/matrix/csr_kernels.cu
Original file line number Diff line number Diff line change
Expand Up @@ -40,8 +40,8 @@ OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
#include <ginkgo/core/base/math.hpp>
#include <ginkgo/core/matrix/coo.hpp>
#include <ginkgo/core/matrix/dense.hpp>
#include <ginkgo/core/matrix/sellp.hpp>
#include <ginkgo/core/matrix/ell.hpp>
#include <ginkgo/core/matrix/sellp.hpp>


#include "core/matrix/dense_kernels.hpp"
Expand Down Expand Up @@ -1138,7 +1138,7 @@ __global__ __launch_bounds__(default_block_size) void fill_in_ell(
}


} // namespace kernel
} // namespace kernel


template <typename ValueType, typename IndexType>
Expand Down Expand Up @@ -1187,10 +1187,11 @@ __global__ __launch_bounds__(default_block_size) void reduce_max_nnz_per_slice(
constexpr auto warp_size = cuda_config::warp_size;
const auto warpid = tidx / warp_size;
const auto tid_in_warp = tidx % warp_size;
const auto slice_num = ceildiv(num_rows, slice_size);

size_type thread_result = 0;
for (auto i = tid_in_warp; i < slice_size; i += warp_size) {
if (warpid * warp_size + i < num_rows) {
if (warpid * slice_size + i < num_rows) {
thread_result =
max(thread_result, nnz_per_row[warpid * slice_size + i]);
}
Expand All @@ -1202,7 +1203,7 @@ __global__ __launch_bounds__(default_block_size) void reduce_max_nnz_per_slice(
warp_tile, thread_result,
[](const size_type &a, const size_type &b) { return max(a, b); });

if (tid_in_warp == 0) {
if (tid_in_warp == 0 && warpid < slice_num) {
result[warpid] = ceildiv(warp_result, stride_factor) * stride_factor;
}
}
Expand Down Expand Up @@ -1250,6 +1251,7 @@ void calculate_total_cols(std::shared_ptr<const CudaExecutor> exec,
as_cuda_type(nnz_per_row.get_const_data()),
as_cuda_type(max_nnz_per_slice.get_data()));

grid_dim = ceildiv(slice_num, default_block_size);
auto block_results = Array<size_type>(exec, grid_dim);

kernel::reduce_total_cols<<<grid_dim, default_block_size,
Expand Down
11 changes: 7 additions & 4 deletions cuda/matrix/dense_kernels.cu
Original file line number Diff line number Diff line change
Expand Up @@ -878,12 +878,13 @@ __global__ __launch_bounds__(default_block_size) void reduce_max_nnz_per_slice(
constexpr auto warp_size = cuda_config::warp_size;
const auto warpid = tidx / warp_size;
const auto tid_in_warp = tidx % warp_size;
const auto slice_num = ceildiv(num_rows, slice_size);

size_type thread_result = 0;
for (auto i = tid_in_warp; i < slice_size; i += warp_size) {
if (warpid * warp_size + i < num_rows) {
if (warpid * slice_size + i < num_rows) {
thread_result =
max(thread_result, nnz_per_row[warpid * warp_size + i]);
max(thread_result, nnz_per_row[warpid * slice_size + i]);
}
}

Expand All @@ -893,7 +894,7 @@ __global__ __launch_bounds__(default_block_size) void reduce_max_nnz_per_slice(
warp_tile, thread_result,
[](const size_type &a, const size_type &b) { return max(a, b); });

if (tid_in_warp == 0) {
if (tid_in_warp == 0 && warpid < slice_num) {
result[warpid] = ceildiv(warp_result, stride_factor) * stride_factor;
}
}
Expand Down Expand Up @@ -933,13 +934,15 @@ void calculate_total_cols(std::shared_ptr<const CudaExecutor> exec,

auto max_nnz_per_slice = Array<size_type>(exec, slice_num);

const auto grid_dim = ceildiv(slice_num, default_block_size);
auto grid_dim =
ceildiv(slice_num * cuda_config::warp_size, default_block_size);

kernel::reduce_max_nnz_per_slice<<<grid_dim, default_block_size>>>(
num_rows, slice_size, stride_factor,
as_cuda_type(nnz_per_row.get_const_data()),
as_cuda_type(max_nnz_per_slice.get_data()));

grid_dim = ceildiv(slice_num, default_block_size);
auto block_results = Array<size_type>(exec, grid_dim);

kernel::reduce_total_cols<<<grid_dim, default_block_size,
Expand Down