FEAT 3
Finite Element Analysis Toolbox
Loading...
Searching...
No Matches
diagonal.cu
1// FEAT3: Finite Element Analysis Toolbox, Version 3
2// Copyright (C) 2010 by Stefan Turek & the FEAT group
3// FEAT3 is released under the GNU General Public License version 3,
4// see the file 'copyright.txt' in the top level directory for details.
5
6// includes, FEAT
7#include <kernel/base_header.hpp>
8#include <kernel/lafem/arch/diagonal.hpp>
9#include <kernel/util/exception.hpp>
10#include <kernel/util/cuda_util.hpp>
11
12namespace FEAT
13{
14 namespace LAFEM
15 {
16 namespace Intern
17 {
18 template <typename IT_>
19 __global__ void cuda_diagonal_csr(IT_ * diag, const IT_ * col_ind, const IT_ * row_ptr, const Index rows)
20 {
21 Index idx = threadIdx.x + blockDim.x * blockIdx.x;
22 if (idx >= rows)
23 return;
24
25 const Index end(row_ptr[idx + 1]);
26 for (Index i(row_ptr[idx]); i < end; ++i)
27 {
28 if (idx == col_ind[i])
29 {
30 diag[idx] = i;
31 return;
32 }
33 }
34 diag[idx] = row_ptr[rows];
35 }
36 }
37 }
38}
39
40
41using namespace FEAT;
42using namespace FEAT::LAFEM;
43using namespace FEAT::LAFEM::Arch;
44
45template <typename IT_>
46void Diagonal::csr_cuda(IT_ * diag, const IT_ * const col_ind, const IT_ * const row_ptr, const Index rows)
47{
48 Index blocksize = Util::cuda_blocksize_axpy;
49 dim3 grid;
50 dim3 block;
51 block.x = (unsigned)blocksize;
52 grid.x = (unsigned)ceil((rows)/(double)(block.x));
53
54 FEAT::LAFEM::Intern::cuda_diagonal_csr<<<grid, block>>>(diag, col_ind, row_ptr, rows);
55
56 cudaDeviceSynchronize();
57#ifdef FEAT_DEBUG_MODE
58 cudaError_t last_error(cudaGetLastError());
59 if (cudaSuccess != last_error)
60 throw InternalError(__func__, __FILE__, __LINE__, "CUDA error occurred in execution!\n" + stringify(cudaGetErrorString(last_error)));
61#endif
62}
63
64template void Diagonal::csr_cuda(std::uint64_t *, const std::uint64_t * const, const std::uint64_t * const, const Index);
65template void Diagonal::csr_cuda(std::uint32_t *, const std::uint32_t * const, const std::uint32_t * const, const Index);