| Title: | Kernel Ridge Regression using 'RcppArmadillo' |
|---|---|
| Description: | Provides core computational operations in C++ via 'RcppArmadillo', enabling faster performance than pure R, improved numerical stability, and parallel execution with OpenMP where available. On systems without OpenMP support, the package automatically falls back to single-threaded execution with no user configuration required. For efficient model selection, it integrates with 'CVST' to provide sequential-testing cross-validation that identifies competitive hyperparameters without exhaustive grid search. The package offers a unified interface for exact kernel ridge regression and three scalable approximations—Nyström, Pivoted Cholesky, and Random Fourier Features—allowing analyses with substantially larger sample sizes than are feasible with exact KRR. It also integrates with the 'tidymodels' ecosystem via the 'parsnip' model specification 'krr_reg', and the S3 method tunable.krr_reg(). To understand the theoretical background, one can refer to Wainwright (2019) <doi:10.1017/9781108627771>. |
| Authors: | Gyeongmin Kim [aut] (Sungshin Women's University), Seyoung Lee [aut] (Sungshin Women's University), Miyoung Jang [aut] (Sungshin Women's University), Kwan-Young Bak [aut, cre, cph] (ORCID: <https://orcid.org/0000-0002-4541-160X>, Sungshin Women's University) |
| Maintainer: | Kwan-Young Bak <[email protected]> |
| License: | GPL (>= 2) |
| Version: | 0.1.2 |
| Built: | 2026-05-14 08:01:09 UTC |
| Source: | https://github.com/kybak90/fastkrr |
The FastKRR implements its core computational operations in C++ via RcppArmadillo, enabling faster performance than pure R, improved numerical stability, and parallel execution with OpenMP where available. On systems without OpenMP support, the package automatically falls back to single-threaded execution with no user configuration required. For efficient model selection, it integrates with CVST to provide sequential-testing cross-validation that identifies competitive hyperparameters without exhaustive grid search. The package offers a unified interface for exact kernel ridge regression and three widely used scalable approximations—Nyström, Pivoted Cholesky, and Random Fourier Features—allowing analyses with substantially larger sample sizes than are feasible with exact KRR while retaining strong predictive performance. This combination of a compiled backend and scalable algorithms addresses limitations of packages that rely solely on exact computation, which is often impractical for large n. It also integrates with the tidymodels ecosystem via the parsnip model specification krr_reg, and the S3 method tunable.krr_reg() (exposes tunable parameters to dials/tune); see their help pages for usage.
R/: High-level R functions and user-facing API
src/: C++ sources (kernel computation, fitting, prediction)
This package links against Rcpp and RcppArmadillo (via LinkingTo). It uses CVST, parsnip, and the tidymodels ecosystem through their public R APIs.
Maintainer: Kwan-Young Bak [email protected] (ORCID) (Sungshin Women's University) [copyright holder]
Authors:
Gyeongmin Kim [email protected] (Sungshin Women's University)
Seyoung Lee [email protected] (Sungshin Women's University)
Miyoung Jang [email protected] (Sungshin Women's University)
CVST, Rcpp, RcppArmadillo, parsnip, tidymodels
Computes low-rank kernel approximation using three methods:
Nyström approximation, Pivoted Cholesky decomposition, and
Random Fourier Features (RFF).
approx_kernel( K = NULL, X = NULL, opt = c("nystrom", "pivoted", "rff"), kernel = c("gaussian", "laplace"), m = NULL, d, rho, eps = 1e-06, W = NULL, b = NULL, n_threads = 4 )approx_kernel( K = NULL, X = NULL, opt = c("nystrom", "pivoted", "rff"), kernel = c("gaussian", "laplace"), m = NULL, d, rho, eps = 1e-06, W = NULL, b = NULL, n_threads = 4 )
K |
Exact Kernel matrix |
X |
Design matrix |
opt |
Method for constructing or approximating :
|
kernel |
Kernel type either "gaussian"or "laplace". |
m |
Approximation rank (number of random features) for the low-rank kernel approximation. If not specified, the recommended choice is
where |
d |
Design matrix's dimension ( |
rho |
Scaling parameter of the kernel ( |
eps |
Tolerance parameter used only in |
W |
Random frequency matrix |
b |
Random phase vector |
n_threads |
Number of parallel threads.
The default is 4. If the system does not support 4 threads,
it automatically falls back to 1 thread. It is applied only for |
Requirements and what to supply:
Common
d and rho must be provided (non-NULL).
nystrom / pivoted
Require a precomputed kernel matrix K; error if K is NULL.
If m is NULL, use .
For "pivoted", a tolerance eps is used; the decomposition stops early
when the next pivot (residual diagonal) drops below eps.
rff
K must be NULL (not used) and X must be provided with d = ncol(X).
The function automatically generates
W (random frequency matrix ) and
b (random phase vector ).
If the user provides them manually, both W and b must be specified and their dimensions must be compatible.
An S3 object of class "approx_kernel" containing the results of the
kernel approximation:
call: The matched function call used to create the object.
opt: The kernel approximation method actually used ("nystrom", "pivoted", "rff").
K_approx: approximated kernel matrix.
m: Kernel approximation degree.
Additional components depend on the value of opt:
nystrom
n_threads: Number of threads used in the computation.
pivoted
eps: Numerical tolerance used for early stopping in the
pivoted Cholesky decomposition.
rff
d: Input design matrix's dimension.
rho: Scaling parameter of the kernel.
W: Random frequency matrix.
b: Random phase –vector.
used_supplied_Wb: Logical; TRUE if user-supplied
W, b were used, FALSE otherwise.
n_threads: Number of threads used in the computation.
# Data setting set.seed(1) d = 1 n = 1000 m = 50 X = matrix(runif(n*d, 0, 1), nrow = n, ncol = d) y = as.vector(sin(2*pi*rowMeans(X)^3) + rnorm(n, 0, 0.1)) K = make_kernel(X, kernel = "gaussian", rho = 1) # Example: RFF approximation K_rff = approx_kernel(X = X, opt = "rff", kernel = "gaussian", m = m, d = d, rho = 1, n_threads = 1) # Exapmle: Nystrom approximation K_nystrom = approx_kernel(K = K, opt = "nystrom", m = m, d = d, rho = 1, n_threads = 1) # Example: Pivoted Cholesky approximation K_pivoted = approx_kernel(K = K, opt = "pivoted", m = m, d = d, rho = 1)# Data setting set.seed(1) d = 1 n = 1000 m = 50 X = matrix(runif(n*d, 0, 1), nrow = n, ncol = d) y = as.vector(sin(2*pi*rowMeans(X)^3) + rnorm(n, 0, 0.1)) K = make_kernel(X, kernel = "gaussian", rho = 1) # Example: RFF approximation K_rff = approx_kernel(X = X, opt = "rff", kernel = "gaussian", m = m, d = d, rho = 1, n_threads = 1) # Exapmle: Nystrom approximation K_nystrom = approx_kernel(K = K, opt = "nystrom", m = m, d = d, rho = 1, n_threads = 1) # Example: Pivoted Cholesky approximation K_pivoted = approx_kernel(K = K, opt = "pivoted", m = m, d = d, rho = 1)
Displays the main coefficient information from a fitted Kernel Ridge Regression (KRR) model,
including the original function call and the first few estimated coefficients.
The type of coefficient reported depends on the kernel approximation method:
for opt = "exact", "nystrom", or "pivoted", the coefficients represent
; for opt = "rff", they represent the coefficient ().
## S3 method for class 'krr' coef(object, ...)## S3 method for class 'krr' coef(object, ...)
object |
An S3 object of class |
... |
Additional arguments (currently ignored). |
A human-readable summary of the fitted KRR model to the console.
# Data setting set.seed(1) lambda = 1e-4 d = 1 n = 50 rho = 1 X = matrix(runif(n*d, 0, 1), nrow = n, ncol = d) y = as.vector(sin(2*pi*rowMeans(X)^3) + rnorm(n, 0, 0.1)) # Example: exact model = fastkrr(X, y, kernel = "gaussian", opt = "exact", rho = rho, lambda = 1e-4) class(model) coef(model)# Data setting set.seed(1) lambda = 1e-4 d = 1 n = 50 rho = 1 X = matrix(runif(n*d, 0, 1), nrow = n, ncol = d) y = as.vector(sin(2*pi*rowMeans(X)^3) + rnorm(n, 0, 0.1)) # Example: exact model = fastkrr(X, y, kernel = "gaussian", opt = "exact", rho = rho, lambda = 1e-4) class(model) coef(model)
Generic function for computing model error.
error(x, ...) ## Default S3 method: error(x, ...)error(x, ...) ## Default S3 method: error(x, ...)
x |
An object to compute model error for. |
... |
Additional arguments passed to methods. |
A numeric value or class-specific result.
Computes the model error for kernel ridge regression ("krr") objects.
Returns the mean squared error (MSE) between the observed responses
and the fitted values stored in the object.
## S3 method for class 'krr' error(x, ...)## S3 method for class 'krr' error(x, ...)
x |
An object of class |
... |
Additional arguments (ignored). |
This method computes:
where y and fitted.values are stored in the "krr" object attributes.
A numeric value giving the mean squared error (MSE).
summary.krr, plot.krr, predict.krr
# Data setting set.seed(1) lambda = 1e-4 d = 1 n = 50 rho = 1 X = matrix(runif(n*d, 0, 1), nrow = n, ncol = d) y = as.vector(sin(2*pi*rowMeans(X)^3) + rnorm(n, 0, 0.1)) model = fastkrr(X, y, kernel = "gaussian", lambda = 0.001) error(model)# Data setting set.seed(1) lambda = 1e-4 d = 1 n = 50 rho = 1 X = matrix(runif(n*d, 0, 1), nrow = n, ncol = d) y = as.vector(sin(2*pi*rowMeans(X)^3) + rnorm(n, 0, 0.1)) model = fastkrr(X, y, kernel = "gaussian", lambda = 0.001) error(model)
This function performs kernel ridge regression (KRR) in high-dimensional
settings. The regularization parameter can be selected via the
CVST (Cross-Validation via Sequential Testing) procedure. For scalability,
three different kernel approximation strategies are supported (Nyström approximation,
Pivoted Cholesky decomposition, Random Fourier Features(RFF)), and kernel matrix
can be computed using two methods(Gaussian kernel, Laplace kerenl).
fastkrr( x, y, kernel = "gaussian", opt = "exact", m = NULL, eps = 1e-06, rho = 1, lambda = NULL, fastcv = FALSE, n_threads = 4, verbose = TRUE )fastkrr( x, y, kernel = "gaussian", opt = "exact", m = NULL, eps = 1e-06, rho = 1, lambda = NULL, fastcv = FALSE, n_threads = 4, verbose = TRUE )
x |
Design matrix |
y |
Response variable |
kernel |
Kernel type either "gaussian"or "laplace". |
opt |
Method for constructing or approximating :
|
m |
Approximation rank(number of random features) used for the low-rank kernel approximation. If not provided by the user, it defaults to
where |
eps |
Tolerance parameter used only in |
rho |
Scaling parameter of the kernel(
|
lambda |
Regularization parameter. If |
fastcv |
If |
n_threads |
Number of parallel threads.
The default is 4. If the system does not support 4 threads,
it automatically falls back to 1 thread.
Parallelization (implemented in C++) is one of the main advantages
of this package and is applied only for |
verbose |
If TRUE, detailed progress and cross-validation results are printed to the console. If FALSE, suppresses intermediate output and only returns the final result. |
The function performs several input checks and automatic adjustments:
If x is a vector, it is converted to a one column matrix.
Otherwise, x must be a matrix; otherwise an error is thrown.
y must be a vector, and its length must match nrow(x).
kernel must be either gaussian or laplace.
opt must be one of "exact", "pivoted",
"nystrom", or "rff".
If m is NULL, it defaults to
where and .
Otherwise, m must be a positive integer.
rho must be a positive real number (default is 1).
lambda can be specified in three ways:
A positive numeric scalar, in which case the model is fitted with this single value.
A numeric vector (length >= 3) of positive values used as a tuning grid;
selection is performed by CVST cross-validation (sequential testing if
fastcv = TRUE).
NULL: use a default grid (internal setting) and tune lambda
via CVST cross-validation (sequential testing if fastcv = TRUE).
n_threads: Number of threads for parallel computation.
Default is 4. If the system has <= 3 available processors,
it uses 1.
An S3 object of class "fastkrr", which is a list containing the
results of the fitted Kernel Ridge Regression model.
coefficients: Estimated coefficient vector . Accessible via model$coefficients.
fitted.values: Fitted values . Accessible via model$fitted.values.
opt: Kernel approximation option. One of "exact", "pivoted", "nystrom", "rff".
kernel: Kernel used ("gaussian" or "laplace").
x: Input design matrix.
y: Response vector.
lambda: Regularization parameter. If NULL, tuned by cross-validation via CVST.
rho: Additional user-specified hyperparameter.
n_threads: Number of threads used for parallelization.
Additional components depend on the value of opt:
K: The full kernel matrix.
K: Exact kernel matrix .
m: Kernel approximation degree.
R: The method provides a low-rank approximation to the kernel matrix
obtained via
Nyström approximation; satisfies .
K: Exact kernel matrix .
m: Kernel pproximation degree.
PR: The method provides a low-rank approximation to the kernel matrix
obtained via
Pivoted Cholesky decomposition; satisfies .
eps: Numerical tolerance used for early stopping in the pivoted Cholesky decomposition.
m: Number of random features.
Z: Random Fourier Feature matrix with
so that .
W: Random frequency matrix
(row is ), drawn i.i.d. from the spectral density of the chosen kernel:
Gaussian: (e.g., ).
Laplace: i.i.d.
b Random phase vector , i.i.d. .
# Data setting set.seed(1) lambda = 1e-4 d = 1 rho = 1 n = 50 X = matrix(runif(n*d, 0, 1), nrow = n, ncol = d) y = as.vector(sin(2*pi*rowMeans(X)^3) + rnorm(n, 0, 0.1)) # Exapmle: pivoted cholesky model = fastkrr(X, y, kernel = "gaussian", opt = "pivoted", rho = rho, lambda = 1e-4) # Example: nystrom model = fastkrr(X, y, kernel = "gaussian", opt = "nystrom", rho = rho, lambda = 1e-4) # Example: random fourier features model = fastkrr(X, y, kernel = "gaussian", opt = "rff", rho = rho, lambda = 1e-4) # Example: Laplace kernel model = fastkrr(X, y, kernel = "laplace", opt = "nystrom", n_threads = 1, rho = rho)# Data setting set.seed(1) lambda = 1e-4 d = 1 rho = 1 n = 50 X = matrix(runif(n*d, 0, 1), nrow = n, ncol = d) y = as.vector(sin(2*pi*rowMeans(X)^3) + rnorm(n, 0, 0.1)) # Exapmle: pivoted cholesky model = fastkrr(X, y, kernel = "gaussian", opt = "pivoted", rho = rho, lambda = 1e-4) # Example: nystrom model = fastkrr(X, y, kernel = "gaussian", opt = "nystrom", rho = rho, lambda = 1e-4) # Example: random fourier features model = fastkrr(X, y, kernel = "gaussian", opt = "rff", rho = rho, lambda = 1e-4) # Example: Laplace kernel model = fastkrr(X, y, kernel = "laplace", opt = "nystrom", n_threads = 1, rho = rho)
Defines a Kernel Ridge Regression model specification
for use with the tidymodels ecosystem via parsnip. This spec can be
paired with the "fastkrr" engine implemented in this package to fit
exact or kernel approximation (Nyström, Pivoted Cholesky, Random Fourier Features) within
recipes/workflows pipelines.
krr_reg( mode = "regression", kernel = NULL, opt = NULL, eps = NULL, n_threads = NULL, m = NULL, rho = NULL, penalty = NULL, fastcv = NULL )krr_reg( mode = "regression", kernel = NULL, opt = NULL, eps = NULL, n_threads = NULL, m = NULL, rho = NULL, penalty = NULL, fastcv = NULL )
mode |
A single string; only '"regression"' is supported. |
kernel |
Kernel matrix |
opt |
Method for constructing or approximating :
|
eps |
Tolerance parameter used only in |
n_threads |
Number of parallel threads. It is applied only for
|
m |
Approximation rank(number of random features) used for the low-rank kernel approximation. |
rho |
Scaling parameter of the kernel( |
penalty |
Regularization parameter. |
fastcv |
If |
A parsnip model specification of class "krr_reg".
if (all(vapply( c("parsnip","stats","modeldata"), requireNamespace, quietly = TRUE, FUN.VALUE = logical(1) ))) { library(tidymodels) library(parsnip) library(stats) library(modeldata) # Data analysis data(ames) ames = ames %>% mutate(Sale_Price = log10(Sale_Price)) set.seed(502) ames_split = initial_split(ames, prop = 0.80, strata = Sale_Price) ames_train = training(ames_split) # dim (2342, 74) ames_test = testing(ames_split) # dim (588, 74) # Model spec krr_spec = krr_reg(kernel = "gaussian", opt = "exact", m = 50, eps = 1e-6, n_threads = 4, rho = 1, penalty = tune()) %>% set_engine("fastkrr") %>% set_mode("regression") # Define rec rec = recipe(Sale_Price ~ Longitude + Latitude, data = ames_train) # workflow wf = workflow() %>% add_recipe(rec) %>% add_model(krr_spec) # Define hyper-parameter grid param_grid = grid_regular( dials::penalty(range = c(-10, -3)), levels = 5 ) # CV setting set.seed(123) cv_folds = vfold_cv(ames_train, v = 5, strata = Sale_Price) # Tuning tune_results = tune_grid( wf, resamples = cv_folds, grid = param_grid, metrics = metric_set(rmse), control = control_grid(verbose = TRUE, save_pred = TRUE) ) # Result check collect_metrics(tune_results) # Select best parameter best_params = select_best(tune_results, metric = "rmse") # Finalized model spec using best parameter final_spec = finalize_model(krr_spec, best_params) final_wf = workflow() %>% add_recipe(rec) %>% add_model(final_spec) # Finalized fitting using best parameter final_fit = final_wf %>% fit(data = ames_train) # Prediction predict(final_fit, new_data = ames_test) print(best_params) }if (all(vapply( c("parsnip","stats","modeldata"), requireNamespace, quietly = TRUE, FUN.VALUE = logical(1) ))) { library(tidymodels) library(parsnip) library(stats) library(modeldata) # Data analysis data(ames) ames = ames %>% mutate(Sale_Price = log10(Sale_Price)) set.seed(502) ames_split = initial_split(ames, prop = 0.80, strata = Sale_Price) ames_train = training(ames_split) # dim (2342, 74) ames_test = testing(ames_split) # dim (588, 74) # Model spec krr_spec = krr_reg(kernel = "gaussian", opt = "exact", m = 50, eps = 1e-6, n_threads = 4, rho = 1, penalty = tune()) %>% set_engine("fastkrr") %>% set_mode("regression") # Define rec rec = recipe(Sale_Price ~ Longitude + Latitude, data = ames_train) # workflow wf = workflow() %>% add_recipe(rec) %>% add_model(krr_spec) # Define hyper-parameter grid param_grid = grid_regular( dials::penalty(range = c(-10, -3)), levels = 5 ) # CV setting set.seed(123) cv_folds = vfold_cv(ames_train, v = 5, strata = Sale_Price) # Tuning tune_results = tune_grid( wf, resamples = cv_folds, grid = param_grid, metrics = metric_set(rmse), control = control_grid(verbose = TRUE, save_pred = TRUE) ) # Result check collect_metrics(tune_results) # Select best parameter best_params = select_best(tune_results, metric = "rmse") # Finalized model spec using best parameter final_spec = finalize_model(krr_spec, best_params) final_wf = workflow() %>% add_recipe(rec) %>% add_model(final_spec) # Finalized fitting using best parameter final_fit = final_wf %>% fit(data = ames_train) # Prediction predict(final_fit, new_data = ames_test) print(best_params) }
construction for given datasetsConstructs a kernel matrix given two
datasets and ,
where and denote the
i-th and j-th rows of and , respectively, and
for a user-specified kernel.
Implemented in C++ via RcppArmadillo.
X |
Design matrix |
X_new |
Second matrix |
kernel |
Kernel type; one of |
rho |
Kernel width parameter ( |
n_threads |
Number of parallel threads.
The default is 4. If the system does not support 4 threads,
it automatically falls back to 1 thread.
Parallelization (implemented in C++) is one of the main advantages
of this package and is applied only for |
Gaussian:
Laplace:
An S3 object of class "kernel_matrix" that represents the computed
kernel matrix. If X_new is NULL, the result is a symmetric matrix
, with .
Otherwise, the result is a rectangular matrix
, with .
# Data setting set.seed(1) d = 1 rho = 1 n = 1000 X = matrix(runif(n*d, 0, 1), nrow = n, ncol = d) # New design matrix new_n = 1500 new_X = matrix(runif(new_n*d, 0, 1), nrow = new_n, ncol = d) # Make kernel : Gaussian kernel K = make_kernel(X, kernel = "gaussian", rho = rho) ## symmetric matrix new_K = make_kernel(X, new_X, kernel = "gaussian", rho = rho) ## rectangular matrix # Make kernel : Laplace kernel K = make_kernel(X, kernel = "laplace", rho = rho, n_threads = 1) ## symmetric matrix new_K = make_kernel(X, new_X, kernel = "laplace", rho = rho, n_threads = 1) ## rectangular matrix# Data setting set.seed(1) d = 1 rho = 1 n = 1000 X = matrix(runif(n*d, 0, 1), nrow = n, ncol = d) # New design matrix new_n = 1500 new_X = matrix(runif(new_n*d, 0, 1), nrow = new_n, ncol = d) # Make kernel : Gaussian kernel K = make_kernel(X, kernel = "gaussian", rho = rho) ## symmetric matrix new_K = make_kernel(X, new_X, kernel = "gaussian", rho = rho) ## rectangular matrix # Make kernel : Laplace kernel K = make_kernel(X, kernel = "laplace", rho = rho, n_threads = 1) ## symmetric matrix new_K = make_kernel(X, new_X, kernel = "laplace", rho = rho, n_threads = 1) ## rectangular matrix
"param()" is a generic S3 function that displays (and invisibly returns)
model hyperparameters. Methods are provided for "krr" objects.
param(x, ...) ## Default S3 method: param(x, ...)param(x, ...) ## Default S3 method: param(x, ...)
x |
An object. |
... |
Additional arguments passed to methods. |
A named list of hyperparameters (invisibly); may print side effects.
Displays (and invisibly returns) the hyperparameters actually used by a
fitted object. For krr objects returned by fastkrr,
this prints a concise hyperparameter panel (e.g., rho, lambda,
m, eps, n_threads, d).
## S3 method for class 'krr' param(x, ...)## S3 method for class 'krr' param(x, ...)
x |
An object of class |
... |
Additional arguments. |
Pivoted approximation note: When opt = "pivoted", the
effective number of pivots m used during the approximation may be
smaller than the user-specified m because the algorithm can stop
early based on eps. If you want to confirm the initial m
that you set, please see the printed Call (the original function
call shows your input arguments).
Prints a human-readable panel and returns (invisibly) a named list of hyperparameters.
# Data setting set.seed(1) lambda = 1e-4 d = 1 n = 50 rho = 1 X = matrix(runif(n*d, 0, 1), nrow = n, ncol = d) y = as.vector(sin(2*pi*rowMeans(X)^3) + rnorm(n, 0, 0.1)) model = fastkrr(X, y, kernel="gaussian", opt="nystrom", rho=1, lambda=1e-4, m=200, n_threads=4, fastcv=FALSE) class(model) param(model)# Data setting set.seed(1) lambda = 1e-4 d = 1 n = 50 rho = 1 X = matrix(runif(n*d, 0, 1), nrow = n, ncol = d) y = as.vector(sin(2*pi*rowMeans(X)^3) + rnorm(n, 0, 0.1)) model = fastkrr(X, y, kernel="gaussian", opt="nystrom", rho=1, lambda=1e-4, m=200, n_threads=4, fastcv=FALSE) class(model) param(model)
Visualizes fitted results from a Kernel Ridge Regression (KRR) model. Automatically generates predictions on a regular grid (120% of training sample size) and overlays them with training data.
## S3 method for class 'krr' plot(x, show_points = TRUE, ...)## S3 method for class 'krr' plot(x, show_points = TRUE, ...)
x |
A fitted KRR model (class |
show_points |
Logical; if |
... |
Additional arguments (currently ignored). |
For multivariate inputs (), visualization requires fixing
all but one variable. For example, in 2D, one can plot
to examine the effect of
while holding at its mean.
A ggplot object showing the fitted regression curve.
set.seed(1) n = 1000 rho = 1 X = runif(n, 0, 1) y = sin(2*pi*X^3) + rnorm(n, 0, 0.1) model_exact = fastkrr(X, y, kernel = "gaussian", rho = rho, opt = "exact", verbose = FALSE) plot(model_exact)set.seed(1) n = 1000 rho = 1 X = runif(n, 0, 1) y = sin(2*pi*X^3) + rnorm(n, 0, 0.1) model_exact = fastkrr(X, y, kernel = "gaussian", rho = rho, opt = "exact", verbose = FALSE) plot(model_exact)
Generates predictions from a fitted Kernel Ridge Regression (KRR) model for new data.
## S3 method for class 'krr' predict(object, newdata, ...)## S3 method for class 'krr' predict(object, newdata, ...)
object |
A S3 object of class |
newdata |
New design matrix or data frame containing new observations
for which predictions are to be made. If |
... |
Additional arguments (currently ignored). |
A numeric vector of predicted values corresponding to newdata or fitted values.
# Data setting n = 30 d = 1 X = matrix(runif(n*d, 0, 1), nrow = n, ncol = d) y = as.vector(sin(2*pi*rowMeans(X)^3) + rnorm(n, 0, 0.1)) lambda = 1e-4 rho = 1 # Fitting model: pivoted model = fastkrr(X, y, kernel = "gaussian", rho = rho, lambda = lambda, opt = "pivoted") # Predict new_n = 50 new_x = matrix(runif(new_n*d, 0, 1), nrow = new_n, ncol = d) new_y = as.vector(sin(2*pi*rowMeans(new_x)^3) + rnorm(new_n, 0, 0.1)) pred = predict(model, new_x) crossprod(pred - new_y) / new_n predict(model) == attributes(model)$fitted.values# Data setting n = 30 d = 1 X = matrix(runif(n*d, 0, 1), nrow = n, ncol = d) y = as.vector(sin(2*pi*rowMeans(X)^3) + rnorm(n, 0, 0.1)) lambda = 1e-4 rho = 1 # Fitting model: pivoted model = fastkrr(X, y, kernel = "gaussian", rho = rho, lambda = lambda, opt = "pivoted") # Predict new_n = 50 new_x = matrix(runif(new_n*d, 0, 1), nrow = new_n, ncol = d) new_y = as.vector(sin(2*pi*rowMeans(new_x)^3) + rnorm(new_n, 0, 0.1)) pred = predict(model, new_x) crossprod(pred - new_y) / new_n predict(model) == attributes(model)$fitted.values
Displays the approximated kernel matrix and key options used to construct it.
## S3 method for class 'approx_kernel' print(x, ...)## S3 method for class 'approx_kernel' print(x, ...)
x |
An S3 object created by |
... |
Additional arguments (currently ignored). |
The function prints the stored approximated kernel matrix (top-left 6x6)
and summarizes options such as the approximation method (opt), approximaion degree (m),
numerical tolerance (eps), and number of threads used (n_threads).
An approximated kernel matrix and its associated options.
approx_kernel, print.krr, print.kernel_matrix
# Data setting set.seed(1) d = 1 n = 1000 m = 50 rho = 1 X = matrix(runif(n*d, 0, 1), nrow = n, ncol = d) y = as.vector(sin(2*pi*rowMeans(X)^3) + rnorm(n, 0, 0.1)) K = make_kernel(X, kernel = "gaussian", rho = rho) # Example: nystrom K_nystrom = approx_kernel(K = K, opt = "nystrom", m = m, d = d, rho = rho, n_threads = 1) class(K_nystrom) print(K_nystrom)# Data setting set.seed(1) d = 1 n = 1000 m = 50 rho = 1 X = matrix(runif(n*d, 0, 1), nrow = n, ncol = d) y = as.vector(sin(2*pi*rowMeans(X)^3) + rnorm(n, 0, 0.1)) K = make_kernel(X, kernel = "gaussian", rho = rho) # Example: nystrom K_nystrom = approx_kernel(K = K, opt = "nystrom", m = m, d = d, rho = rho, n_threads = 1) class(K_nystrom) print(K_nystrom)
Displays the top-left 6×6 portion of a kernel or approximated kernel matrix for quick inspection.
## S3 method for class 'kernel_matrix' print(x, ...)## S3 method for class 'kernel_matrix' print(x, ...)
x |
An object of class |
... |
Additional arguments (currently ignored). |
A top-left 6x6 block of the kernel matrix to the console.
approx_kernel, fastkrr, print.approx_kernel, print.krr
# data setting set.seed(1) n = 1000 ; d = 1 m = 100 rho = 1 X = matrix(runif(n*d, 0, 1), nrow = n, ncol = d) y = as.vector(sin(2*pi*X^3) + rnorm(n, 0, 0.1)) # Example for fastkrr fit_pivoted = fastkrr(X, y, kernel = "gaussian", opt = "pivoted", m = 100, fastcv = TRUE, verbose = FALSE) class(attr(fit_pivoted, "K")) print(class(attr(fit_pivoted, "K"))) class(attr(fit_pivoted, "K_approx")) print(class(attr(fit_pivoted, "K_approx"))) # Example for make_kernel K = make_kernel(X, kernel = "gaussian", rho = rho) class(K) print(K) # Example for make_kernel K_rff = approx_kernel(X = X, opt = "rff", kernel = "gaussian", d = d, rho = rho, n_threads = 1, m = 100) class(attr(K_rff, "K_approx")) print(attr(K_rff, "K_approx"))# data setting set.seed(1) n = 1000 ; d = 1 m = 100 rho = 1 X = matrix(runif(n*d, 0, 1), nrow = n, ncol = d) y = as.vector(sin(2*pi*X^3) + rnorm(n, 0, 0.1)) # Example for fastkrr fit_pivoted = fastkrr(X, y, kernel = "gaussian", opt = "pivoted", m = 100, fastcv = TRUE, verbose = FALSE) class(attr(fit_pivoted, "K")) print(class(attr(fit_pivoted, "K"))) class(attr(fit_pivoted, "K_approx")) print(class(attr(fit_pivoted, "K_approx"))) # Example for make_kernel K = make_kernel(X, kernel = "gaussian", rho = rho) class(K) print(K) # Example for make_kernel K_rff = approx_kernel(X = X, opt = "rff", kernel = "gaussian", d = d, rho = rho, n_threads = 1, m = 100) class(attr(K_rff, "K_approx")) print(attr(K_rff, "K_approx"))
Displays key information from a fitted Kernel Ridge Regression (KRR) model, including the original call, first few coefficients, a 6×6 block of the kernel (or approximated kernel) matrix, and the main kernel options.
## S3 method for class 'krr' print(x, ...)## S3 method for class 'krr' print(x, ...)
x |
An S3 object of class |
... |
Additional arguments (currently ignored). |
A human-readable summary of the fitted KRR model to the console.
fastkrr, print.approx_kernel,
print.kernel_matrix
# Data setting set.seed(1) lambda = 1e-4 d = 1 n = 50 rho = 1 X = matrix(runif(n*d, 0, 1), nrow = n, ncol = d) y = as.vector(sin(2*pi*rowMeans(X)^3) + rnorm(n, 0, 0.1)) # Example: exact model = fastkrr(X, y, kernel = "gaussian", opt = "exact", rho = rho, lambda = 1e-4) class(model) print(model)# Data setting set.seed(1) lambda = 1e-4 d = 1 n = 50 rho = 1 X = matrix(runif(n*d, 0, 1), nrow = n, ncol = d) y = as.vector(sin(2*pi*rowMeans(X)^3) + rnorm(n, 0, 0.1)) # Example: exact model = fastkrr(X, y, kernel = "gaussian", opt = "exact", rho = rho, lambda = 1e-4) class(model) print(model)
Displays key information from a fitted Kernel Ridge Regression (KRR) model, including the original call, first few coefficients, a 6×6 block of the kernel (or approximated kernel) matrix, and the main kernel options.
## S3 method for class 'krr' summary(object, ...)## S3 method for class 'krr' summary(object, ...)
object |
An S3 object of class |
... |
Additional arguments (currently ignored). |
A human-readable summary of the fitted KRR model to the console.
# Data setting set.seed(1) lambda = 1e-4 d = 1 n = 50 rho = 1 X = matrix(runif(n*d, 0, 1), nrow = n, ncol = d) y = as.vector(sin(2*pi*rowMeans(X)^3) + rnorm(n, 0, 0.1)) # Example: exact model = fastkrr(X, y, kernel = "gaussian", opt = "exact", rho = rho, lambda = 1e-4) class(model) summary(model)# Data setting set.seed(1) lambda = 1e-4 d = 1 n = 50 rho = 1 X = matrix(runif(n*d, 0, 1), nrow = n, ncol = d) y = as.vector(sin(2*pi*rowMeans(X)^3) + rnorm(n, 0, 0.1)) # Example: exact model = fastkrr(X, y, kernel = "gaussian", opt = "exact", rho = rho, lambda = 1e-4) class(model) summary(model)
"krr_reg"
Supplies a tibble of tunable arguments for "krr_reg()".
## S3 method for class 'krr_reg' tunable(x, ...)## S3 method for class 'krr_reg' tunable(x, ...)
x |
A |
... |
Not used; included for S3 method compatibility. |
A tibble (one row per tunable parameter)
with columns "name", "call_info", "source",
"component", and "component_id".