12 #ifndef PLSSVM_BACKENDS_HIP_PREDICT_KERNEL_HPP_
13 #define PLSSVM_BACKENDS_HIP_PREDICT_KERNEL_HPP_
16 #include "hip/hip_runtime.h"
17 #include "hip/hip_runtime_api.h"
34 template <
typename real_type>
37 real_type temp{ 0.0 };
38 if (index < num_features) {
42 temp += alpha_d[num_data_points - 1] * data_last_d[index];
63 template <
typename real_type>
64 __global__
void device_kernel_predict_polynomial(real_type *out_d,
const real_type *data_d,
const real_type *data_last_d,
const real_type *alpha_d,
const kernel_index_type num_data_points,
const real_type *points,
const kernel_index_type num_predict_points,
const kernel_index_type num_features,
const int degree,
const real_type gamma,
const real_type coef0) {
66 const kernel_index_type predict_point_index = blockIdx.y * blockDim.y + threadIdx.y;
68 real_type temp{ 0.0 };
69 if (predict_point_index < num_predict_points) {
70 for (
kernel_index_type feature_index = 0; feature_index < num_features; ++feature_index) {
71 if (data_point_index == num_data_points - 1) {
78 temp = alpha_d[data_point_index] * pow(gamma * temp + coef0, degree);
80 atomicAdd(&out_d[predict_point_index], temp);
98 template <
typename real_type>
100 const kernel_index_type data_point_index = blockIdx.x * blockDim.x + threadIdx.x;
101 const kernel_index_type predict_point_index = blockIdx.y * blockDim.y + threadIdx.y;
103 real_type temp{ 0.0 };
104 if (predict_point_index < num_predict_points) {
105 for (
kernel_index_type feature_index = 0; feature_index < num_features; ++feature_index) {
106 if (data_point_index == num_data_points - 1) {
113 temp = alpha_d[data_point_index] * exp(-gamma * temp);
115 atomicAdd(&out_d[predict_point_index], temp);
__device__ __forceinline__ double atomicAdd(double *addr, const double val)
Atomically add the double precision val to the value denoted by addr.
Definition: atomics.cuh:24
Global type definitions and compile-time constants.
Namespace containing the C-SVM using the HIP backend.
Definition: csvm.hpp:34
__global__ void device_kernel_w_linear(real_type *w_d, const real_type *data_d, const real_type *data_last_d, const real_type *alpha_d, const kernel_index_type num_data_points, const kernel_index_type num_features)
Calculate the w vector to speed up the prediction of the labels for data points using the linear kern...
Definition: predict_kernel.hip.hpp:35
__global__ void device_kernel_predict_polynomial(real_type *out_d, const real_type *data_d, const real_type *data_last_d, const real_type *alpha_d, const kernel_index_type num_data_points, const real_type *points, const kernel_index_type num_predict_points, const kernel_index_type num_features, const int degree, const real_type gamma, const real_type coef0)
Predicts the labels for data points using the polynomial kernel function.
Definition: predict_kernel.hip.hpp:64
__global__ void device_kernel_predict_rbf(real_type *out_d, const real_type *data_d, const real_type *data_last_d, const real_type *alpha_d, const kernel_index_type num_data_points, const real_type *points, const kernel_index_type num_predict_points, const kernel_index_type num_features, const real_type gamma)
Predicts the labels for data points using the radial basis functions kernel function.
Definition: predict_kernel.hip.hpp:99
constexpr kernel_index_type THREAD_BLOCK_SIZE
Global compile-time constant used for internal caching. May be changed during the CMake configuration...
Definition: constants.hpp:25
int kernel_index_type
Integer type used inside kernels.
Definition: constants.hpp:19
constexpr kernel_index_type INTERNAL_BLOCK_SIZE
Global compile-time constant used for internal caching. May be changed during the CMake configuration...
Definition: constants.hpp:32