PLSSVM - Parallel Least Squares Support Vector Machine
2.0.0
A Least Squares Support Vector Machine implementation using different backends.
|
Namespace containing the C-SVM using the CUDA backend. More...
Namespaces | |
detail | |
Namespace containing CUDA backend specific implementation details. Should not directly be used by users. | |
Classes | |
class | csvm |
A C-SVM implementation using CUDA as backend. More... | |
class | backend_exception |
Exception type thrown if a problem with the CUDA backend occurs. More... | |
Functions | |
template<typename real_type > | |
__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 kernel function. More... | |
template<typename real_type > | |
__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. More... | |
template<typename real_type > | |
__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. More... | |
template<typename real_type > | |
__global__ void | device_kernel_q_linear (real_type *q, const real_type *data_d, const real_type *data_last, const kernel_index_type num_rows, const kernel_index_type feature_range) |
Calculates the q vector using the linear C-SVM kernel. More... | |
template<typename real_type > | |
__global__ void | device_kernel_q_polynomial (real_type *q, const real_type *data_d, const real_type *data_last, const kernel_index_type num_rows, const kernel_index_type num_cols, const int degree, const real_type gamma, const real_type coef0) |
Calculates the q vector using the polynomial C-SVM kernel. More... | |
template<typename real_type > | |
__global__ void | device_kernel_q_rbf (real_type *q, const real_type *data_d, const real_type *data_last, const kernel_index_type num_rows, const kernel_index_type num_cols, const real_type gamma) |
Calculates the q vector using the radial basis functions C-SVM kernel. More... | |
template<typename real_type > | |
__global__ void | device_kernel_linear (const real_type *q, real_type *ret, const real_type *d, const real_type *data_d, const real_type QA_cost, const real_type cost, const kernel_index_type num_rows, const kernel_index_type feature_range, const real_type add, const kernel_index_type id) |
Calculates the C-SVM kernel using the linear kernel function. More... | |
template<typename real_type > | |
__global__ void | device_kernel_polynomial (const real_type *q, real_type *ret, const real_type *d, const real_type *data_d, const real_type QA_cost, const real_type cost, const kernel_index_type num_rows, const kernel_index_type num_cols, const real_type add, const int degree, const real_type gamma, const real_type coef0) |
Calculates the C-SVM kernel using the polynomial kernel function. More... | |
template<typename real_type > | |
__global__ void | device_kernel_rbf (const real_type *q, real_type *ret, const real_type *d, const real_type *data_d, const real_type QA_cost, const real_type cost, const kernel_index_type num_rows, const kernel_index_type num_cols, const real_type add, const real_type gamma) |
Calculates the C-SVM kernel using the radial basis function kernel function. More... | |
Namespace containing the C-SVM using the CUDA backend.
__global__ void plssvm::cuda::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 kernel function.
Supports multi-GPU execution.
real_type | the type of the data |
[out] | w_d | the w vector to assemble |
[in] | data_d | the one-dimension support vector matrix |
[in] | data_last_d | the last row of the support vector matrix |
[in] | alpha_d | the previously calculated weight for each data point |
[in] | num_data_points | the total number of support vectors |
[in] | num_features | the number of features per support vector |
__global__ void plssvm::cuda::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.
Currently only single GPU execution is supported.
real_type | the type of the data |
[out] | out_d | the calculated predictions |
[in] | data_d | the one-dimension support vector matrix |
[in] | data_last_d | the last row of the support vector matrix |
[in] | alpha_d | the previously calculated weight for each data point |
[in] | num_data_points | the total number of support vectors |
[in] | points | the data points to predict |
[in] | num_predict_points | the total number of data points to predict |
[in] | num_features | the number of features per support vector and point to predict |
[in] | degree | the degree parameter used in the polynomial kernel function |
[in] | gamma | the gamma parameter used in the polynomial kernel function |
[in] | coef0 | the coef0 parameter used in the polynomial kernel function |
__global__ void plssvm::cuda::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.
Currently only single GPU execution is supported.
real_type | the type of the data |
[out] | out_d | the calculated predictions |
[in] | data_d | the one-dimension support vector matrix |
[in] | data_last_d | the last row of the support vector matrix |
[in] | alpha_d | the previously calculated weight for each data point |
[in] | num_data_points | the total number of support vectors |
[in] | points | the data points to predict |
[in] | num_predict_points | the total number of data points to predict |
[in] | num_features | the number of features per support vector and point to predict |
[in] | gamma | the gamma parameter used in the rbf kernel function |
__global__ void plssvm::cuda::device_kernel_q_linear | ( | real_type * | q, |
const real_type * | data_d, | ||
const real_type * | data_last, | ||
const kernel_index_type | num_rows, | ||
const kernel_index_type | feature_range | ||
) |
Calculates the q
vector using the linear C-SVM kernel.
Supports multi-GPU execution.
real_type | the type of the data |
[out] | q | the calculated q vector |
[in] | data_d | the one-dimensional data matrix |
[in] | data_last | the last row in the data matrix |
[in] | num_rows | the number of rows in the data matrix |
[in] | feature_range | number of features used for the calculation |
__global__ void plssvm::cuda::device_kernel_q_polynomial | ( | real_type * | q, |
const real_type * | data_d, | ||
const real_type * | data_last, | ||
const kernel_index_type | num_rows, | ||
const kernel_index_type | num_cols, | ||
const int | degree, | ||
const real_type | gamma, | ||
const real_type | coef0 | ||
) |
Calculates the q
vector using the polynomial C-SVM kernel.
Currently only single GPU execution is supported.
real_type | the type of the data |
[out] | q | the calculated q vector |
[in] | data_d | the one-dimensional data matrix |
[in] | data_last | the last row in the data matrix |
[in] | num_rows | the number of rows in the data matrix |
[in] | num_cols | the number of columns in the data matrix |
[in] | degree | the degree parameter used in the polynomial kernel function |
[in] | gamma | the gamma parameter used in the polynomial kernel function |
[in] | coef0 | the coef0 parameter used in the polynomial kernel function |
__global__ void plssvm::cuda::device_kernel_q_rbf | ( | real_type * | q, |
const real_type * | data_d, | ||
const real_type * | data_last, | ||
const kernel_index_type | num_rows, | ||
const kernel_index_type | num_cols, | ||
const real_type | gamma | ||
) |
Calculates the q
vector using the radial basis functions C-SVM kernel.
Currently only single GPU execution is supported.
real_type | the type of the data |
[out] | q | the calculated q vector |
[in] | data_d | the one-dimensional data matrix |
[in] | data_last | the last row in the data matrix |
[in] | num_rows | the number of rows in the data matrix |
[in] | num_cols | the number of columns in the data matrix |
[in] | gamma | the gamma parameter used in the rbf kernel function |
__global__ void plssvm::cuda::device_kernel_linear | ( | const real_type * | q, |
real_type * | ret, | ||
const real_type * | d, | ||
const real_type * | data_d, | ||
const real_type | QA_cost, | ||
const real_type | cost, | ||
const kernel_index_type | num_rows, | ||
const kernel_index_type | feature_range, | ||
const real_type | add, | ||
const kernel_index_type | id | ||
) |
Calculates the C-SVM kernel using the linear kernel function.
Supports multi-GPU execution.
real_type | the type of the data |
[in] | q | the q vector |
[out] | ret | the result vector |
[in] | d | the right-hand side of the equation |
[in] | data_d | the one-dimension data matrix |
[in] | QA_cost | the bottom right matrix entry multiplied by cost |
[in] | cost | 1 / the cost parameter in the C-SVM |
[in] | num_rows | the number of columns in the data matrix |
[in] | feature_range | number of features used for the calculation on the device id |
[in] | add | denotes whether the values are added or subtracted from the result vector |
[in] | id | the id of the current device |
__global__ void plssvm::cuda::device_kernel_polynomial | ( | const real_type * | q, |
real_type * | ret, | ||
const real_type * | d, | ||
const real_type * | data_d, | ||
const real_type | QA_cost, | ||
const real_type | cost, | ||
const kernel_index_type | num_rows, | ||
const kernel_index_type | num_cols, | ||
const real_type | add, | ||
const int | degree, | ||
const real_type | gamma, | ||
const real_type | coef0 | ||
) |
Calculates the C-SVM kernel using the polynomial kernel function.
Currently only single GPU execution is supported.
real_type | the type of the data |
[in] | q | the q vector |
[out] | ret | the result vector |
[in] | d | the right-hand side of the equation |
[in] | data_d | the one-dimension data matrix |
[in] | QA_cost | he bottom right matrix entry multiplied by cost |
[in] | cost | 1 / the cost parameter in the C-SVM |
[in] | num_rows | the number of columns in the data matrix |
[in] | num_cols | the number of rows in the data matrix |
[in] | add | denotes whether the values are added or subtracted from the result vector |
[in] | degree | the degree parameter used in the polynomial kernel function |
[in] | gamma | the gamma parameter used in the polynomial kernel function |
[in] | coef0 | the coef0 parameter used in the polynomial kernel function |
__global__ void plssvm::cuda::device_kernel_rbf | ( | const real_type * | q, |
real_type * | ret, | ||
const real_type * | d, | ||
const real_type * | data_d, | ||
const real_type | QA_cost, | ||
const real_type | cost, | ||
const kernel_index_type | num_rows, | ||
const kernel_index_type | num_cols, | ||
const real_type | add, | ||
const real_type | gamma | ||
) |
Calculates the C-SVM kernel using the radial basis function kernel function.
Currently only single GPU execution is supported.
real_type | the type of the data |
[in] | q | the q vector |
[out] | ret | the result vector |
[in] | d | the right-hand side of the equation |
[in] | data_d | the one-dimension data matrix |
[in] | QA_cost | he bottom right matrix entry multiplied by cost |
[in] | cost | 1 / the cost parameter in the C-SVM |
[in] | num_rows | the number of columns in the data matrix |
[in] | num_cols | the number of rows in the data matrix |
[in] | add | denotes whether the values are added or subtracted from the result vector |
[in] | gamma | the gamma parameter used in the rbf kernel function |