FEAT 3
Finite Element Analysis Toolbox
Loading...
Searching...
No Matches
min_abs_index.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/min_abs_index.hpp>
9#include <kernel/util/exception.hpp>
10#include <kernel/util/memory_pool.hpp>
11
12// includes, CUDA
13#include <cublas_v2.h>
14
15namespace FEAT
16{
17 namespace LAFEM
18 {
19 namespace Intern
20 {
21
22 Index cuda_min_abs_index(const float * x, const Index size)
23 {
24 int result;
25 cublasStatus_t status;
26 status = cublasIsamin(Util::Intern::cublas_handle, int(size), x, 1, &result);
27 if (status != CUBLAS_STATUS_SUCCESS)
28 throw InternalError(__func__, __FILE__, __LINE__, "cuda error: " + stringify(cublasGetStatusString(status)));
29 return (Index)result - 1;
30 }
31
32 Index cuda_min_abs_index(const double * x, const Index size)
33 {
34 int result;
35 cublasStatus_t status;
36 status = cublasIdamin(Util::Intern::cublas_handle, int(size), x, 1, &result);
37 if (status != CUBLAS_STATUS_SUCCESS)
38 throw InternalError(__func__, __FILE__, __LINE__, "cuda error: " + stringify(cublasGetStatusString(status)));
39 return (Index)result - 1;
40 }
41 }
42 }
43}
44
45using namespace FEAT;
46using namespace FEAT::LAFEM;
47using namespace FEAT::LAFEM::Arch;
48
49template <typename DT_>
50Index MinAbsIndex::value_cuda(const DT_ * const x, const Index size)
51{
52 Index result = Intern::cuda_min_abs_index(x, size);
53
54 cudaDeviceSynchronize();
55#ifdef FEAT_DEBUG_MODE
56 cudaError_t last_error(cudaGetLastError());
57 if (cudaSuccess != last_error)
58 throw InternalError(__func__, __FILE__, __LINE__, "CUDA error occurred in execution!\n" + stringify(cudaGetErrorString(last_error)));
59#endif
60 return result;
61}
62
63template Index MinAbsIndex::value_cuda(const float * const, const Index);
64template Index MinAbsIndex::value_cuda(const double * const, const Index);