Package 'FastKRR'

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

Help Index


Kernel Ridge Regression using the RcppArmadillo Package

Description

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.

Directory structure

  • 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.

Author(s)

Maintainer: Kwan-Young Bak [email protected] (ORCID) (Sungshin Women's University) [copyright holder]

Authors:

See Also

CVST, Rcpp, RcppArmadillo, parsnip, tidymodels


Compute low-rank approximations(Nyström, Pivoted Cholesky, RFF)

Description

Computes low-rank kernel approximation K~Rn×n\tilde{K} \in \mathbb{R}^{n \times n}using three methods: Nyström approximation, Pivoted Cholesky decomposition, and Random Fourier Features (RFF).

Usage

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
)

Arguments

K

Exact Kernel matrix KRn×nK \in \mathbb{R}^{n \times n}. Used in "nystrom" and "pivoted".

X

Design matrix XRn×dX \in \mathbb{R}^{n \times d}. Only required for "rff".

opt

Method for constructing or approximating :

"nystrom"

Construct a low-rank approximation of the kernel matrix KRn×nK \in \mathbb{R}^{n \times n} using the Nyström approximation.

"pivoted"

Construct a low-rank approximation of the kernel matrix KRn×nK \in \mathbb{R}^{n \times n} using Pivoted Cholesky decomposition.

"rff"

Construct a low-rank approximation of the kernel matrix KRn×nK \in \mathbb{R}^{n \times n} using Random Fourier Features (RFF).

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

nlog(d+5)/10\lceil n \cdot \log(d + 5) / 10 \rceil

where XX is design matrix, n=nrow(X)n = nrow(X) and d=ncol(X)d = ncol(X).

d

Design matrix's dimension (d=ncol(X)d = ncol(X)).

rho

Scaling parameter of the kernel (ρ\rho), specified by the user.

eps

Tolerance parameter used only in "pivoted" for stopping criterion of the Pivoted Cholesky decomposition.

W

Random frequency matrix ωRm×d\omega \in \mathbb{R}^{m \times d}

b

Random phase vector bRmb \in \mathbb{R}^m, i.i.d. Unif[0,2π]\mathrm{Unif} [ 0,\,2\pi ].

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 opt = "nystrom" or opt = "rff" , and for the Laplace kernel (kernel = "laplace").

Details

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 nlog(d+5)/10\lceil n \cdot \log(d + 5) / 10 \rceil.

  • 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 ωRm×d\omega \in \mathbb{R}^{m \times d}) and b (random phase vector bRmb \in \mathbb{R}^{m}).

  • If the user provides them manually, both W and b must be specified and their dimensions must be compatible.

Value

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: n×nn \times n 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: m×dm \times d Random frequency matrix.

  • b: Random phase mm–vector.

  • used_supplied_Wb: Logical; TRUE if user-supplied W, b were used, FALSE otherwise.

  • n_threads: Number of threads used in the computation.

Examples

# 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)

Coef method for fitted Kernel Ridge Regression models

Description

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 α\alpha; for opt = "rff", they represent the coefficient (β\beta).

Usage

## S3 method for class 'krr'
coef(object, ...)

Arguments

object

An S3 object of class krr, typically returned by fastkrr.

...

Additional arguments (currently ignored).

Value

A human-readable summary of the fitted KRR model to the console.

See Also

fastkrr

Examples

# 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)

Compute Model Error (Generic)

Description

Generic function for computing model error.

Usage

error(x, ...)

## Default S3 method:
error(x, ...)

Arguments

x

An object to compute model error for.

...

Additional arguments passed to methods.

Value

A numeric value or class-specific result.


Compute Model Error for Kernel Ridge Regression Models

Description

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.

Usage

## S3 method for class 'krr'
error(x, ...)

Arguments

x

An object of class "krr", typically returned by fastkrr.

...

Additional arguments (ignored).

Details

This method computes:

MSE=1ni(yiy^i)2\text{MSE} = \frac{1}{n} \sum_i (y_i - \hat{y}_i)^2

where y and fitted.values are stored in the "krr" object attributes.

Value

A numeric value giving the mean squared error (MSE).

See Also

summary.krr, plot.krr, predict.krr

Examples

# 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)

Fit kernel ridge regression using exact or approximate methods

Description

This function performs kernel ridge regression (KRR) in high-dimensional settings. The regularization parameter λ\lambda 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).

Usage

fastkrr(
  x,
  y,
  kernel = "gaussian",
  opt = "exact",
  m = NULL,
  eps = 1e-06,
  rho = 1,
  lambda = NULL,
  fastcv = FALSE,
  n_threads = 4,
  verbose = TRUE
)

Arguments

x

Design matrix XRn×dX \in \mathbb{R}^{n\times d}.

y

Response variable yRny \in \mathbb{R}^{n}.

kernel

Kernel type either "gaussian"or "laplace".

opt

Method for constructing or approximating :

"exact"

Construct the full kernel matrix KRn×nK \in \mathbb{R}^{n\times n} using design martix XX.

"nystrom"

Construct a low-rank approximation of the kernel matrix KRn×nK \in \mathbb{R}^{n \times n} using the Nyström approximation.

"pivoted"

Construct a low-rank approximation of the kernel matrix KRn×nK \in \mathbb{R}^{n \times n} using Pivoted Cholesky decomposition.

"rff"

Use Random Fourier Features to construct a feature map ZRn×mZ \in \mathbb{R}^{n \times m} (with mm random features) so that KZZK \approx Z Z^\top. Here, mm is the number of features.

m

Approximation rank(number of random features) used for the low-rank kernel approximation. If not provided by the user, it defaults to

nlog(d+5)10,\lceil n \cdot \frac{\log(d + 5)}{10} \rceil,

where n=nrow(X)n = nrow(X) and d=ncol(X)d = ncol(X).

eps

Tolerance parameter used only in "pivoted" for stopping criterion of the Pivoted Cholesky decomposition.

rho

Scaling parameter of the kernel(ρ\rho), specified by the user. Defaults to 1.

Gaussian kernel : K(x,x)=exp(ρxx22)\text{Gaussian kernel : } \mathcal{K}(x, x') = \exp(-\rho \| x - x'\|^2_2)

Laplace kernel : K(x,x)=exp(ρxx1)\text{Laplace kernel : } \mathcal{K}(x, x') = \exp(-\rho \| x - x'\|_1)

lambda

Regularization parameter. If NULL, the penalty parameter is chosen automatically via CVST package. If not provided, the argument is set to a kernel-specific grid of 100 values: [1010,103][10^{-10}, 10^{-3}] for Gaussian, [105,102][10^{-5}, 10^{-2}] for Laplace.

fastcv

If TRUE, accelerated cross-validation is performed via sequential testing (early stopping) as implemented in the CVST package. The default is FALSE.

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 opt = "nystrom" or opt = "rff", and for the Laplace kernel (kernel = "laplace").

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.

Details

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

    nlog(d+5)/10\lceil n \cdot \log(d + 5) / 10 \rceil

    where n=nrow(X)n = nrow(X) and d=ncol(X)d = ncol(X). Otherwise, m must be a positive integer.

  • rho must be a positive real number (default is 1).

  • lambda can be specified in three ways:

    1. A positive numeric scalar, in which case the model is fitted with this single value.

    2. 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).

    3. 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.

Value

An S3 object of class "fastkrr", which is a list containing the results of the fitted Kernel Ridge Regression model.

  • coefficients: Estimated coefficient vector Rn\mathbb{R}^{n}. Accessible via model$coefficients.

  • fitted.values: Fitted values Rn\mathbb{R}^{n}. 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:

opt = “exact”

  • K: The full kernel matrix.

opt = “nystrom”

  • K: Exact kernel matrix KRn×nK \in \mathbb{R}^{n \times n}.

  • m: Kernel approximation degree.

  • R: The method provides a low-rank approximation to the kernel matrix RRn×mR \in \mathbb{R}^{n \times m} obtained via Nyström approximation; satisfies KRRK \approx R R^\top.

opt = “pivoted”

  • K: Exact kernel matrix KRn×nK \in \mathbb{R}^{n \times n}.

  • m: Kernel pproximation degree.

  • PR: The method provides a low-rank approximation to the kernel matrix PRRn×mPR \in \mathbb{R}^{ n \times m} obtained via Pivoted Cholesky decomposition; satisfies KPR(PR)K \approx PR\,(PR)^\top.

  • eps: Numerical tolerance used for early stopping in the pivoted Cholesky decomposition.

opt = “rff”

  • m: Number of random features.

  • Z: Random Fourier Feature matrix ZRn×mZ \in \mathbb{R}^{n \times m} with Zij=zj(xi)=2/mcos(ωjxi+bj),j=1,,m,Z_{ij} = z_j(x_i) = \sqrt{2/m}\cos(\omega_j^\top x_i + b_j), \quad j = 1, \cdots, m, so that KZZK \approx Z Z^\top.

  • W: Random frequency matrix ωRm×d\omega \in \mathbb{R}^{m \times d} (row jj is ωjRd\omega_j^\top \in \mathbb{R}^d), drawn i.i.d. from the spectral density of the chosen kernel:

    • Gaussian: ωjkN(0,2γ)\omega_{jk} \sim \mathcal{N}(0, 2\gamma) (e.g., γ=1/2\gamma=1/\ell^2).

    • Laplace: ωjkCauchy(0,1/σ)\omega_{jk} \sim \mathrm{Cauchy}(0, 1/\sigma) i.i.d.

  • b Random phase vector bRmb \in \mathbb{R}^m, i.i.d. Unif[0,2π]\mathrm{Unif}[0,\,2\pi].

Examples

# 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)

Kernel Ridge Regression

Description

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.

Usage

krr_reg(
  mode = "regression",
  kernel = NULL,
  opt = NULL,
  eps = NULL,
  n_threads = NULL,
  m = NULL,
  rho = NULL,
  penalty = NULL,
  fastcv = NULL
)

Arguments

mode

A single string; only '"regression"' is supported.

kernel

Kernel matrix KK has two kinds of Kernel ("gaussian", "laplace").

opt

Method for constructing or approximating :

"exact"

Construct the full kernel matrix KRn×nK \in \mathbb{R}^{n\times n} using design matrix XX.

"nystrom"

Construct a low-rank approximation of the kernel matrix KRn×nK \in \mathbb{R}^{n \times n} using the Nyström approximation.

"pivoted"

Construct a low-rank approximation of the kernel matrix KRn×nK \in \mathbb{R}^{n \times n} using Pivoted Cholesky decomposition.

"rff"

Use Random Fourier Features to construct a feature map ZRn×mZ \in \mathbb{R}^{n \times m} (with mm random features) so that KZZK \approx Z Z^\top. Here, mm is the number of features.

eps

Tolerance parameter used only in "pivoted" for stopping criterion of the Pivoted Cholesky decomposition.

n_threads

Number of parallel threads. It is applied only for opt = "nystrom" or opt = "rff", and for the Laplace kernel (kernel = "laplace").

m

Approximation rank(number of random features) used for the low-rank kernel approximation.

rho

Scaling parameter of the kernel(ρ\rho).

penalty

Regularization parameter.

fastcv

If TRUE, accelerated cross-validation is performed via sequential testing (early stopping) as implemented in the CVST package.

Value

A parsnip model specification of class "krr_reg".

Examples

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)

}

Kernel matrix KK construction for given datasets

Description

Constructs a kernel matrix KRn×nK \in \mathbb{R}^{n \times n'} given two datasets XRn×dX \in \mathbb{R}^{n \times d} and XRn×dX' \in \mathbb{R}^{n' \times d}, where xiRdx_i \in \mathbb{R}^d and xjRdx'_j \in \mathbb{R}^d denote the i-th and j-th rows of XX and XX', respectively, and Kij=K(xi,xj)K_{ij}=\mathcal{K}(x_i, x'_j) for a user-specified kernel. Implemented in C++ via RcppArmadillo.

Arguments

X

Design matrix XRn×dX \in \mathbb{R}^{n \times d} (rows xiRdx_i \in \mathbb{R}^d).

X_new

Second matrix XRn×dX' \in \mathbb{R}^{n' \times d} (rows xjRdx'_j \in \mathbb{R}^d). If omitted, X=XX' = X and n=nn' = n.

kernel

Kernel type; one of "gaussian" or "laplace".

rho

Kernel width parameter (ρ>0\rho > 0).

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 "laplace" kernels.

Details

Gaussian:

K(xi,xj)=exp ⁣(ρxixj22)\mathcal{K}(x_i,x_j)=\exp\!\big(-\rho\|x_i-x_j\|_2^2\big)

Laplace:

K(xi,xj)=exp ⁣(ρxixj1)\mathcal{K}(x_i,x_j)=\exp\!\big(-\rho\|x_i-x_j\|_1\big)

Value

An S3 object of class "kernel_matrix" that represents the computed kernel matrix. If X_new is NULL, the result is a symmetric matrix Kij=K(xi,xj)K_{ij} = \mathcal{K}(x_i, x_j), with KRn×nK \in \mathbb{R}^{n \times n}. Otherwise, the result is a rectangular matrix Kij=K(xi,xj)K'_{ij} = \mathcal{K}(x_i, x'_j), with KRn×nK' \in \mathbb{R}^{n \times n'}.

Examples

# 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

Extract/print hyperparameters of fitted models

Description

"param()" is a generic S3 function that displays (and invisibly returns) model hyperparameters. Methods are provided for "krr" objects.

Usage

param(x, ...)

## Default S3 method:
param(x, ...)

Arguments

x

An object.

...

Additional arguments passed to methods.

Value

A named list of hyperparameters (invisibly); may print side effects.


Param method for fitted Kernel Ridge Regression models

Description

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).

Usage

## S3 method for class 'krr'
param(x, ...)

Arguments

x

An object of class "krr", typically returned by fastkrr.

...

Additional arguments.

Details

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).

Value

Prints a human-readable panel and returns (invisibly) a named list of hyperparameters.

See Also

fastkrr

Examples

# 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)

Plot method for fitted Kernel Ridge Regression (KRR) models

Description

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.

Usage

## S3 method for class 'krr'
plot(x, show_points = TRUE, ...)

Arguments

x

A fitted KRR model (class "krr") returned by fastkrr.

show_points

Logical; if TRUE, displays the training data points. Default = TRUE.

...

Additional arguments (currently ignored).

Details

For multivariate inputs (d2d \ge 2), visualization requires fixing all but one variable. For example, in 2D, one can plot f(x1,x2=xˉ2)f(x_1, x_2 = \bar{x}_2) to examine the effect of x1x_1 while holding x2x_2 at its mean.

Value

A ggplot object showing the fitted regression curve.

See Also

fastkrr, predict.krr

Examples

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)

Predict responses for new data using fitted KRR model

Description

Generates predictions from a fitted Kernel Ridge Regression (KRR) model for new data.

Usage

## S3 method for class 'krr'
predict(object, newdata, ...)

Arguments

object

A S3 object of class krr created by fastkrr.

newdata

New design matrix or data frame containing new observations for which predictions are to be made. If newdata is missing, the function returns fitted values.

...

Additional arguments (currently ignored).

Value

A numeric vector of predicted values corresponding to newdata or fitted values.

See Also

fastkrr, make_kernel

Examples

# 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

Print method for approximated kernel matrices

Description

Displays the approximated kernel matrix and key options used to construct it.

Usage

## S3 method for class 'approx_kernel'
print(x, ...)

Arguments

x

An S3 object created by approx_kernel.

...

Additional arguments (currently ignored).

Details

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).

Value

An approximated kernel matrix and its associated options.

See Also

approx_kernel, print.krr, print.kernel_matrix

Examples

# 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)

Print method for kernel matrices

Description

Displays the top-left 6×6 portion of a kernel or approximated kernel matrix for quick inspection.

Usage

## S3 method for class 'kernel_matrix'
print(x, ...)

Arguments

x

An object of class kernel_matrix, which may represent either an exact kernel matrix (from make_kernel or fastkrr) or an approximated kernel matrix (from approx_kernel).

...

Additional arguments (currently ignored).

Value

A top-left 6x6 block of the kernel matrix to the console.

See Also

approx_kernel, fastkrr, print.approx_kernel, print.krr

Examples

# 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"))

Print method for fitted Kernel Ridge Regression models

Description

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.

Usage

## S3 method for class 'krr'
print(x, ...)

Arguments

x

An S3 object of class krr, typically returned by fastkrr.

...

Additional arguments (currently ignored).

Value

A human-readable summary of the fitted KRR model to the console.

See Also

fastkrr, print.approx_kernel, print.kernel_matrix

Examples

# 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)

Summary method for fitted Kernel Ridge Regression models

Description

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.

Usage

## S3 method for class 'krr'
summary(object, ...)

Arguments

object

An S3 object of class krr, typically returned by fastkrr.

...

Additional arguments (currently ignored).

Value

A human-readable summary of the fitted KRR model to the console.

See Also

fastkrr

Examples

# 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)

Expose tunable parameters for "krr_reg"

Description

Supplies a tibble of tunable arguments for "krr_reg()".

Usage

## S3 method for class 'krr_reg'
tunable(x, ...)

Arguments

x

A "krr_reg" model specification.

...

Not used; included for S3 method compatibility.

Value

A tibble (one row per tunable parameter) with columns "name", "call_info", "source", "component", and "component_id".